RasterProcessTool/Toolbox/SimulationSARTool/SimulationSAR/TBPImageAlgCls.cpp

620 lines
23 KiB
C++
Raw Blame History

This file contains ambiguous Unicode characters!

This file contains ambiguous Unicode characters that may be confused with others in your current locale. If your use case is intentional and legitimate, you can safely ignore this warning. Use the Escape button to highlight these characters.

#include "stdafx.h"
#include "TBPImageAlgCls.h"
#include <QDateTime>
#include <QDebug>
#include <QString>
#include <cmath>
#include <QProgressDialog>
#include <QMessageBox>
#include "GPUTool.cuh"
#include "GPUTBPImage.cuh"
#include "ImageOperatorBase.h"
#include "BPBasic0_CUDA.cuh"
#include "BaseTool.h"
void CreatePixelXYZ(std::shared_ptr<EchoL0Dataset> 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<double> 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<double> Pxs((double*)mallocCUDAHost(sizeof(double)*prfcount), FreeCUDAHost);
std::shared_ptr<double> Pys((double*)mallocCUDAHost(sizeof(double) * prfcount), FreeCUDAHost);
std::shared_ptr<double> Pzs((double*)mallocCUDAHost(sizeof(double) * prfcount), FreeCUDAHost);
std::shared_ptr<double> AntDirectX((double*)mallocCUDAHost(sizeof(double) * prfcount), FreeCUDAHost);
std::shared_ptr<double> AntDirectY((double*)mallocCUDAHost(sizeof(double) * prfcount), FreeCUDAHost);
std::shared_ptr<double> AntDirectZ((double*)mallocCUDAHost(sizeof(double) * prfcount), FreeCUDAHost);
{
std::shared_ptr<double> 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<double> d_Pxs((double*)mallocCUDADevice(sizeof(double) * prfcount), FreeCUDADevice);
std::shared_ptr<double> d_Pys((double*)mallocCUDADevice(sizeof(double) * prfcount), FreeCUDADevice);
std::shared_ptr<double> d_Pzs((double*)mallocCUDADevice(sizeof(double) * prfcount), FreeCUDADevice);
std::shared_ptr<double> d_AntDirectX((double*)mallocCUDADevice(sizeof(double) * prfcount), FreeCUDADevice);
std::shared_ptr<double> d_AntDirectY((double*)mallocCUDADevice(sizeof(double) * prfcount), FreeCUDADevice);
std::shared_ptr<double> 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<double> demx = readDataArr<double>(xyzRaster, 0, startcolidx, prfcount, tempechocol, 1, GDALREADARRCOPYMETHOD::VARIABLEMETHOD);
std::shared_ptr<double> demy = readDataArr<double>(xyzRaster, 0, startcolidx, prfcount, tempechocol, 2, GDALREADARRCOPYMETHOD::VARIABLEMETHOD);
std::shared_ptr<double> demz = readDataArr<double>(xyzRaster, 0, startcolidx, prfcount, tempechocol, 3, GDALREADARRCOPYMETHOD::VARIABLEMETHOD);
std::shared_ptr<double> h_demx((double*)mallocCUDAHost(sizeof(double) * prfcount*tempechocol), FreeCUDAHost);
std::shared_ptr<double> h_demy((double*)mallocCUDAHost(sizeof(double) * prfcount*tempechocol), FreeCUDAHost);
std::shared_ptr<double> 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<double> d_demx((double*)mallocCUDADevice(sizeof(double) * prfcount * tempechocol), FreeCUDADevice);
std::shared_ptr<double> d_demy((double*)mallocCUDADevice(sizeof(double) * prfcount * tempechocol), FreeCUDADevice);
std::shared_ptr<double> 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<EchoL0Dataset> 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<EchoL0Dataset> inL0ds)
{
this->L0ds = inL0ds;
}
void TBPImageAlgCls::setImageL1(std::shared_ptr<SARSimulationImageL1Dataset> inL1ds)
{
this->L1ds = inL1ds;
}
std::shared_ptr<EchoL0Dataset> TBPImageAlgCls::getEchoL1()
{
return this->L0ds;
}
std::shared_ptr<SARSimulationImageL1Dataset> 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);
return this->ProcessWithGridNet(num_thread, outRasterXYZ);
}
ErrorCode TBPImageAlgCls::ProcessWithGridNet(long num_thread,QString xyzRasterPath)
{
this->outRasterXYZPath = xyzRasterPath;
//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<std::complex<double>> 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<double>(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 freqpoint = this->L0ds->getPlusePoints();
long PRFCount = this->L0ds->getPluseCount();
// 图像范围
long imWidth = this->L1ds->getcolCount();
long imHeight = this->L1ds->getrowCount();
double refRange = this->L0ds->getRefPhaseRange();
// 内存分配大小
long echoBlockline = Memory1GB / 8 / 2 / freqpoint * 4; //2GB
echoBlockline = echoBlockline < 1 ? 1 : echoBlockline;
long imageBlockline = Memory1GB / 8 / 2 / imHeight * 4; //2GB
imageBlockline = imageBlockline < 1 ? 1 : imageBlockline;
qDebug() << "echo block rows: " << echoBlockline;
qDebug() << "image block rows: " << imageBlockline;
// 天线坐标
std::shared_ptr<double> Pxs(new double[this->L0ds->getPluseCount()], delArrPtr);
std::shared_ptr<double> Pys(new double[this->L0ds->getPluseCount()], delArrPtr);
std::shared_ptr<double> Pzs(new double[this->L0ds->getPluseCount()], delArrPtr);
{
std::shared_ptr<double> antpos = this->L0ds->getAntPos();
double time = 0;
double Px = 0;
double Py = 0;
double Pz = 0;
for (long i = 0; i < PRFCount; 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 startFreq=this->L0ds->getCenterFreq() - this->L0ds->getBandwidth() / 2;
gdalImage imageXYZ(this->outRasterXYZPath); // 图像坐标
for (long img_rid = 0; img_rid < imHeight; img_rid = img_rid + imageBlockline) {
// 获取坐标范围
long imrowcount = imageBlockline;
long imcolcount = imWidth;
std::shared_ptr<double> img_x = readDataArr<double>(imageXYZ, img_rid, 0, imrowcount, imcolcount, 1, GDALREADARRCOPYMETHOD::VARIABLEMETHOD);
std::shared_ptr<double> img_y = readDataArr<double>(imageXYZ, img_rid, 0, imrowcount, imcolcount, 2, GDALREADARRCOPYMETHOD::VARIABLEMETHOD);
std::shared_ptr<double> img_z = readDataArr<double>(imageXYZ, img_rid, 0, imrowcount, imcolcount, 3, GDALREADARRCOPYMETHOD::VARIABLEMETHOD);
std::shared_ptr<std::complex<double>> imgArr = this->L1ds->getImageRaster(img_rid, imrowcount); // 回波值
double img_x_min=0, img_x_max=0;
double img_y_min=0, img_y_max=0;
double img_z_min=0, img_z_max=0;
minValueInArr<double>(img_x.get(), imrowcount * imcolcount, img_x_min);
maxValueInArr<double>(img_x.get(), imrowcount * imcolcount, img_x_max);
minValueInArr<double>(img_y.get(), imrowcount * imcolcount, img_y_min);
maxValueInArr<double>(img_y.get(), imrowcount * imcolcount, img_y_max);
minValueInArr<double>(img_z.get(), imrowcount * imcolcount, img_z_min);
maxValueInArr<double>(img_z.get(), imrowcount * imcolcount, img_z_max);
qDebug() << "imgX:\t" << img_x_min << " ~ " << img_x_max;
qDebug() << "imgY:\t" << img_y_min << " ~ " << img_y_max;
qDebug() << "imgZ:\t" << img_z_min << " ~ " << img_z_max;
//shared_complexPtrToHostCuComplex(std::complex<double>*src, cuComplex * dst, long len)
// 处理
GPUDATA h_data;
h_data.Nfft = freqpoint;
h_data.K = freqpoint;
h_data.deltaF = this->L0ds->getBandwidth() / (freqpoint - 1);
// 计算maxWr需要先计算deltaF
double deltaF = h_data.deltaF; // 从输入参数获取
double maxWr = 299792458.0f / (2.0f * deltaF);
// 生成r_vec主机端
double* r_vec_host = (double*)mallocCUDAHost(sizeof(double) * h_data.Nfft);// new double[data.Nfft];
const double step = maxWr / h_data.Nfft;
const double start = -1*h_data.Nfft / 2.0f * step;
printf("nfft=%d\n", h_data.Nfft);
for (int i = 0; i < h_data.Nfft; ++i) {
r_vec_host[i] = start + i * step;
}
h_data.r_vec = r_vec_host;
h_data.x_mat = img_x.get();//地面
h_data.y_mat = img_y.get();
h_data.z_mat = img_z.get();
h_data.nx = imcolcount;
h_data.ny = imrowcount;
minValueInArr<double>(h_data.x_mat, h_data.nx * h_data.ny, img_x_min);
maxValueInArr<double>(h_data.x_mat, h_data.nx * h_data.ny, img_x_max);
minValueInArr<double>(h_data.y_mat, h_data.nx * h_data.ny, img_y_min);
maxValueInArr<double>(h_data.y_mat, h_data.nx * h_data.ny, img_y_max);
minValueInArr<double>(h_data.z_mat, h_data.nx * h_data.ny, img_z_min);
maxValueInArr<double>(h_data.z_mat, h_data.nx * h_data.ny, img_z_max);
qDebug() << "imgX:\t" << img_x_min << " ~ " << img_x_max;
qDebug() << "imgY:\t" << img_y_min << " ~ " << img_y_max;
qDebug() << "imgZ:\t" << img_z_min << " ~ " << img_z_max;
qDebug() << "imgXYZ" << h_data.x_mat[5016600] << " , " << h_data.y_mat[5016600] << " , " << h_data.z_mat[5016600] << " , ";
h_data.R0 = refRange;// 参考斜距
h_data.im_final = (cuComplex*)mallocCUDAHost(sizeof(cuComplex) * imrowcount * imcolcount);
shared_complexPtrToHostCuComplex(imgArr.get(), h_data.im_final, imrowcount * imcolcount);
// 保存fft 结果
this->TimeEchoDataPath = JoinPath(this->L1ds->getoutFolderPath(), this->L0ds->getSimulationTaskName() + "_Timeecho.bin");
gdalImageComplex outTimeEchoImg = CreategdalImageComplexNoProj(this->TimeEchoDataPath,PRFCount, freqpoint, 1);
for (long prfid = 0; prfid < PRFCount; prfid = prfid + echoBlockline) {
long ehcoprfcount = echoBlockline;
long echofreqpoint = freqpoint;
std::shared_ptr<std::complex<double>> echoArr = this->L0ds->getEchoArr(prfid, ehcoprfcount);
// 复制天线方向图
std::shared_ptr<double> antpx(new double[ehcoprfcount], delArrPtr);
std::shared_ptr<double> antpy(new double[ehcoprfcount], delArrPtr);
std::shared_ptr<double> antpz(new double[ehcoprfcount], delArrPtr);
for (long anti = 0; anti < ehcoprfcount; anti++) {
if (anti + prfid < PRFCount) {
antpx.get()[anti] = Pxs.get()[anti + prfid];
antpy.get()[anti] = Pys.get()[anti + prfid];
antpz.get()[anti] = Pzs.get()[anti + prfid];
//if (abs(antpx.get()[anti]) < 10 || abs(antpy.get()[anti]) < 10 || abs(antpz.get()[anti]) < 10) {
// qDebug() << anti << ":" << antpx.get()[anti] << "," << antpy.get()[anti] << "," << antpz.get()[anti];
//}
}
}
// 起始频率
h_data.minF = (double*)mallocCUDAHost(sizeof(double) * ehcoprfcount);
h_data.Freq = (double*)mallocCUDAHost(sizeof(double) * ehcoprfcount);
for (long anti = 0; anti < ehcoprfcount; anti++) {
h_data.minF[anti] = startFreq;
h_data.Freq[anti] = startFreq;
}
h_data.AntX = antpx.get(); // 天线
h_data.AntY = antpy.get();
h_data.AntZ = antpz.get();
// checkout
if (!this->checkZeros(h_data.AntX, ehcoprfcount) ||
!this->checkZeros(h_data.AntY, ehcoprfcount) ||
!this->checkZeros(h_data.AntZ, ehcoprfcount)
) {
printf("the ant position is not zeros!!!");
}
h_data.Np = ehcoprfcount;
h_data.phdata= (cuComplex*)mallocCUDAHost(sizeof(cuComplex) * ehcoprfcount * echofreqpoint);
qDebug() << ehcoprfcount <<":\t"<<h_data.Np;
shared_complexPtrToHostCuComplex(echoArr.get(), h_data.phdata, ehcoprfcount * echofreqpoint);
{
// BP
GPUDATA d_data;
initGPUData(h_data, d_data);
qDebug() << "------------------------------------------------------";
double min_r_vec = 0, max_r_cev = 0;
minValueInArr<double>(h_data.r_vec, h_data.Nfft, min_r_vec);
maxValueInArr<double>(h_data.r_vec, h_data.Nfft, max_r_cev);
qDebug() << "r_vec:\t" << min_r_vec << " ~ " << max_r_cev;
minValueInArr<double>(h_data.Freq, h_data.Nfft, min_r_vec);
maxValueInArr<double>(h_data.Freq, h_data.Nfft, max_r_cev);
qDebug() << "Freq:\t" << min_r_vec << " ~ " << max_r_cev;
// 打印参数
qDebug() << "R0\t" << h_data.R0;
qDebug() << "deltaF\t" << h_data.deltaF;
qDebug() << "Nfft\t" << h_data.Nfft;
qDebug() << "R0\t" << h_data.R0;
qDebug() << "K\t" << h_data.K;
qDebug() << "Np\t" << h_data.Np;
qDebug() << "nx\t" << h_data.nx;
qDebug() << "ny\t" << h_data.ny;
qDebug() << "------------------------------------------------------";
bpBasic0CUDA(d_data, 0);
printf("BP finished!!!\n");
DeviceToHost(h_data.im_final, d_data.im_final, sizeof(cuComplex) * h_data.nx * h_data.ny);
DeviceToHost(h_data.phdata, d_data.phdata, sizeof(cuComplex) * h_data.Np * h_data.Nfft);
gdalImageComplex outTimeEchoImg = CreategdalImageComplexNoProj(this->TimeEchoDataPath, h_data.Np, h_data.Nfft, 1);
std::shared_ptr<std::complex<double>> data_time(new std::complex<double>[h_data.Np * h_data.Nfft], delArrPtr);
for (long i = 0; i < h_data.Np * h_data.Nfft; i++) {
data_time.get()[i] = std::complex<double>(h_data.phdata[i].x, h_data.phdata[i].y);
}
outTimeEchoImg.saveImage(data_time, 0, 0, h_data.Np, h_data.Nfft, 1);
HostCuComplexToshared_complexPtr(h_data.im_final, imgArr.get(), imrowcount * imcolcount);
freeGPUData(d_data);
}
this->L1ds->saveImageRaster(imgArr, img_rid, imrowcount);
}
freeHostData(h_data);
}
this->L1ds->saveToXml();
return ErrorCode::SUCCESS;
}
ErrorCode TBPImageAlgCls::BPProcessBlockGPU()
{
return ErrorCode::SUCCESS;
}
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<std::complex<double>> echoArr = this->L0ds->getEchoArr(startechoid, tempechoBlockline);
std::shared_ptr<std::complex<double>> IFFTArr = outTimeEchoImg.getDataComplexSharePtr(startechoid, 0, tempechoBlockline, outColCount, 1);
std::shared_ptr<cuComplex> host_echoArr((cuComplex*)mallocCUDAHost(sizeof(cuComplex) * tempechoBlockline * outColCount), FreeCUDAHost);
std::shared_ptr<cuComplex> 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<cuComplex> device_echoArr((cuComplex*)mallocCUDADevice(sizeof(cuComplex) * tempechoBlockline * inColCount), FreeCUDADevice);
std::shared_ptr<cuComplex> 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<double>(host_IFFTechoArr.get()[ii * outColCount + jj].x, host_IFFTechoArr.get()[ii * outColCount + jj].y);
//IFFTArr.get()[ii * outColCount + jj] = std::complex<double>(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;
}
bool TBPImageAlgCls::checkZeros(double* data, long long len)
{
bool flag = true;
#pragma omp parallel for
for (long long i = 0; i < len; i++) {
if (abs(data[i]) < 1e-7) {
flag= false;
break;
}
}
return flag;
}