#include "stdafx.h" #include "TBPImageAlgCls.h" #include #include #include #include #include #include #include "GPUTool.cuh" #include "GPUTBPImage.cuh" #include "ImageOperatorBase.h" void CreatePixelXYZ(std::shared_ptr echoL0ds, QString outPixelXYZPath) { long bandwidth = echoL0ds->getBandwidth(); double Rnear = echoL0ds->getNearRange(); double Rfar = echoL0ds->getFarRange(); double refRange = echoL0ds->getRefPhaseRange(); double dx = LIGHTSPEED / 2.0 / bandwidth; // c/2b // 创建坐标系统 long prfcount = echoL0ds->getPluseCount(); long freqcount = echoL0ds->getPlusePoints(); Eigen::MatrixXd gt = Eigen::MatrixXd::Zero(2, 3); gt(0, 0) = 0; gt(0, 1) = 1; gt(0, 2) = 0; gt(1, 0) = 0; gt(1, 1) = 0; gt(1, 2) = 1; gdalImage xyzRaster = CreategdalImage(outPixelXYZPath, prfcount, freqcount, 3, gt, QString(""), false, true,true,GDT_Float64); std::shared_ptr antpos = echoL0ds->getAntPos(); dx = (echoL0ds->getFarRange()-echoL0ds->getNearRange())/(echoL0ds->getPlusePoints()-1); Rnear = echoL0ds->getNearRange(); double Rref = echoL0ds->getRefPhaseRange(); double centerInc = echoL0ds->getCenterAngle()*d2r; long echocol = Memory1GB * 1.0 / 8 / 4 / prfcount * 6; qDebug() << "echocol:\t " << echocol ; echocol = echocol < 1 ? 1: echocol; std::shared_ptr Pxs((double*)mallocCUDAHost(sizeof(double)*prfcount), FreeCUDAHost); std::shared_ptr Pys((double*)mallocCUDAHost(sizeof(double) * prfcount), FreeCUDAHost); std::shared_ptr Pzs((double*)mallocCUDAHost(sizeof(double) * prfcount), FreeCUDAHost); std::shared_ptr AntDirectX((double*)mallocCUDAHost(sizeof(double) * prfcount), FreeCUDAHost); std::shared_ptr AntDirectY((double*)mallocCUDAHost(sizeof(double) * prfcount), FreeCUDAHost); std::shared_ptr AntDirectZ((double*)mallocCUDAHost(sizeof(double) * prfcount), FreeCUDAHost); { std::shared_ptr antpos = echoL0ds->getAntPos(); double time = 0; double Px = 0; double Py = 0; double Pz = 0; for (long i = 0; i < prfcount; i++) { Pxs.get()[i] = antpos.get()[i * 19 + 1]; // 卫星坐标 Pys.get()[i] = antpos.get()[i * 19 + 2]; Pzs.get()[i] = antpos.get()[i * 19 + 3]; AntDirectX.get()[i] = antpos.get()[i * 19 + 13];// zero doppler AntDirectY.get()[i] = antpos.get()[i * 19 + 14]; AntDirectZ.get()[i] = antpos.get()[i * 19 + 15]; double NormAnt = std::sqrt(AntDirectX.get()[i] * AntDirectX.get()[i] + AntDirectY.get()[i] * AntDirectY.get()[i] + AntDirectZ.get()[i] * AntDirectZ.get()[i]); AntDirectX.get()[i] = AntDirectX.get()[i] / NormAnt; AntDirectY.get()[i] = AntDirectY.get()[i] / NormAnt; AntDirectZ.get()[i] = AntDirectZ.get()[i] / NormAnt;// 归一化 } antpos.reset(); } std::shared_ptr d_Pxs((double*)mallocCUDADevice(sizeof(double) * prfcount), FreeCUDADevice); std::shared_ptr d_Pys((double*)mallocCUDADevice(sizeof(double) * prfcount), FreeCUDADevice); std::shared_ptr d_Pzs((double*)mallocCUDADevice(sizeof(double) * prfcount), FreeCUDADevice); std::shared_ptr d_AntDirectX((double*)mallocCUDADevice(sizeof(double) * prfcount), FreeCUDADevice); std::shared_ptr d_AntDirectY((double*)mallocCUDADevice(sizeof(double) * prfcount), FreeCUDADevice); std::shared_ptr d_AntDirectZ((double*)mallocCUDADevice(sizeof(double) * prfcount), FreeCUDADevice); HostToDevice(Pxs.get(), d_Pxs.get(), sizeof(double) * prfcount); HostToDevice(Pys.get(), d_Pys.get(), sizeof(double) * prfcount); HostToDevice(Pzs.get(), d_Pzs.get(), sizeof(double) * prfcount); HostToDevice(AntDirectX.get(), d_AntDirectX.get(), sizeof(double) * prfcount); HostToDevice(AntDirectY.get(), d_AntDirectY.get(), sizeof(double) * prfcount); HostToDevice(AntDirectZ.get(), d_AntDirectZ.get(), sizeof(double) * prfcount); for (long startcolidx = 0; startcolidx < freqcount; startcolidx = startcolidx + echocol) { long tempechocol = echocol; if (startcolidx + tempechocol >= freqcount) { tempechocol = freqcount - startcolidx; } qDebug() << "\r[" << QDateTime::currentDateTime().toString("yyyy-MM-dd hh:mm:ss.zzz") << "] imgxyz :\t" << startcolidx << "\t-\t" << startcolidx + tempechocol << " / " << freqcount ; std::shared_ptr demx = readDataArr(xyzRaster, 0, startcolidx, prfcount, tempechocol, 1, GDALREADARRCOPYMETHOD::VARIABLEMETHOD); std::shared_ptr demy = readDataArr(xyzRaster, 0, startcolidx, prfcount, tempechocol, 2, GDALREADARRCOPYMETHOD::VARIABLEMETHOD); std::shared_ptr demz = readDataArr(xyzRaster, 0, startcolidx, prfcount, tempechocol, 3, GDALREADARRCOPYMETHOD::VARIABLEMETHOD); std::shared_ptr h_demx((double*)mallocCUDAHost(sizeof(double) * prfcount*tempechocol), FreeCUDAHost); std::shared_ptr h_demy((double*)mallocCUDAHost(sizeof(double) * prfcount*tempechocol), FreeCUDAHost); std::shared_ptr h_demz((double*)mallocCUDAHost(sizeof(double) * prfcount*tempechocol), FreeCUDAHost); #pragma omp parallel for for (long ii = 0; ii < prfcount; ii++) { for (long jj = 0; jj < tempechocol; jj++) { h_demx.get()[ii*tempechocol+jj]=demx.get()[ii*tempechocol+jj]; h_demy.get()[ii*tempechocol+jj]=demy.get()[ii*tempechocol+jj]; h_demz.get()[ii*tempechocol+jj]=demz.get()[ii*tempechocol+jj]; } } std::shared_ptr d_demx((double*)mallocCUDADevice(sizeof(double) * prfcount * tempechocol), FreeCUDADevice); std::shared_ptr d_demy((double*)mallocCUDADevice(sizeof(double) * prfcount * tempechocol), FreeCUDADevice); std::shared_ptr d_demz((double*)mallocCUDADevice(sizeof(double) * prfcount * tempechocol), FreeCUDADevice); HostToDevice(h_demx.get(), d_demx.get(), sizeof(double) * prfcount * tempechocol); HostToDevice(h_demy.get(), d_demy.get(), sizeof(double) * prfcount * tempechocol); HostToDevice(h_demz.get(), d_demz.get(), sizeof(double) * prfcount * tempechocol); TIMEBPCreateImageGrid( d_Pxs.get(), d_Pys.get(), d_Pzs.get(), d_AntDirectX.get(), d_AntDirectY.get(), d_AntDirectZ.get(), d_demx.get(), d_demy.get(), d_demz.get(), prfcount, tempechocol, 1000, Rnear+dx* startcolidx, dx, refRange // 更新最近修读 ); DeviceToHost(h_demx.get(), d_demx.get(), sizeof(double) * prfcount * tempechocol); DeviceToHost(h_demy.get(), d_demy.get(), sizeof(double) * prfcount * tempechocol); DeviceToHost(h_demz.get(), d_demz.get(), sizeof(double) * prfcount * tempechocol); #pragma omp parallel for for (long ii = 0; ii < prfcount; ii++) { for (long jj = 0; jj < tempechocol; jj++) { demx.get()[ii * tempechocol + jj]=h_demx.get()[ii * tempechocol + jj] ; demy.get()[ii * tempechocol + jj]=h_demy.get()[ii * tempechocol + jj] ; demz.get()[ii * tempechocol + jj]=h_demz.get()[ii * tempechocol + jj] ; } } xyzRaster.saveImage(demx, 0, startcolidx,prfcount,tempechocol, 1); xyzRaster.saveImage(demy, 0, startcolidx,prfcount,tempechocol, 2); xyzRaster.saveImage(demz, 0, startcolidx,prfcount,tempechocol, 3); } } void TBPImageProcess(QString echofile, QString outImageFolder, QString imagePlanePath,long num_thread) { std::shared_ptr echoL0ds(new EchoL0Dataset); echoL0ds->Open(echofile); std::shared_ptr< SARSimulationImageL1Dataset> imagL1(new SARSimulationImageL1Dataset); imagL1->setCenterAngle(echoL0ds->getCenterAngle()); imagL1->setCenterFreq(echoL0ds->getCenterFreq()); imagL1->setNearRange(echoL0ds->getNearRange()); imagL1->setRefRange((echoL0ds->getNearRange() + echoL0ds->getFarRange()) / 2); imagL1->setFarRange(echoL0ds->getFarRange()); imagL1->setFs(echoL0ds->getFs()); imagL1->setLookSide(echoL0ds->getLookSide()); imagL1->OpenOrNew(outImageFolder, echoL0ds->getSimulationTaskName(), echoL0ds->getPluseCount(), echoL0ds->getPlusePoints()); TBPImageAlgCls TBPimag; TBPimag.setEchoL0(echoL0ds); TBPimag.setImageL1(imagL1); TBPimag.setImagePlanePath(imagePlanePath); TBPimag.Process(num_thread); } void TBPImageAlgCls::setImagePlanePath(QString INimagePlanePath) { this->imagePlanePath = INimagePlanePath; } QString TBPImageAlgCls::getImagePlanePath() { return this->imagePlanePath; } void TBPImageAlgCls::setEchoL0(std::shared_ptr inL0ds) { this->L0ds = inL0ds; } void TBPImageAlgCls::setImageL1(std::shared_ptr inL1ds) { this->L1ds = inL1ds; } std::shared_ptr TBPImageAlgCls::getEchoL1() { return this->L0ds; } std::shared_ptr TBPImageAlgCls::getImageL0() { return this->L1ds; } ErrorCode TBPImageAlgCls::Process(long num_thread) { qDebug() << u8"开始成像"; qDebug() << u8"创建成像平面的XYZ"; QString outRasterXYZ = JoinPath(this->L1ds->getoutFolderPath(), this->L0ds->getSimulationTaskName() + "_xyz.bin"); CreatePixelXYZ(this->L0ds, outRasterXYZ); this->outRasterXYZPath = outRasterXYZ; qDebug() << u8"频域回波-> 时域回波"; this->TimeEchoDataPath = JoinPath(this->L1ds->getoutFolderPath(), this->L0ds->getSimulationTaskName() + "_Timeecho.bin"); this->EchoFreqToTime(); // 初始化Raster qDebug() << u8"初始化影像"; long imageheight = this->L1ds->getrowCount(); long imagewidth = this->L1ds->getcolCount(); long blokline = Memory1GB / 8 / 4 / imageheight * 32; blokline = blokline < 1000 ? 1000 : blokline; long startline = 0; for (startline = 0; startline < imageheight; startline = startline + blokline) { long templine = blokline; if (startline + templine >= imageheight) { templine = imageheight - startline; } qDebug() << "\r[" << QDateTime::currentDateTime().toString("yyyy-MM-dd hh:mm:ss.zzz") << "] imgxyz :\t" << startline << "\t-\t" << startline + templine << " / " << imageheight ; std::shared_ptr> imageRaster = this->L1ds->getImageRaster(startline, templine); for (long i = 0; i < templine; i++) { for (long j = 0; j < imagewidth; j++) { imageRaster.get()[i * imagewidth + j] = std::complex(0,0); } } this->L1ds->saveImageRaster(imageRaster, startline,templine); } qDebug() << u8"频域回波-> 时域回波 结束"; if (GPURUN) { return this->ProcessGPU(); } else { QMessageBox::information(nullptr,u8"提示",u8"目前只支持显卡"); return ErrorCode::FAIL; } return ErrorCode::SUCCESS; } ErrorCode TBPImageAlgCls::ProcessGPU() { // 常用参数 long rowCount = this->L1ds->getrowCount(); long colCount = this->L1ds->getcolCount(); long pixelCount = rowCount * colCount; long PRFCount = this->L0ds->getPluseCount(); long PlusePoints = this->L0ds->getPlusePoints(); long bandwidth = this->L0ds->getBandwidth(); double Rnear = this->L1ds->getNearRange(); double Rfar = this->L1ds->getFarRange(); double refRange = this->L0ds->getRefPhaseRange(); double dx = LIGHTSPEED / 2.0 / bandwidth; // c/2b Rfar = Rnear + dx * (PlusePoints - 1); // 更新斜距距离 float freq = this->L1ds->getCenterFreq(); double factorj = freq * 4 * M_PI / LIGHTSPEED ; qDebug() << "------------------------------------------------"; qDebug() << "TBP params:"; qDebug() << "Rnear:\t" << Rnear; qDebug() << "Rfar:\t" << Rfar; qDebug() << "refRange:\t" << this->getEchoL1()->getRefPhaseRange(); qDebug() << "dx:\t" << dx; qDebug() << "freq:\t" << freq; qDebug() << "rowCount:\t" << rowCount; qDebug() << "colCount:\t" << colCount; qDebug() << "PRFCount:\t" << PRFCount; qDebug() << "PlusePoints:\t" << PlusePoints; qDebug() << "bandwidth:\t" << bandwidth; // 反方向计算起始相位 double deltaF = bandwidth / (PlusePoints - 1); double startfreq = freq - bandwidth / 2; double startlamda = LIGHTSPEED / startfreq; qDebug() << "deltaF:\t" << deltaF; std::shared_ptr Pxs (new double[this->L0ds->getPluseCount()],delArrPtr); std::shared_ptr Pys (new double[this->L0ds->getPluseCount()],delArrPtr); std::shared_ptr Pzs (new double[this->L0ds->getPluseCount()],delArrPtr); { std::shared_ptr antpos = this->L0ds->getAntPos(); double time = 0; double Px = 0; double Py = 0; double Pz = 0; for (long i = 0; i < rowCount; i++) { time = antpos.get()[i *19 + 0]; Px = antpos.get()[i *19 + 1]; Py = antpos.get()[i *19 + 2]; Pz = antpos.get()[i *19 + 3]; Pxs.get()[i] = Px; Pys.get()[i] = Py; Pzs.get()[i] = Pz; } antpos.reset(); } // 计算成像的基本参数 // 距离向分辨率 double dr = sqrt(pow(Pxs.get()[2]- Pxs.get()[1],2)+pow(Pys.get()[2] - Pys.get()[1],2)+pow(Pzs.get()[2] - Pzs.get()[1],2)); qDebug() << "------- resolution ----------------------------------"; qDebug() << "Range Resolution (m):\t" << dx ; qDebug() << "Cross Resolution (m):\t" << dr; qDebug() << "Range Range (m):\t" << dx*PlusePoints; qDebug() << "Cross Range (m):\t" << dr*PRFCount; qDebug() << "start Freq (Hz):\t" << startfreq; qDebug() << "start lamda (m):\t" << startlamda; qDebug() << "rowCount:\t" << rowCount; qDebug() << "colCount:\t" << colCount; qDebug() << "PRFCount:\t" << PRFCount; qDebug() << "PlusePoints:\t" << PlusePoints; qDebug() << "Rnear:\t" << Rnear; qDebug() << "Rfar:\t" << Rfar; qDebug() << "refRange:\t" << refRange; // 方位向分辨率 long echoBlockline = Memory1GB / 8 / 2 / PlusePoints * 2; //2GB echoBlockline = echoBlockline < 1 ? 1 : echoBlockline; long imageBlockline = Memory1GB / 8 / 2 / colCount * 2; //2GB imageBlockline = imageBlockline < 1 ? 1 : imageBlockline; gdalImage imageXYZ(this->outRasterXYZPath); gdalImageComplex imagetimeimg(this->TimeEchoDataPath); long startimgrowid = 0; for (startimgrowid = 0; startimgrowid < rowCount; startimgrowid = startimgrowid + imageBlockline) { long tempimgBlockline = imageBlockline; if (startimgrowid + imageBlockline >= rowCount) { tempimgBlockline = rowCount - startimgrowid; } qDebug() << "\r image create Row Range :\t"< img_x = readDataArr(imageXYZ,startimgrowid,0,tempimgBlockline,colCount,1,GDALREADARRCOPYMETHOD::VARIABLEMETHOD); std::shared_ptr img_y = readDataArr(imageXYZ,startimgrowid,0,tempimgBlockline,colCount,2,GDALREADARRCOPYMETHOD::VARIABLEMETHOD); std::shared_ptr img_z = readDataArr(imageXYZ,startimgrowid,0,tempimgBlockline,colCount,3,GDALREADARRCOPYMETHOD::VARIABLEMETHOD); std::shared_ptr> imgArr = this->L1ds->getImageRaster(startimgrowid, tempimgBlockline); // 获取回波 long startechoid = 0; long iffeechoLen = PlusePoints; for (long startechoid = 0; startechoid < PRFCount; startechoid = startechoid + echoBlockline) { long tempechoBlockline = echoBlockline; if (startechoid + tempechoBlockline >= PRFCount) { tempechoBlockline = PRFCount - startechoid; } std::shared_ptr> echoArr = readDataArrComplex < std::complex>(imagetimeimg,startechoid,long(0), tempechoBlockline, iffeechoLen, 1, GDALREADARRCOPYMETHOD::VARIABLEMETHOD);//; this->L0ds->getEchoArr(startechoid, tempechoBlockline); std::shared_ptr antpx(new double[tempechoBlockline],delArrPtr); std::shared_ptr antpy(new double[tempechoBlockline], delArrPtr); std::shared_ptr antpz(new double[tempechoBlockline], delArrPtr); for (long anti = 0; anti < tempechoBlockline; anti++) { // 天线坐标 antpx.get()[anti] = Pxs.get()[anti + startechoid]; antpy.get()[anti] = Pys.get()[anti + startechoid]; antpz.get()[anti] = Pzs.get()[anti + startechoid]; } TBPImageGPUAlg2( antpx, antpy, antpz, // 坐标 img_x, img_y, img_z, // 图像坐标 echoArr, // 回波 imgArr, // 图像 startfreq, dx, Rnear, Rfar, refRange, tempimgBlockline, colCount, tempechoBlockline, PlusePoints ); qDebug() << QString(" image block PRF:[%1] \t").arg((startechoid + tempechoBlockline) * 100.0 / PRFCount) << startechoid << "\t-\t" << startechoid + tempechoBlockline; } this->L1ds->saveImageRaster(imgArr, startimgrowid, tempimgBlockline); } qDebug() << "\r[" << QDateTime::currentDateTime().toString("yyyy-MM-dd hh:mm:ss.zzz") << "] image writing:\t" << this->L1ds->getxmlFilePath(); this->L1ds->saveToXml(); return ErrorCode::SUCCESS; } void TBPImageGPUAlg2(std::shared_ptr antPx, std::shared_ptr antPy, std::shared_ptr antPz, std::shared_ptr img_x, std::shared_ptr img_y, std::shared_ptr img_z, std::shared_ptr> echoArr, std::shared_ptr> img_arr, double freq, double dx, double Rnear, double Rfar, double refRange, long rowcount, long colcount, long prfcount, long freqcount ) { // 天线 std::shared_ptr h_antPx ((double*)mallocCUDAHost(sizeof(double) * prfcount),FreeCUDAHost); std::shared_ptr h_antPy ((double*)mallocCUDAHost(sizeof(double) * prfcount),FreeCUDAHost); std::shared_ptr h_antPz ((double*)mallocCUDAHost(sizeof(double) * prfcount),FreeCUDAHost); std::shared_ptr d_antPx ((double*)mallocCUDADevice(sizeof(double) * prfcount),FreeCUDADevice); std::shared_ptr d_antPy ((double*)mallocCUDADevice(sizeof(double) * prfcount),FreeCUDADevice); std::shared_ptr d_antPz ((double*)mallocCUDADevice(sizeof(double) * prfcount),FreeCUDADevice); // 网格坐标 std::shared_ptr h_imgx((double*)mallocCUDAHost(sizeof(double) * rowcount * colcount),FreeCUDAHost); std::shared_ptr h_imgy((double*)mallocCUDAHost(sizeof(double) * rowcount * colcount),FreeCUDAHost); std::shared_ptr h_imgz((double*)mallocCUDAHost(sizeof(double) * rowcount * colcount),FreeCUDAHost); std::shared_ptr d_imgx ((double*)mallocCUDADevice(sizeof(double) * rowcount * colcount),FreeCUDADevice); std::shared_ptr d_imgy ((double*)mallocCUDADevice(sizeof(double) * rowcount * colcount),FreeCUDADevice); std::shared_ptr d_imgz ((double*)mallocCUDADevice(sizeof(double) * rowcount * colcount),FreeCUDADevice); // 回波范围 std::shared_ptr h_echoArr((cuComplex*)mallocCUDAHost(sizeof(cuComplex) * prfcount * freqcount), FreeCUDAHost); std::shared_ptr d_echoArr((cuComplex*)mallocCUDADevice(sizeof(cuComplex) * prfcount * freqcount), FreeCUDADevice); // 成像结果 std::shared_ptr h_imgArr((cuComplex*)mallocCUDAHost(sizeof(cuComplex) * rowcount * colcount), FreeCUDAHost); std::shared_ptr d_imgArr((cuComplex*)mallocCUDADevice(sizeof(cuComplex) * rowcount * colcount), FreeCUDADevice); #pragma omp parallel for for (long i = 0; i < prfcount; i++) { for (long j = 0; j < freqcount; j++) { h_echoArr.get()[i * freqcount + j] = make_cuComplex(echoArr.get()[i * freqcount + j].real(), echoArr.get()[i * freqcount + j].imag()); } } // 天线位置 for (long i = 0; i < prfcount; i++) { h_antPx.get()[i] = antPx.get()[i]; h_antPy.get()[i] = antPy.get()[i]; h_antPz.get()[i] = antPz.get()[i]; } #pragma omp parallel for for (long i = 0; i < rowcount; i++) { for (long j = 0; j < colcount; j++) { h_imgx.get()[i * colcount + j] = img_x.get()[i * colcount + j]; h_imgy.get()[i * colcount + j] = img_y.get()[i * colcount + j]; h_imgz.get()[i * colcount + j] = img_z.get()[i * colcount + j]; } } #pragma omp parallel for for (long i = 0; i < rowcount; i++) { for (long j = 0; j < colcount; j++) { h_imgArr.get()[i * colcount + j].x = img_arr.get()[i * colcount + j].real(); h_imgArr.get()[i * colcount + j].y = img_arr.get()[i * colcount + j].imag(); } } HostToDevice(h_imgArr.get(), d_imgArr.get(), sizeof(cuComplex) * rowcount * colcount); HostToDevice(h_echoArr.get(), d_echoArr.get(), sizeof(cuComplex) * prfcount * freqcount); HostToDevice(h_antPx.get(), d_antPx.get(), sizeof(double) * prfcount); HostToDevice(h_antPy.get(), d_antPy.get(), sizeof(double) * prfcount); HostToDevice(h_antPz.get(), d_antPz.get(), sizeof(double) * prfcount); HostToDevice(h_imgx.get(), d_imgx.get(), sizeof(double) * rowcount * colcount); HostToDevice(h_imgy.get(), d_imgy.get(), sizeof(double) * rowcount * colcount); HostToDevice(h_imgz.get(), d_imgz.get(), sizeof(double) * rowcount * colcount); // 直接使用 double startlamda = LIGHTSPEED / freq; TimeBPImage( d_antPx.get(), d_antPy.get(), d_antPz.get(), d_imgx.get(), d_imgy.get(), d_imgz.get(), d_echoArr.get(), prfcount, freqcount, d_imgArr.get(), rowcount, colcount, startlamda, Rnear, dx, refRange,Rfar ); // Device -> Host DeviceToHost(h_imgArr.get(), d_imgArr.get(), sizeof(cuComplex)* rowcount* colcount); #pragma omp parallel for for (long i = 0; i < rowcount; i++) { for (long j = 0; j < colcount; j++) { img_arr.get()[i * colcount + j] = std::complex(h_imgArr.get()[i * colcount + j].x, h_imgArr.get()[i * colcount + j].y); } } } void TBPImageAlgCls::setGPU(bool flag) { this->GPURUN = flag; } bool TBPImageAlgCls::getGPU( ) { return this->GPURUN; } void TBPImageAlgCls::EchoFreqToTime() { // 读取数据 long PRFCount = this->L0ds->getPluseCount(); long inColCount = this->L0ds->getPlusePoints(); long outColCount = inColCount;// nextpow2(inColCount); this->TimeEchoRowCount = PRFCount; this->TimeEchoColCount = outColCount; qDebug() << "IFFT : " << this->TimeEchoDataPath; qDebug() << "PRF Count:\t" << PRFCount; qDebug() << "inColCount:\t" << inColCount; qDebug() << "outColCount:\t" << outColCount; // 创建二进制文件 gdalImageComplex outTimeEchoImg = CreategdalImageComplexNoProj(this->TimeEchoDataPath, this->TimeEchoRowCount, this->TimeEchoColCount, 1); // 分块 long echoBlockline = Memory1GB / 8 / 2 / outColCount * 3; //1GB echoBlockline = echoBlockline < 1 ? 1 : echoBlockline; long startechoid = 0; for (long startechoid = 0; startechoid < PRFCount; startechoid = startechoid + echoBlockline) { long tempechoBlockline = echoBlockline; if (startechoid + tempechoBlockline >= PRFCount) { tempechoBlockline = PRFCount - startechoid; } std::shared_ptr> echoArr = this->L0ds->getEchoArr(startechoid, tempechoBlockline); std::shared_ptr> IFFTArr = outTimeEchoImg.getDataComplexSharePtr(startechoid, 0, tempechoBlockline, outColCount, 1); std::shared_ptr host_echoArr((cuComplex*)mallocCUDAHost(sizeof(cuComplex) * tempechoBlockline * outColCount), FreeCUDAHost); std::shared_ptr host_IFFTechoArr((cuComplex*)mallocCUDAHost(sizeof(cuComplex) * tempechoBlockline * outColCount), FreeCUDAHost); memset(host_echoArr.get(), 0, sizeof(cuComplex) * tempechoBlockline * outColCount); #pragma omp parallel for for (long ii = 0; ii < tempechoBlockline; ii++) { for (long jj = 0; jj < inColCount; jj++) { host_echoArr.get()[ii * outColCount + jj] = make_cuComplex(echoArr.get()[ii * inColCount + jj].real(), echoArr.get()[ii * inColCount + jj].imag()); } } #pragma omp parallel for for (long ii = 0; ii < tempechoBlockline * outColCount; ii++) { host_IFFTechoArr.get()[ii] = make_cuComplex(0, 0); } std::shared_ptr device_echoArr((cuComplex*)mallocCUDADevice(sizeof(cuComplex) * tempechoBlockline * inColCount), FreeCUDADevice); std::shared_ptr device_IFFTechoArr((cuComplex*)mallocCUDADevice(sizeof(cuComplex) * tempechoBlockline * outColCount), FreeCUDADevice); HostToDevice(host_echoArr.get(), device_echoArr.get(), sizeof(cuComplex) * tempechoBlockline * inColCount); HostToDevice(host_IFFTechoArr.get(), device_IFFTechoArr.get(), sizeof(cuComplex) * tempechoBlockline * outColCount); CUDAIFFT(device_echoArr.get(), device_IFFTechoArr.get(), tempechoBlockline, outColCount, outColCount); FFTShift1D(device_IFFTechoArr.get(), tempechoBlockline, outColCount); DeviceToHost(host_IFFTechoArr.get(), device_IFFTechoArr.get(), sizeof(cuComplex) * tempechoBlockline * outColCount); #pragma omp parallel for for (long ii = 0; ii < tempechoBlockline; ii++) { for (long jj = 0; jj < outColCount; jj++) { IFFTArr.get()[ii * outColCount + jj] = std::complex(host_IFFTechoArr.get()[ii * outColCount + jj].x, host_IFFTechoArr.get()[ii * outColCount + jj].y); //IFFTArr.get()[ii * outColCount + jj] = std::complex(host_echoArr.get()[ii * outColCount + jj].x, host_echoArr.get()[ii * outColCount + jj].y); } } outTimeEchoImg.saveImage(IFFTArr, startechoid, 0, tempechoBlockline, outColCount, 1); qDebug() << QString(" image block PRF:[%1] \t").arg((startechoid + tempechoBlockline) * 100.0 / PRFCount) << startechoid << "\t-\t" << startechoid + tempechoBlockline; } return; }