#include "stdafx.h" #include "RTPCProcessCls.h" #include "BaseConstVariable.h" #include "SARSatelliteSimulationAbstractCls.h" #include "SARSimulationTaskSetting.h" #include "SatelliteOribtModel.h" #include #include "ImageOperatorBase.h" #include "GeoOperator.h" #include "EchoDataFormat.h" #include #include #include #include #include #ifdef DEBUGSHOWDIALOG #include "ImageShowDialogClass.h" #endif #ifdef __CUDANVCC___ #include "GPUTool.cuh" #endif // __CUDANVCC___ RTPCProcessCls::RTPCProcessCls() { this->PluseCount = 0; this->PlusePoint = 0; this->TaskSetting = nullptr; this->EchoSimulationData = nullptr; this->DEMTiffPath = ""; this->LandCoverPath = ""; this->HHSigmaPath = ""; this->HVSigmaPath = ""; this->VHSigmaPath = ""; this->VVSigmaPath = ""; this->OutEchoPath = ""; this->DEMTiffPath.clear(); this->LandCoverPath.clear(); this->HHSigmaPath.clear(); this->HVSigmaPath.clear(); this->VHSigmaPath.clear(); this->VVSigmaPath.clear(); this->OutEchoPath.clear(); this->SigmaDatabasePtr = std::shared_ptr(new SigmaDatabase); } RTPCProcessCls::~RTPCProcessCls() { } void RTPCProcessCls::setTaskSetting(std::shared_ptr < AbstractSARSatelliteModel> TaskSetting) { this->TaskSetting = std::shared_ptr < AbstractSARSatelliteModel>(TaskSetting); qDebug() << "RTPCProcessCls::setTaskSetting"; } void RTPCProcessCls::setEchoSimulationDataSetting(std::shared_ptr EchoSimulationData) { this->EchoSimulationData = std::shared_ptr(EchoSimulationData); qDebug() << "RTPCProcessCls::setEchoSimulationDataSetting"; } void RTPCProcessCls::setTaskFileName(QString EchoFileName) { this->TaskFileName = EchoFileName; } void RTPCProcessCls::setDEMTiffPath(QString DEMTiffPath) { this->DEMTiffPath = DEMTiffPath; } void RTPCProcessCls::setLandCoverPath(QString LandCoverPath) { this->LandCoverPath = LandCoverPath; } void RTPCProcessCls::setHHSigmaPath(QString HHSigmaPath) { this->HHSigmaPath = HHSigmaPath; } void RTPCProcessCls::setHVSigmaPath(QString HVSigmaPath) { this->HVSigmaPath = HVSigmaPath; } void RTPCProcessCls::setVHSigmaPath(QString VHSigmaPath) { this->VHSigmaPath = VHSigmaPath; } void RTPCProcessCls::setVVSigmaPath(QString VVSigmaPath) { this->VVSigmaPath = VVSigmaPath; } void RTPCProcessCls::setOutEchoPath(QString OutEchoPath) { this->OutEchoPath = OutEchoPath; } ErrorCode RTPCProcessCls::Process(long num_thread) { // RTPC 算法 qDebug() << u8"params init ...."; ErrorCode stateCode = this->InitParams(); if (stateCode != ErrorCode::SUCCESS) { return stateCode; } else {} qDebug() << "DEMMainProcess"; stateCode = this->DEMPreprocess(); if (stateCode != ErrorCode::SUCCESS) { return stateCode; } else {} qDebug() << "RTPCMainProcess"; stateCode = this->RTPCMainProcess(num_thread); if (stateCode != ErrorCode::SUCCESS) { return stateCode; } else {} return ErrorCode::SUCCESS; } ErrorCode RTPCProcessCls::InitParams() { if (nullptr == this->TaskSetting || this->DEMTiffPath.isEmpty() || this->LandCoverPath.isEmpty() || this->HHSigmaPath.isEmpty() || this->HVSigmaPath.isEmpty() || this->VHSigmaPath.isEmpty() || this->VVSigmaPath.isEmpty()) { return ErrorCode::RTPC_PARAMSISEMPTY; } else { } // 归一化绝对路径 this->OutEchoPath = QDir(this->OutEchoPath).absolutePath(); // 回波大小 double imgStart_end = this->TaskSetting->getSARImageEndTime() - this->TaskSetting->getSARImageStartTime(); this->PluseCount = ceil(imgStart_end * this->TaskSetting->getPRF()); double rangeTimeSample = (this->TaskSetting->getFarRange() - this->TaskSetting->getNearRange()) * 2.0 / LIGHTSPEED; this->PlusePoint = ceil(rangeTimeSample * this->TaskSetting->getFs()); // 初始化回波存放位置 qDebug() << "--------------Echo Data Setting ---------------------------------------"; this->EchoSimulationData = std::shared_ptr(new EchoL0Dataset); this->EchoSimulationData->setCenterFreq(this->TaskSetting->getCenterFreq()); this->EchoSimulationData->setNearRange(this->TaskSetting->getNearRange()); this->EchoSimulationData->setFarRange(this->TaskSetting->getFarRange()); this->EchoSimulationData->setFs(this->TaskSetting->getFs()); this->EchoSimulationData->setCenterAngle(this->TaskSetting->getCenterLookAngle()); this->EchoSimulationData->setLookSide(this->TaskSetting->getIsRightLook() ? "R" : "L"); this->EchoSimulationData->OpenOrNew(OutEchoPath, TaskFileName, PluseCount, PlusePoint); QString tmpfolderPath = QDir(OutEchoPath).filePath("tmp"); if (QDir(tmpfolderPath).exists() == false) { QDir(OutEchoPath).mkpath(tmpfolderPath); } this->tmpfolderPath = tmpfolderPath; return ErrorCode::SUCCESS; } ErrorCode RTPCProcessCls::DEMPreprocess() { this->demxyzPath = QDir(tmpfolderPath).filePath("demxyz.tif"); gdalImage demds(this->DEMTiffPath); gdalImage demxyz = CreategdalImage(demxyzPath, demds.height, demds.width, 3, demds.gt, demds.projection, true, true);// X,Y,Z // 分块计算并转换为XYZ Eigen::MatrixXd demArr = demds.getData(0, 0, demds.height, demds.width, 1); Eigen::MatrixXd demR = demArr; Landpoint LandP{ 0,0,0 }; Point3 GERpoint{ 0,0,0 }; double R = 0; double dem_row = 0, dem_col = 0, dem_alt = 0; long line_invert = 1000; double rowidx = 0; double colidx = 0; for (int max_rows_ids = 0; max_rows_ids < demds.height; max_rows_ids = max_rows_ids + line_invert) { Eigen::MatrixXd demdata = demds.getData(max_rows_ids, 0, line_invert, demds.width, 1); Eigen::MatrixXd xyzdata_x = demdata.array() * 0; Eigen::MatrixXd xyzdata_y = demdata.array() * 0; Eigen::MatrixXd xyzdata_z = demdata.array() * 0; int datarows = demdata.rows(); int datacols = demdata.cols(); for (int i = 0; i < datarows; i++) { for (int j = 0; j < datacols; j++) { rowidx = i + max_rows_ids; colidx = j; demds.getLandPoint(rowidx, colidx, demdata(i, j), LandP); // 获取地理坐标 LLA2XYZ(LandP, GERpoint); // 经纬度转换为地心坐标系 xyzdata_x(i, j) = GERpoint.x; xyzdata_y(i, j) = GERpoint.y; xyzdata_z(i, j) = GERpoint.z; } } demxyz.saveImage(xyzdata_x, max_rows_ids, 0, 1); demxyz.saveImage(xyzdata_y, max_rows_ids, 0, 2); demxyz.saveImage(xyzdata_z, max_rows_ids, 0, 3); } // 计算坡向角 this->demsloperPath = QDir(tmpfolderPath).filePath("demsloper.tif"); this->demmaskPath = QDir(tmpfolderPath).filePath("demmask.tif"); gdalImage demsloperxyz = CreategdalImage(this->demsloperPath, demds.height, demds.width, 4, demds.gt, demds.projection, true, true);// X,Y,Z,cosangle gdalImage demmask = CreategdalImage(this->demmaskPath, demxyz.height, demxyz.width, 1, demxyz.gt, demxyz.projection, true, true);// X,Y,Z line_invert = 1000; long start_ids = 0; long dem_rows = 0, dem_cols = 0; for (start_ids = 1; start_ids < demds.height; start_ids = start_ids + line_invert) { Eigen::MatrixXd demdata = demds.getData(start_ids - 1, 0, line_invert + 2, demxyz.width, 1); long startlineid = start_ids; Eigen::MatrixXd maskdata = demmask.getData(start_ids - 1, 0, line_invert + 2, demxyz.width, 1); Eigen::MatrixXd demsloper_x = demsloperxyz.getData(start_ids - 1, 0, line_invert + 2, demxyz.width, 1); Eigen::MatrixXd demsloper_y = demsloperxyz.getData(start_ids - 1, 0, line_invert + 2, demxyz.width, 2); Eigen::MatrixXd demsloper_z = demsloperxyz.getData(start_ids - 1, 0, line_invert + 2, demxyz.width, 3); Eigen::MatrixXd demsloper_angle = demsloperxyz.getData(start_ids - 1, 0, line_invert + 2, demxyz.width, 4); maskdata = maskdata.array() * 0; Landpoint p0, p1, p2, p3, p4, pslopeVector, pp; Vector3D slopeVector; dem_rows = maskdata.rows(); dem_cols = maskdata.cols(); double sloperAngle = 0; Vector3D Zaxis = { 0,0,1 }; double rowidx = 0, colidx = 0; for (long i = 1; i < dem_rows - 1; i++) { for (long j = 1; j < dem_cols - 1; j++) { rowidx = i + startlineid; colidx = j; demds.getLandPoint(rowidx, colidx, demdata(i, j), p0); demds.getLandPoint(rowidx - 1, colidx, demdata(i - 1, j), p1); demds.getLandPoint(rowidx, colidx - 1, demdata(i, j - 1), p2); demds.getLandPoint(rowidx + 1, colidx, demdata(i + 1, j), p3); demds.getLandPoint(rowidx, colidx + 1, demdata(i, j + 1), p4); pslopeVector = getSlopeVector(p0, p1, p2, p3, p4); // 地面坡向矢量 slopeVector = { pslopeVector.lon,pslopeVector.lat,pslopeVector.ati }; pp = LLA2XYZ(p0); Zaxis.x = pp.lon; Zaxis.y = pp.lat; Zaxis.z = pp.ati; sloperAngle = getCosAngle(slopeVector, Zaxis); // 地面坡向角 demsloper_x(i, j) = slopeVector.x; demsloper_y(i, j) = slopeVector.y; demsloper_z(i, j) = slopeVector.z; demsloper_angle(i, j) = sloperAngle; maskdata(i, j)++; } } demmask.saveImage(maskdata, start_ids - 1, 0, 1); demsloperxyz.saveImage(demsloper_x, start_ids - 1, 0, 1); demsloperxyz.saveImage(demsloper_y, start_ids - 1, 0, 2); demsloperxyz.saveImage(demsloper_z, start_ids - 1, 0, 3); demsloperxyz.saveImage(demsloper_angle, start_ids - 1, 0, 4); } return ErrorCode::SUCCESS; } ErrorCode RTPCProcessCls::RTPCMainProcess(long num_thread) { omp_set_num_threads(num_thread);// 设置openmp 线程数量 double widthSpace = LIGHTSPEED / 2 / this->TaskSetting->getFs(); double prf_time = 0; double dt = 1 / this->TaskSetting->getPRF();// 获取每次脉冲的时间间隔 bool antflag = true; // 计算天线方向图 Landpoint LandP{ 0,0,0 }; Point3 GERpoint{ 0,0,0 }; double R = 0; double dem_row = 0, dem_col = 0, dem_alt = 0; long double imageStarttime = 0; imageStarttime = this->TaskSetting->getSARImageStartTime(); //std::vector sateOirbtNodes(this->PluseCount); std::shared_ptr sateOirbtNodes(new SatelliteOribtNode[this->PluseCount], delArrPtr); { // 姿态计算不同 // 计算姿态 std::shared_ptr antpos = this->EchoSimulationData->getAntPos(); double dAt = 1e-6; double prf_time_dt = 0; Landpoint InP{ 0,0,0 }, outP{ 0,0,0 }; for (long prf_id = 0; prf_id < this->PluseCount; prf_id++) { prf_time = dt * prf_id; prf_time_dt = prf_time + dAt; SatelliteOribtNode sateOirbtNode; SatelliteOribtNode sateOirbtNode_dAt; this->TaskSetting->getSatelliteOribtNode(prf_time, sateOirbtNode, antflag); this->TaskSetting->getSatelliteOribtNode(prf_time_dt, sateOirbtNode_dAt, antflag); sateOirbtNode.AVx = (sateOirbtNode_dAt.Vx - sateOirbtNode.Vx) / dAt; // 加速度 sateOirbtNode.AVy = (sateOirbtNode_dAt.Vy - sateOirbtNode.Vy) / dAt; sateOirbtNode.AVz = (sateOirbtNode_dAt.Vz - sateOirbtNode.Vz) / dAt; InP.lon = sateOirbtNode.Px; InP.lat = sateOirbtNode.Py; InP.ati = sateOirbtNode.Pz; outP = XYZ2LLA(InP); antpos.get()[prf_id * 19 + 0] = prf_time + imageStarttime; antpos.get()[prf_id * 19 + 1] = sateOirbtNode.Px; antpos.get()[prf_id * 19 + 2] = sateOirbtNode.Py; antpos.get()[prf_id * 19 + 3] = sateOirbtNode.Pz; antpos.get()[prf_id * 19 + 4] = sateOirbtNode.Vx; antpos.get()[prf_id * 19 + 5] = sateOirbtNode.Vy; antpos.get()[prf_id * 19 + 6] = sateOirbtNode.Vz; antpos.get()[prf_id * 19 + 7] = sateOirbtNode.AntDirecX; antpos.get()[prf_id * 19 + 8] = sateOirbtNode.AntDirecY; antpos.get()[prf_id * 19 + 9] = sateOirbtNode.AntDirecZ; antpos.get()[prf_id * 19 + 10] = sateOirbtNode.AVx; antpos.get()[prf_id * 19 + 11] = sateOirbtNode.AVy; antpos.get()[prf_id * 19 + 12] = sateOirbtNode.AVz; antpos.get()[prf_id * 19 + 13] = sateOirbtNode.zeroDopplerDirectX; antpos.get()[prf_id * 19 + 14] = sateOirbtNode.zeroDopplerDirectY; antpos.get()[prf_id * 19 + 15] = sateOirbtNode.zeroDopplerDirectZ; antpos.get()[prf_id * 19 + 16] = outP.lon; antpos.get()[prf_id * 19 + 17] = outP.lat; antpos.get()[prf_id * 19 + 18] = outP.ati; sateOirbtNodes[prf_id] = sateOirbtNode; } this->EchoSimulationData->saveAntPos(antpos); antpos.reset(); qDebug() << "Ant position finished sucessfully !!!"; } // 回波 long echoIdx = 0; double NearRange = this->EchoSimulationData->getNearRange(); // 近斜据 double FarRange = this->EchoSimulationData->getFarRange(); double TimgNearRange = 2 * NearRange / LIGHTSPEED; double TimgFarRange = 2 * FarRange / LIGHTSPEED; double Fs = this->TaskSetting->getFs(); // 距离向采样率 double Pt = this->TaskSetting->getPt() * this->TaskSetting->getGri();// 发射电压 1v //double GainAntLen = -3;// -3dB 为天线半径 long pluseCount = this->PluseCount; double lamda = this->TaskSetting->getCenterLamda(); // 波长 // 天线方向图 std::shared_ptr TransformPattern = this->TaskSetting->getTransformRadiationPattern(); // 发射天线方向图 std::shared_ptr ReceivePattern = this->TaskSetting->getReceiveRadiationPattern(); // 接收天线方向图 std::shared_ptr> echo = this->EchoSimulationData->getEchoArr(); long PlusePoint = this->EchoSimulationData->getPlusePoints(); // 初始化 为 0 for (long i = 0; i < pluseCount * PlusePoint; i++) { echo.get()[i] = std::complex(0, 0); } this->EchoSimulationData->saveEchoArr(echo, 0, PluseCount); POLARTYPEENUM polartype = this->TaskSetting->getPolarType(); #ifndef __CUDANVCC___ QMessageBox::information(this, u8"程序提示", u8"请确定安装了CUDA库"); #else // RTPC CUDA版本 if (pluseCount * 4 * 18 > Memory1MB * 100) { long max = Memory1MB * 100 / 4 / 20 / PluseCount; QMessageBox::warning(nullptr, u8"仿真场景太大了", u8"当前频点数下,脉冲数量最多为:"+QString::number(max)); } gdalImage demxyz(this->demxyzPath);// 地面点坐标 gdalImage demlandcls(this->LandCoverPath);// 地表覆盖类型 gdalImage demsloperxyz(this->demsloperPath);// 地面坡向 long demRow = demxyz.height; long demCol = demxyz.width; long blokline = 100; // 每块 250MB*16 = 4GB // dem 9 // blokline = Memory1MB * 500 / 8 / demCol; blokline = blokline < 1 ? 1 : blokline; bool bloklineflag = false; Eigen::MatrixXd dem_x = demxyz.getData(0, 0, blokline, demxyz.width, 1); // 地面坐标 long tempDemRows = dem_x.rows(); long tempDemCols = dem_x.cols(); Eigen::MatrixXd dem_y = Eigen::MatrixXd::Zero(tempDemRows, tempDemCols); Eigen::MatrixXd dem_z = Eigen::MatrixXd::Zero(tempDemRows, tempDemCols); Eigen::MatrixXd demsloper_x = Eigen::MatrixXd::Zero(tempDemRows, tempDemCols); Eigen::MatrixXd demsloper_y = Eigen::MatrixXd::Zero(tempDemRows, tempDemCols); Eigen::MatrixXd demsloper_z = Eigen::MatrixXd::Zero(tempDemRows, tempDemCols); Eigen::MatrixXd sloperAngle = Eigen::MatrixXd::Zero(tempDemRows, tempDemCols); float* h_dem_x; float* h_dem_y; float* h_dem_z; float* h_demsloper_x; float* h_demsloper_y; float* h_demsloper_z; float* h_demsloper_angle; float* d_dem_x; float* d_dem_y; float* d_dem_z; float* d_demsloper_x; float* d_demsloper_y; float* d_demsloper_z; float* d_demsloper_angle; mallocCUDAHost((void*)h_dem_x, sizeof(float) * blokline * tempDemCols); mallocCUDAHost((void*)h_dem_y, sizeof(float) * blokline * tempDemCols); mallocCUDAHost((void*)h_dem_z, sizeof(float) * blokline * tempDemCols); mallocCUDAHost((void*)h_demsloper_x, sizeof(float) * blokline * tempDemCols); mallocCUDAHost((void*)h_demsloper_y, sizeof(float) * blokline * tempDemCols); mallocCUDAHost((void*)h_demsloper_z, sizeof(float) * blokline * tempDemCols); mallocCUDAHost((void*)h_demsloper_angle, sizeof(float) * blokline * tempDemCols); mallocCUDADevice((void*)d_dem_x, sizeof(float) * blokline * tempDemCols); // 7 mallocCUDADevice((void*)d_dem_y, sizeof(float) * blokline * tempDemCols); mallocCUDADevice((void*)d_dem_z, sizeof(float) * blokline * tempDemCols); mallocCUDADevice((void*)d_demsloper_x, sizeof(float) * blokline * tempDemCols); mallocCUDADevice((void*)d_demsloper_y, sizeof(float) * blokline * tempDemCols); mallocCUDADevice((void*)d_demsloper_z, sizeof(float) * blokline * tempDemCols); mallocCUDADevice((void*)d_demsloper_angle, sizeof(float) * blokline * tempDemCols); float* h_dem_theta; // 天线方向图 float* h_dem_phi; float* d_dem_theta; float* d_dem_phi; mallocCUDAHost((void*)h_dem_theta, sizeof(float) * blokline * tempDemCols); mallocCUDAHost((void*)h_dem_phi, sizeof(float) * blokline * tempDemCols); mallocCUDADevice((void*)d_dem_theta, sizeof(float) * blokline * tempDemCols);// 9 mallocCUDADevice((void*)d_dem_phi, sizeof(float) * blokline * tempDemCols); HostToDevice((void*)h_dem_theta, (void*)d_dem_theta, sizeof(float) * blokline * tempDemCols); HostToDevice((void*)h_dem_phi, (void*)d_dem_phi, sizeof(float) * blokline * tempDemCols); // 提前声明参数变量 float* h_R;// 辐射方向 float* h_localangle;//入射角 float* d_R;// 辐射方向 float* d_localangle;//入射角 mallocCUDAHost((void*)h_R, sizeof(float) * blokline * tempDemCols); mallocCUDAHost((void*)h_localangle, sizeof(float) * blokline * tempDemCols); // 11 mallocCUDADevice((void*)d_R, sizeof(float) * blokline * tempDemCols); mallocCUDADevice((void*)d_localangle, sizeof(float) * blokline * tempDemCols); float* h_RstX; float* h_RstY; float* h_RstZ; float* d_RstX; float* d_RstY; float* d_RstZ; mallocCUDAHost((void*)h_RstX, sizeof(float) * blokline * tempDemCols); mallocCUDAHost((void*)h_RstY, sizeof(float) * blokline * tempDemCols); mallocCUDAHost((void*)h_RstZ, sizeof(float) * blokline * tempDemCols); // 14 mallocCUDADevice((void*)d_RstX, sizeof(float)* blokline* tempDemCols); mallocCUDADevice((void*)d_RstY, sizeof(float)* blokline* tempDemCols); mallocCUDADevice((void*)d_RstZ, sizeof(float)* blokline* tempDemCols); float* h_sigma0; float* h_TransAnt; float* h_ReciveAnt; float* d_sigma0; float* d_TransAnt; float* d_ReciveAnt; mallocCUDAHost((void*)h_sigma0 , sizeof(float)* blokline* tempDemCols); mallocCUDAHost((void*)h_TransAnt , sizeof(float)* blokline* tempDemCols); mallocCUDAHost((void*)h_ReciveAnt, sizeof(float)* blokline* tempDemCols); // 17 mallocCUDADevice((void*)d_sigma0, sizeof(float)* blokline* tempDemCols); mallocCUDADevice((void*)d_TransAnt, sizeof(float)* blokline* tempDemCols); mallocCUDADevice((void*)d_ReciveAnt, sizeof(float)* blokline* tempDemCols); // 回波 cuComplex* h_echoAmp; cuComplex* d_echoAmp; mallocCUDAHost((void*)h_echoAmp, sizeof(cuComplex)* blokline* tempDemCols); mallocCUDADevice((void*)d_echoAmp, sizeof(cuComplex)* blokline* tempDemCols); //19 long* h_FreqID; long* d_FreqID; mallocCUDAHost((void*)h_FreqID, sizeof(long)* blokline* tempDemCols); mallocCUDADevice((void*)d_FreqID, sizeof(long)* blokline* tempDemCols); //21 Eigen::MatrixXd landcover= Eigen::MatrixXd::Zero(blokline, tempDemCols);// 地面覆盖类型 long startline = 0; for (startline = 0; startline < demRow; startline = startline + blokline) { long newblokline = blokline; if ((startline + blokline) < demRow) { newblokline = demRow - startline; bloklineflag = true; } dem_x = demxyz.getData(startline, 0, newblokline, demxyz.width, 1); // 地面坐标 dem_y = demxyz.getData(startline, 0, newblokline, demxyz.width, 2); dem_z = demxyz.getData(startline, 0, newblokline, demxyz.width, 3); demsloper_x = demsloperxyz.getData(startline, 0, newblokline, demsloperxyz.width, 1); demsloper_y = demsloperxyz.getData(startline, 0, newblokline, demsloperxyz.width, 2); demsloper_z = demsloperxyz.getData(startline, 0, newblokline, demsloperxyz.width, 3); sloperAngle = demsloperxyz.getData(startline, 0, newblokline, demsloperxyz.width, 4); landcover = demlandcls.getData(startline, 0, newblokline, demlandcls.width, 1); if (bloklineflag) { FreeCUDAHost(h_dem_x); FreeCUDADevice(d_dem_x); FreeCUDAHost(h_dem_y); FreeCUDADevice(d_dem_y); FreeCUDAHost(h_dem_z); FreeCUDADevice(d_dem_z); FreeCUDAHost(h_demsloper_x); FreeCUDADevice(d_demsloper_x); FreeCUDAHost(h_demsloper_y); FreeCUDADevice(d_demsloper_y); FreeCUDAHost(h_demsloper_z); FreeCUDADevice(d_demsloper_z); //6 FreeCUDAHost(h_demsloper_angle); FreeCUDADevice(d_demsloper_angle);//7 FreeCUDAHost(h_dem_theta); FreeCUDADevice(d_dem_theta); FreeCUDAHost(h_dem_phi); FreeCUDADevice(d_dem_phi); //9 FreeCUDAHost(h_R); FreeCUDADevice(d_R); FreeCUDAHost(h_localangle); FreeCUDADevice(d_localangle); //11 FreeCUDAHost(h_RstX); FreeCUDADevice(d_RstX); FreeCUDAHost(h_RstY); FreeCUDADevice(d_RstY); FreeCUDAHost(h_RstZ); FreeCUDADevice(d_RstZ); //14 FreeCUDAHost(h_sigma0 ); FreeCUDADevice(d_sigma0); FreeCUDAHost(h_TransAnt ); FreeCUDADevice(d_TransAnt); FreeCUDAHost(h_ReciveAnt); FreeCUDADevice(d_ReciveAnt); //17 FreeCUDAHost(h_echoAmp); FreeCUDADevice(d_echoAmp);//19 FreeCUDAHost(h_FreqID); FreeCUDADevice(d_FreqID);//20 mallocCUDAHost((void*)h_dem_x, sizeof(float) * newblokline * tempDemCols); mallocCUDAHost((void*)h_dem_y, sizeof(float) * newblokline * tempDemCols); mallocCUDAHost((void*)h_dem_z, sizeof(float) * newblokline * tempDemCols); mallocCUDAHost((void*)h_demsloper_x, sizeof(float) * newblokline * tempDemCols); mallocCUDAHost((void*)h_demsloper_y, sizeof(float) * newblokline * tempDemCols); mallocCUDAHost((void*)h_demsloper_z, sizeof(float) * newblokline * tempDemCols); mallocCUDAHost((void*)h_demsloper_angle, sizeof(float) * blokline * tempDemCols);//7 mallocCUDAHost((void*)h_dem_theta, sizeof(float) * newblokline * tempDemCols); mallocCUDAHost((void*)h_dem_phi, sizeof(float) * newblokline * tempDemCols); //9 mallocCUDAHost((void*)h_R, sizeof(float) * newblokline * tempDemCols); mallocCUDAHost((void*)h_localangle, sizeof(float) * newblokline * tempDemCols);//11 mallocCUDAHost((void*)h_RstX, sizeof(float) * newblokline * tempDemCols); mallocCUDAHost((void*)h_RstY, sizeof(float) * newblokline * tempDemCols); mallocCUDAHost((void*)h_RstZ, sizeof(float) * newblokline * tempDemCols);//14 mallocCUDAHost((void*)h_sigma0, sizeof(float) * newblokline * tempDemCols); mallocCUDAHost((void*)h_TransAnt, sizeof(float) * newblokline * tempDemCols); mallocCUDAHost((void*)h_ReciveAnt, sizeof(float) * newblokline * tempDemCols);//17 mallocCUDAHost((void*)h_echoAmp, sizeof(cuComplex) * newblokline * tempDemCols);//19 mallocCUDAHost((void*)h_FreqID, sizeof(long) * newblokline * tempDemCols);//20 mallocCUDADevice((void*)d_dem_x, sizeof(float) * newblokline * tempDemCols); mallocCUDADevice((void*)d_dem_y, sizeof(float) * newblokline * tempDemCols); mallocCUDADevice((void*)d_dem_z, sizeof(float) * newblokline * tempDemCols); mallocCUDADevice((void*)d_demsloper_x, sizeof(float) * newblokline * tempDemCols); mallocCUDADevice((void*)d_demsloper_y, sizeof(float) * newblokline * tempDemCols); mallocCUDADevice((void*)d_demsloper_z, sizeof(float) * newblokline * tempDemCols);//6 mallocCUDADevice((void*)d_demsloper_angle, sizeof(float) * newblokline * tempDemCols);//7 mallocCUDADevice((void*)d_dem_theta, sizeof(float) * newblokline * tempDemCols); mallocCUDADevice((void*)d_dem_phi, sizeof(float) * newblokline * tempDemCols);// 9 mallocCUDADevice((void*)d_R, sizeof(float) * newblokline * tempDemCols); mallocCUDADevice((void*)d_localangle, sizeof(float) * newblokline * tempDemCols);//11 mallocCUDADevice((void*)d_RstX, sizeof(float) * newblokline * tempDemCols); mallocCUDADevice((void*)d_RstY, sizeof(float) * newblokline * tempDemCols); mallocCUDADevice((void*)d_RstZ, sizeof(float) * newblokline * tempDemCols);//14 mallocCUDADevice((void*)d_sigma0, sizeof(float) * newblokline * tempDemCols); mallocCUDADevice((void*)d_TransAnt, sizeof(float) * newblokline * tempDemCols); mallocCUDADevice((void*)d_ReciveAnt, sizeof(float) * newblokline * tempDemCols);//17 mallocCUDADevice((void*)d_echoAmp, sizeof(cuComplex) * newblokline * tempDemCols); //19 mallocCUDADevice((void*)d_FreqID, sizeof(long) * newblokline * tempDemCols); //20 } for (long i = 0; i < newblokline; i++) { for (long j = 0; j < demxyz.width; j++) { h_dem_x[i * demxyz.width + j] = dem_x(i, j); h_dem_y[i * demxyz.width + j] = dem_y(i, j); h_dem_z[i * demxyz.width + j] = dem_z(i, j); h_demsloper_x[i * demxyz.width + j] = demsloper_x(i, j); h_demsloper_y[i * demxyz.width + j] = demsloper_y(i, j); h_demsloper_z[i * demxyz.width + j] = demsloper_z(i, j); h_demsloper_angle[i * demxyz.width + j] = sloperAngle(i, j); } } HostToDevice((void*)h_dem_x, (void*)d_dem_x, sizeof(float) * newblokline * tempDemCols); // 复制 机器 -> GPU HostToDevice((void*)h_dem_y, (void*)d_dem_y, sizeof(float) * newblokline * tempDemCols); HostToDevice((void*)h_dem_z, (void*)d_dem_z, sizeof(float) * newblokline * tempDemCols); HostToDevice((void*)h_demsloper_x, (void*)d_demsloper_x, sizeof(float) * newblokline * tempDemCols); HostToDevice((void*)h_demsloper_y, (void*)d_demsloper_y, sizeof(float) * newblokline * tempDemCols); HostToDevice((void*)h_demsloper_z, (void*)d_demsloper_z, sizeof(float) * newblokline * tempDemCols); HostToDevice((void*)h_demsloper_angle, (void*)d_demsloper_angle, sizeof(float) * newblokline * tempDemCols); HostToDevice((void*)h_dem_theta, (void*)d_dem_theta, sizeof(float) * blokline * tempDemCols); HostToDevice((void*)h_dem_phi, (void*)d_dem_phi, sizeof(float) * blokline * tempDemCols); long pixelcount = newblokline * tempDemCols; for (long prfid = 0; prfid < pluseCount; prfid++) { {// 计算 // 天线位置 float antpx = sateOirbtNodes[prfid].Px; float antpy = sateOirbtNodes[prfid].Py; float antpz = sateOirbtNodes[prfid].Pz; float antvx = sateOirbtNodes[prfid].Vx; float antvy = sateOirbtNodes[prfid].Vy; float antvz = sateOirbtNodes[prfid].Vz; //6 float antdirectx = sateOirbtNodes[prfid].AntDirecX; float antdirecty = sateOirbtNodes[prfid].AntDirecY; float antdirectz = sateOirbtNodes[prfid].AntDirecZ; // 9 天线指向 float antXaxisX = sateOirbtNodes[prfid].AntXaxisX; float antXaxisY = sateOirbtNodes[prfid].AntXaxisY; float antXaxisZ = sateOirbtNodes[prfid].AntXaxisZ;//12 天线坐标系 float antYaxisX = sateOirbtNodes[prfid].AntYaxisX; float antYaxisY = sateOirbtNodes[prfid].AntYaxisY; float antYaxisZ = sateOirbtNodes[prfid].AntYaxisZ;//15 float antZaxisX = sateOirbtNodes[prfid].AntZaxisX; float antZaxisY = sateOirbtNodes[prfid].AntZaxisY; float antZaxisZ = sateOirbtNodes[prfid].AntZaxisZ;//18 make_VectorA_B(antpx, antpy, antpz, d_dem_x, d_dem_y, d_dem_z, d_RstX, d_RstY, d_RstZ, pixelcount); // Rst = Rs - Rt; Norm_Vector(d_RstX, d_RstY, d_RstZ,d_R,pixelcount); // R cosAngle_VA_AB(d_RstX, d_RstY, d_RstZ, d_demsloper_x, d_demsloper_y, d_demsloper_z,d_localangle ,pixelcount); // 局部入射角 SatelliteAntDirectNormal(d_RstX, d_RstY, d_RstZ,antXaxisX, antXaxisY, antXaxisZ,antYaxisX, antYaxisY, antYaxisZ,antZaxisX, antZaxisY, antZaxisZ,antdirectx, antdirecty, antdirectz,d_dem_theta, d_dem_phi , pixelcount);// 计算角度 DeviceToHost(h_dem_theta, d_dem_theta, sizeof(float)*pixelcount); // 从GPU -> 主机 DeviceToHost(h_dem_phi, d_dem_phi, sizeof(float)*pixelcount); DeviceToHost(h_localangle, d_localangle, sizeof(float)*pixelcount); for (long ii =0; ii getGainLearThetaPhi(h_dem_theta[ii*tempDemCols+ jj], h_dem_phi[ii * tempDemCols + jj] ); } } for (long ii = 0; ii < blokline; ii++) { // 计算接收天线方向图 for (long jj = 0; jj < tempDemCols; jj++) { h_ReciveAnt[ii * tempDemCols + jj] = ReceivePattern->getGainLearThetaPhi(h_dem_theta[ii * tempDemCols + jj], h_dem_phi[ii * tempDemCols + jj]); } } for (long ii = 0; ii < blokline; ii++) {// 后向散射图 for (long jj = 0; jj < tempDemCols; jj++) { h_sigma0[ii * tempDemCols + jj] = this->SigmaDatabasePtr->getAmp(landcover(ii,jj), h_localangle[ii*tempDemCols+jj] * r2d, polartype); } } HostToDevice((void*)h_sigma0, (void*)d_sigma0, sizeof(float)* pixelcount); HostToDevice((void*)h_TransAnt, (void*)d_TransAnt, sizeof(float)* pixelcount); HostToDevice((void*)h_ReciveAnt, (void*)d_ReciveAnt, sizeof(float)* pixelcount); // 计算回波 calculationEcho(d_sigma0, d_TransAnt, d_ReciveAnt,d_localangle, d_R, d_demsloper_angle,NearRange, Fs, Pt, lamda, PlusePoint,d_echoAmp, d_FreqID, pixelcount); DeviceToHost(h_echoAmp, d_echoAmp, sizeof(float)* pixelcount); DeviceToHost(h_FreqID, d_FreqID, sizeof(long)* pixelcount); // 回波存档 for (long i = 0; i < pixelcount; i++) { echo.get()[i*PlusePoint+ h_FreqID[i]] =std::complex(h_echoAmp[i].x, h_echoAmp[i].y); } std::cout << "\r"<<"dem:\t"<EchoSimulationData->saveEchoArr(echo, 0, PluseCount); this->EchoSimulationData->saveToXml(); return ErrorCode::SUCCESS; } void RTPCProcessMain(long num_thread, QString TansformPatternFilePath, QString ReceivePatternFilePath, QString simulationtaskName, QString OutEchoPath, QString GPSXmlPath, QString TaskXmlPath, QString demTiffPath, QString LandCoverPath, QString HHSigmaPath, QString HVSigmaPath, QString VHSigmaPath, QString VVSigmaPath) { std::vector TansformPatternGainpoints = ReadGainFile(TansformPatternFilePath); std::shared_ptr TansformPatternGainPtr = CreateAbstractRadiationPattern(TansformPatternGainpoints); std::vector ReceivePatternGainpoints = ReadGainFile(ReceivePatternFilePath); std::shared_ptr ReceivePatternGainPtr = CreateAbstractRadiationPattern(ReceivePatternGainpoints); std::shared_ptr < AbstractSARSatelliteModel> task = ReadSimulationSettingsXML(TaskXmlPath); if (nullptr == task) { return; } else { // 打印参数 qDebug() << "--------------Task Seting ---------------------------------------"; qDebug() << "SARImageStartTime: " << task->getSARImageStartTime(); qDebug() << "SARImageEndTime: " << task->getSARImageEndTime(); qDebug() << "BandWidth: " << task->getBandWidth(); qDebug() << "CenterFreq: " << task->getCenterFreq(); qDebug() << "PRF: " << task->getPRF(); qDebug() << "Fs: " << task->getFs(); qDebug() << "POLAR: " << task->getPolarType(); qDebug() << "NearRange: " << task->getNearRange(); qDebug() << "FarRange: " << task->getFarRange(); qDebug() << (task->getFarRange() - task->getNearRange()) * 2 / LIGHTSPEED * task->getFs(); qDebug() << "\n\n"; } // 1.2 设置天线方向图 task->setTransformRadiationPattern(TansformPatternGainPtr); task->setReceiveRadiationPattern(ReceivePatternGainPtr); //2. 读取GPS节点 std::vector nodes; ErrorCode stateCode = ReadSateGPSPointsXML(GPSXmlPath, nodes); if (stateCode != ErrorCode::SUCCESS) { qWarning() << QString::fromStdString(errorCode2errInfo(stateCode)); return; } else {} std::shared_ptr SatelliteOribtModel = CreataPolyfitSatelliteOribtModel(nodes, task->getSARImageStartTime(), 3); // 以成像开始时间作为 时间参考起点 SatelliteOribtModel->setbeamAngle(task->getCenterLookAngle(), task->getIsRightLook()); // 设置天线方向图 if (nullptr == SatelliteOribtModel) { return; } else { task->setSatelliteOribtModel(SatelliteOribtModel); } qDebug() << "-------------- RTPC init ---------------------------------------"; RTPCProcessCls rtpc; rtpc.setTaskSetting(task); //qDebug() << "setTaskSetting"; rtpc.setTaskFileName(simulationtaskName); //qDebug() << "setTaskFileName"; rtpc.setDEMTiffPath(demTiffPath); //qDebug() << "setDEMTiffPath"; rtpc.setLandCoverPath(LandCoverPath); //qDebug() << "setLandCoverPath"; rtpc.setHHSigmaPath(HHSigmaPath); //qDebug() << "setHHSigmaPath"; rtpc.setHVSigmaPath(HVSigmaPath); //qDebug() << "setHVSigmaPath"; rtpc.setVHSigmaPath(VHSigmaPath); //qDebug() << "setVHSigmaPath"; rtpc.setVVSigmaPath(VVSigmaPath); //qDebug() << "setVVSigmaPath"; rtpc.setOutEchoPath(OutEchoPath); //qDebug() << "setOutEchoPath"; qDebug() << "-------------- RTPC start---------------------------------------"; rtpc.Process(num_thread); // 处理程序 qDebug() << "-------------- RTPC end---------------------------------------"; }