You can not select more than 25 topics
Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.
800 lines
25 KiB
800 lines
25 KiB
#include <omp.h> // For OpenMP parallelization
|
|
#include <mutex> // For thread-safe access
|
|
#include <sstream> // For building file names
|
|
#include <iomanip> // For zero-padded formatting
|
|
#include <filesystem> // For file existence check
|
|
#include <iostream> // For console output
|
|
#include <vector> // For dynamic arrays
|
|
#include <array> // For fixed-size field vectors
|
|
#include <string> // For file and path strings
|
|
#include <cstring> // required for strncpy
|
|
|
|
#include "vtr.h" // real vector
|
|
#include "cvtr.h" // complex vector
|
|
#include "Constants.h"
|
|
|
|
#include <sys/wait.h>
|
|
#include <sys/times.h>
|
|
|
|
#include "rapidcsv.h" // For easy CSV reading
|
|
#include "OutputSurfMesh.h" // For tri file operations
|
|
|
|
#include <complex> // For complex numbers
|
|
|
|
|
|
#include <fstream> // For file operations
|
|
#include <cstdlib> // For std::exit
|
|
#include <algorithm> // For std::copy
|
|
#include <tuple> // For std::tuple
|
|
|
|
|
|
#include <fftw3.h> // FFT header file
|
|
|
|
|
|
|
|
namespace fs = std::filesystem;
|
|
|
|
// Define types for 3D field vectors and time series containers
|
|
using FieldSeries3D = std::vector<std::vector<vtr>>; // [probe][time][field polarization]
|
|
using FFTSeries3D = std::vector<std::vector<cVtr>>; // [probe][freq][field polarization]
|
|
using current3D = std::vector<std::vector<std::vector<cVtr>>>; // [freq][traingle][nodes]
|
|
using current2D = std::vector<std::vector<cVtr>>; // [triangle][nodes]
|
|
|
|
// Count number of probe files in the folder
|
|
int countCSVFiles(const std::string& folder)
|
|
{
|
|
int count = 0;
|
|
for (const auto& entry : fs::directory_iterator(folder))
|
|
{
|
|
if (!entry.is_regular_file())
|
|
continue;
|
|
|
|
if (entry.path().extension() == ".csv")
|
|
++count;
|
|
}
|
|
return count;
|
|
}
|
|
|
|
// Function to load electric and magnetic field data from CSV files over time
|
|
void loadSeparatedProbeFields(
|
|
const std::string& folder, // Directory containing CSV files
|
|
const std::string& baseName, // Base file name (e.g., "Probes_sphere")
|
|
int startTimeStep, // Index of the first time step (e.g., 0)
|
|
int timeStepIncrement, // Step size between time steps (e.g., 1)
|
|
int nFiles, // Number of time steps to read (inclusive)
|
|
std::vector<std::vector<vtr>>& E_fieldSeries, // Output: E-field [probe][time]
|
|
std::vector<std::vector<vtr>>& H_fieldSeries // Output: H-field [probe][time]
|
|
)
|
|
{
|
|
std::mutex mtx;
|
|
int numProbes = 0;
|
|
bool initialized = false;
|
|
|
|
|
|
|
|
#pragma omp parallel for
|
|
for (int t = 0; t < nFiles; ++t)
|
|
{
|
|
int Time = startTimeStep + t * timeStepIncrement;
|
|
|
|
std::stringstream filename;
|
|
if (Time < 10000)
|
|
filename << folder << "/" << baseName << "_"
|
|
<< std::setfill('0') << std::setw(4) << Time << ".csv";
|
|
else
|
|
filename << folder << "/" << baseName << "_" << Time << ".csv";
|
|
|
|
if (!fs::exists(filename.str()))
|
|
{
|
|
#pragma omp critical
|
|
std::cerr << "[Warning] File not found: " << filename.str() << std::endl;
|
|
continue;
|
|
}
|
|
|
|
rapidcsv::Document doc(filename.str(), rapidcsv::LabelParams(0, -1));
|
|
int numRows = doc.GetRowCount();
|
|
|
|
for (int i = 0; i < numRows; ++i)
|
|
{
|
|
std::vector<float> row = doc.GetRow<float>(i);
|
|
|
|
if (row.size() < 6)
|
|
{
|
|
#pragma omp critical
|
|
std::cerr << "[Error] Row " << i << " in " << filename.str()
|
|
<< " has fewer than 6 columns.\n";
|
|
continue;
|
|
}
|
|
|
|
vtr E(row[0], row[1], row[2]);
|
|
vtr H(row[3], row[4], row[5]);
|
|
|
|
// Ordered, thread-safe assignment
|
|
E_fieldSeries[i][t] = E;
|
|
H_fieldSeries[i][t] = H;
|
|
}
|
|
}
|
|
|
|
// Print summary
|
|
if (!E_fieldSeries.empty())
|
|
{
|
|
std::cout << "Loaded E/H fields for " << E_fieldSeries.size()
|
|
<< " probes across " << E_fieldSeries[0].size() << " time steps.\n";
|
|
}
|
|
else
|
|
{
|
|
std::cerr << "[Error] No field data loaded.\n";
|
|
}
|
|
}
|
|
|
|
|
|
/*
|
|
Basic Workflow for FFTW3
|
|
|
|
1. Allocate input/output memory
|
|
--> Use fftwf_malloc() for aligned memory (important for SIMD)
|
|
|
|
2. Create a plan
|
|
A plan tells FFTW how to compute the FFT efficiently for a specific input size and stride.
|
|
--> fftwf_plan plan = fftwf_plan_dft_r2c_1d(N, in, out, FFTW_ESTIMATE);
|
|
|
|
3. Execute the plan
|
|
fftwf_execute(plan);
|
|
|
|
4. Destroy the plan
|
|
fftwf_destroy_plan(plan);
|
|
|
|
5. Free the Memory
|
|
fftwf_free(in); fftwf_free(out);
|
|
|
|
*/
|
|
|
|
#include <iostream>
|
|
#include <iomanip>
|
|
|
|
void computeFieldFFT(const FieldSeries3D& timeSeries, FFTSeries3D& fftSeries)
|
|
{
|
|
int numProbes = timeSeries.size();
|
|
if (numProbes == 0) return;
|
|
|
|
int numTimeSamples = timeSeries[0].size();
|
|
|
|
int Npad = 4 * numTimeSamples; // zero-padding factor (e.g., 4x)
|
|
int numFreq = Npad / 2 + 1;
|
|
|
|
fftSeries.resize(numProbes);
|
|
|
|
fftwf_init_threads();
|
|
fftwf_plan_with_nthreads(1); // Use 1 thread per plan creation
|
|
|
|
for (int p = 0; p < numProbes; ++p)
|
|
{
|
|
fftSeries[p].resize(numFreq);
|
|
|
|
float* inX = (float*)fftwf_malloc(sizeof(float) * Npad);
|
|
float* inY = (float*)fftwf_malloc(sizeof(float) * Npad);
|
|
float* inZ = (float*)fftwf_malloc(sizeof(float) * Npad);
|
|
|
|
fftwf_complex* outX = (fftwf_complex*)fftwf_malloc(sizeof(fftwf_complex) * numFreq);
|
|
fftwf_complex* outY = (fftwf_complex*)fftwf_malloc(sizeof(fftwf_complex) * numFreq);
|
|
fftwf_complex* outZ = (fftwf_complex*)fftwf_malloc(sizeof(fftwf_complex) * numFreq);
|
|
|
|
for (int t = 0; t < numTimeSamples; ++t)
|
|
{
|
|
inX[t] = timeSeries[p][t].getx();
|
|
inY[t] = timeSeries[p][t].gety();
|
|
inZ[t] = timeSeries[p][t].getz();
|
|
}
|
|
for (int t = numTimeSamples; t < Npad; ++t)
|
|
{
|
|
inX[t] = 0.0f;
|
|
inY[t] = 0.0f;
|
|
inZ[t] = 0.0f;
|
|
}
|
|
|
|
fftwf_plan planX = fftwf_plan_dft_r2c_1d(Npad, inX, outX, FFTW_ESTIMATE);
|
|
fftwf_plan planY = fftwf_plan_dft_r2c_1d(Npad, inY, outY, FFTW_ESTIMATE);
|
|
fftwf_plan planZ = fftwf_plan_dft_r2c_1d(Npad, inZ, outZ, FFTW_ESTIMATE);
|
|
|
|
fftwf_execute(planX);
|
|
fftwf_execute(planY);
|
|
fftwf_execute(planZ);
|
|
|
|
for (int f = 0; f < numFreq; ++f)
|
|
{
|
|
cVtr cv;
|
|
cv.setcvtr(
|
|
Complex(outX[f][0], outX[f][1]),
|
|
Complex(outY[f][0], outY[f][1]),
|
|
Complex(outZ[f][0], outZ[f][1])
|
|
);
|
|
fftSeries[p][f] = cv;
|
|
}
|
|
|
|
fftwf_destroy_plan(planX);
|
|
fftwf_destroy_plan(planY);
|
|
fftwf_destroy_plan(planZ);
|
|
fftwf_free(inX); fftwf_free(inY); fftwf_free(inZ);
|
|
fftwf_free(outX); fftwf_free(outY); fftwf_free(outZ);
|
|
|
|
// Update textual progress
|
|
float progress = (float)(p + 1) / numProbes * 100.0f;
|
|
std::cout << "\r[";
|
|
int barWidth = 50;
|
|
int pos = progress / 100.0f * barWidth;
|
|
for (int i = 0; i < barWidth; ++i)
|
|
std::cout << (i < pos ? "=" : (i == pos ? ">" : " "));
|
|
std::cout << "] " << std::fixed << std::setprecision(1) << progress << "%";
|
|
std::cout.flush();
|
|
}
|
|
|
|
std::cout << "\nDone computing FFTs." << std::endl;
|
|
fftwf_cleanup_threads();
|
|
}
|
|
|
|
|
|
|
|
|
|
struct BandAnalysisResult {
|
|
int idx_min;
|
|
int idx_max;
|
|
double fmin_Hz;
|
|
double fmax_Hz;
|
|
int num_points;
|
|
bool valid;
|
|
};
|
|
|
|
BandAnalysisResult analyzeFFTAmplitudeBand(
|
|
const std::vector<cVtr>& fftProbe,
|
|
const std::vector<double>& freqAxis
|
|
) {
|
|
int N = fftProbe.size();
|
|
std::vector<double> magnitudes(N);
|
|
|
|
double peak = 0.0;
|
|
|
|
|
|
|
|
for (int i = 0; i < N; ++i)
|
|
{
|
|
Complex cx = fftProbe[i].getx();
|
|
Complex cy = fftProbe[i].gety();
|
|
Complex cz = fftProbe[i].getz();
|
|
|
|
double mag = std::sqrt(
|
|
std::norm(cx) + std::norm(cy) + std::norm(cz)
|
|
);
|
|
magnitudes[i] = mag;
|
|
if (mag > peak) peak = mag;
|
|
}
|
|
|
|
|
|
std::ofstream fout("fft_spectrum.csv");
|
|
fout << "Frequency_Hz,Magnitude\n";
|
|
for (int i = 0; i < N; ++i)
|
|
{
|
|
fout << freqAxis[i] << "," << magnitudes[i] << "\n";
|
|
}
|
|
fout.close();
|
|
std::cout << "[Info] FFT spectrum saved to fft_spectrum.csv\n";
|
|
|
|
|
|
|
|
double threshold = 0.05 * peak;
|
|
int idx_min = -1, idx_max = -1;
|
|
|
|
for (int i = 0; i < N; ++i) {
|
|
if (magnitudes[i] > threshold) {
|
|
if (idx_min == -1) idx_min = i;
|
|
idx_max = i;
|
|
}
|
|
}
|
|
|
|
BandAnalysisResult result;
|
|
if (idx_min != -1 && idx_max != -1)
|
|
{
|
|
result.idx_min = idx_min;
|
|
result.idx_max = idx_max;
|
|
result.fmin_Hz = freqAxis[idx_min];
|
|
result.fmax_Hz = freqAxis[idx_max];
|
|
result.num_points = idx_max - idx_min + 1;
|
|
result.valid = true;
|
|
|
|
std::cout << "Frequency Band: " << result.fmin_Hz / 1e9 << " GHz to " << result.fmax_Hz / 1e9 << " GHz\n";
|
|
std::cout << "Number of points in band: " << result.num_points << "\n";
|
|
std::cout << "Peak Magnitude: " << peak << "\n";
|
|
|
|
|
|
}
|
|
else
|
|
{
|
|
std::cout << "No significant frequency band found.\n";
|
|
result = { -1, -1, 0.0, 0.0, 0, false };
|
|
}
|
|
|
|
|
|
|
|
return result;
|
|
}
|
|
|
|
|
|
|
|
|
|
std::vector<double> computeFreqAxis(int numTimeSamples, double dt)
|
|
{
|
|
int Npad = 4 * numTimeSamples;
|
|
int numFreq = Npad / 2 + 1;
|
|
std::vector<double> freqAxis(numFreq);
|
|
double fs = 1.0 / dt; // Sampling frequency
|
|
|
|
for (int k = 0; k < numFreq; ++k)
|
|
freqAxis[k] = k * fs / Npad; // f_k = k * (fs / N)
|
|
|
|
return freqAxis;
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
void writeCurrents(const std::string& filename, const current2D& current)
|
|
{
|
|
|
|
int Num_tri = current.size(); // Number of triangles
|
|
|
|
std::ofstream fid(filename);
|
|
for (int t = 0; t < Num_tri; ++t)
|
|
{
|
|
for (int n = 0; n < 3; ++n)
|
|
{
|
|
cVtr cv = current[t][n]; // Get the first node's current
|
|
fid << cv.getx().real() << " " << cv.getx().imag() << "\n"; // Write real and imaginary parts
|
|
fid << cv.gety().real() << " " << cv.gety().imag() << "\n"; // Write real and imaginary parts
|
|
fid << cv.getz().real() << " " << cv.getz().imag() << "\n"; // Write real and imaginary parts
|
|
}
|
|
fid << "\n"; // Newline after each node
|
|
}
|
|
fid.close();
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
std::pair<double, double> process_RCS(
|
|
const current2D& F_curJ,
|
|
const current2D& F_curM,
|
|
double FREQ,
|
|
const std::string& FREQ_STR,
|
|
const std::string& project_name,
|
|
const std::string& outfile,
|
|
double THETA_sc,
|
|
double PHI_sc,
|
|
int writeMode,
|
|
int output_opt)
|
|
{
|
|
|
|
std::string current_dir = fs::current_path();
|
|
std::string freq_dir = current_dir + "/freq/FREQ" + FREQ_STR;
|
|
fs::create_directories(freq_dir);
|
|
std::cout << "Working Folder = " << freq_dir << std::endl;
|
|
fs::current_path(freq_dir);
|
|
|
|
// Write .curJ and .curM
|
|
writeCurrents(outfile + ".curJ", F_curJ);
|
|
writeCurrents(outfile + ".curM", F_curM);
|
|
|
|
// Write .region file
|
|
std::ofstream reg(outfile + ".region");
|
|
reg << "1\n";
|
|
reg << "0 FEBI\n";
|
|
reg << "0 FEBI " << outfile << "\n";
|
|
reg << "Coupling\n";
|
|
reg.close();
|
|
|
|
// Copy geometry and binary
|
|
fs::copy(current_dir + "/" + project_name + "_out.tri", freq_dir + "/" + outfile + ".tri", fs::copy_options::overwrite_existing);
|
|
fs::copy(current_dir + "/n2f_main", freq_dir + "/n2f_main", fs::copy_options::overwrite_existing);
|
|
|
|
// Run far-field computation
|
|
// Projectname freq_in_MHz centerX centerY centerZ numTheta numPhi mag_or_imaginary radiation_or_scattering
|
|
if (writeMode == 0)
|
|
{
|
|
std::string cmd_n2f = "./n2f_main " + outfile + " " + std::to_string(int(FREQ)) + " 0 0 0 180 360 " + std::to_string(output_opt) + " 1";
|
|
std::system(cmd_n2f.c_str());
|
|
} else if (writeMode == 1)
|
|
{
|
|
std::string cmd_n2f = "./n2f_main " + outfile + " " + std::to_string(int(FREQ)) + " 0 0 0 180 4 " + std::to_string(output_opt) + " 1";
|
|
std::system(cmd_n2f.c_str());
|
|
} else if (writeMode == 2)
|
|
{
|
|
std::string cmd_n2f = "./n2f_main " + outfile + " " + std::to_string(int(FREQ)) + " 0 0 0 2 360 " + std::to_string(output_opt) + " 1";
|
|
std::system(cmd_n2f.c_str());
|
|
}
|
|
|
|
std::system(("emsurftranslator -s " + outfile).c_str());
|
|
|
|
std::cout << "\n\n==== Process Farfield Result ====\n\n";
|
|
std::system(("sed -i 's/ / /g' " + outfile + ".cs").c_str());
|
|
|
|
// Parse .cs file
|
|
std::ifstream csfile(outfile + ".cs");
|
|
std::vector<std::tuple<double, double, double, double>> rcs_data;
|
|
std::string line;
|
|
std::getline(csfile, line); // skip header
|
|
while (std::getline(csfile, line)) {
|
|
std::istringstream iss(line);
|
|
double theta, phi, rcs_v, rcs_h;
|
|
if (iss >> theta >> phi >> rcs_v >> rcs_h) {
|
|
rcs_data.emplace_back(theta, phi, rcs_v, rcs_h);
|
|
}
|
|
}
|
|
csfile.close();
|
|
|
|
fs::current_path(current_dir); // Go back to root dir
|
|
|
|
auto it = std::find_if(rcs_data.begin(), rcs_data.end(), [&](auto& row) {
|
|
return std::get<0>(row) == THETA_sc && std::get<1>(row) == PHI_sc;
|
|
});
|
|
|
|
if (it != rcs_data.end()) {
|
|
return {std::get<2>(*it), std::get<3>(*it)};
|
|
} else {
|
|
std::cerr << "Requested (theta, phi) RCS not found.\n";
|
|
return {0.0, 0.0};
|
|
}
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
// Sample usage and test
|
|
int main(int argc, char *argv[])
|
|
{
|
|
|
|
// Computation Time Start and End
|
|
tms tmsStart, tmsEnd;
|
|
times(&tmsStart);
|
|
|
|
// Project name
|
|
char prjName[180];
|
|
|
|
// Scattering directions
|
|
int Theta_sc, Phi_sc;
|
|
|
|
int output_opt; // Output format
|
|
double dt; // Simulation mode
|
|
int step;
|
|
|
|
int writeMode; // writeMode 0: Full | 1: Principal Phi cuts = 0 cut / 90 cut FULL | 2: Theta = 90 cut
|
|
|
|
|
|
|
|
// Read input parameters
|
|
if (argc == 8)
|
|
{
|
|
strcpy( prjName, argv[1] );
|
|
dt = atof(argv[2]);
|
|
step = atoi(argv[3]);
|
|
Theta_sc = atoi(argv[4]);
|
|
Phi_sc = atoi(argv[5]);
|
|
writeMode = atoi(argv[6]);
|
|
output_opt = atoi(argv[7]);
|
|
}
|
|
else
|
|
{
|
|
std::cout << "Usage:\n";
|
|
std::cout << " ./program <ProjectName> <dt> <Theta_sc> <Phi_sc> <writeMode> <output_opt>\n\n";
|
|
|
|
std::cout << "Where:\n";
|
|
std::cout << " <ProjectName> = Name of the project (string)\n";
|
|
std::cout << " <dt> = Time interval in seconds\n";
|
|
std::cout << " <step> = Step in between samples\n";
|
|
std::cout << " <Theta_sc> = Scattering angle in theta (deg)\n";
|
|
std::cout << " <Phi_sc> = Scattering angle in phi (deg)\n";
|
|
std::cout << " <writeMode> = 0 (Full), 1 (Principal Theta cuts), 2 (Principal Phi cuts)\n";
|
|
std::cout << " <output_opt> = 0 (magnitude), 1 (real & imag)\n";
|
|
|
|
exit(1);
|
|
}
|
|
|
|
|
|
// Set center of origin
|
|
vtr origin;
|
|
origin.setvtr( 0.0, 0.0, 0.0 );
|
|
|
|
|
|
|
|
|
|
OutputSurfMesh outputMesh;
|
|
|
|
// Read the mesh from .tri file
|
|
char triName[256];
|
|
sprintf(triName, "./%s_out.tri", prjName);
|
|
std::cout << "--------------------" << std::endl;
|
|
std::cout << "Reading Tri surface mesh " << triName << std::endl;
|
|
outputMesh.readFromFile(triName);
|
|
std::cout << "--------------------" << std::endl;
|
|
std::cout << "Compute Normals " << std::endl;
|
|
outputMesh.computeTriangleNormals();
|
|
std::cout << "--------------------" << std::endl;
|
|
outputMesh.printSummary();
|
|
std::cout << "--------------------" << std::endl;
|
|
|
|
// Store the nodes of the triangles
|
|
int num_nodes = outputMesh.num_nodes;
|
|
int numFiles = countCSVFiles("./PROBES_sc");
|
|
std::cout << "Number of Nodes = " << num_nodes << std::endl;
|
|
std::cout << "Number of Files = " << numFiles << std::endl;
|
|
|
|
|
|
|
|
|
|
// ------------------------------ //
|
|
// READ FIELDS AT THE PROBES //
|
|
// ------------------------------ //
|
|
|
|
// Read Scattered Field
|
|
FieldSeries3D E_fields_sc(num_nodes, std::vector<vtr>(numFiles));
|
|
FieldSeries3D H_fields_sc(num_nodes, std::vector<vtr>(numFiles));
|
|
|
|
loadSeparatedProbeFields(
|
|
"./PROBES_sc", // Folder
|
|
"Probes_sphere", // Base name
|
|
0, // Start time
|
|
step, // Increment
|
|
numFiles, // Number of files
|
|
E_fields_sc,
|
|
H_fields_sc
|
|
);
|
|
|
|
|
|
// Read Incident Field
|
|
FieldSeries3D E_fields_inc(num_nodes, std::vector<vtr>(numFiles));
|
|
FieldSeries3D H_fields_inc(num_nodes, std::vector<vtr>(numFiles));
|
|
loadSeparatedProbeFields(
|
|
"./PROBES_inc", // Folder
|
|
"Probes_sphere", // Base name
|
|
0, // Start time
|
|
step, // Increment
|
|
numFiles, // Number of files
|
|
E_fields_inc,
|
|
H_fields_inc
|
|
);
|
|
|
|
|
|
// -------------------------------------- //
|
|
// Perform FFT of time domain fields //
|
|
// -------------------------------------- //
|
|
FFTSeries3D E_fields_sc_fft, H_fields_sc_fft;
|
|
FFTSeries3D E_fields_inc_fft;
|
|
E_fields_sc_fft.resize(num_nodes);
|
|
H_fields_sc_fft.resize(num_nodes);
|
|
E_fields_inc_fft.resize(num_nodes);
|
|
|
|
std::cout << "Performing FFT on scattered Efields..." << std::endl;
|
|
computeFieldFFT(E_fields_sc, E_fields_sc_fft);
|
|
std::cout << "Performing FFT on scattered Hfields..." << std::endl;
|
|
computeFieldFFT(H_fields_sc, H_fields_sc_fft);
|
|
std::cout << "Performing FFT on Incident Efields..." << std::endl;
|
|
computeFieldFFT(E_fields_inc, E_fields_inc_fft);
|
|
|
|
|
|
int numTime = E_fields_inc[0].size();
|
|
std::ofstream fout("probe0_incident.csv");
|
|
fout << "TimeStep,Ex,Ey,Ez\n";
|
|
for (int t = 0; t < numTime; ++t)
|
|
{
|
|
fp_t inX = E_fields_inc[0][t].getx();
|
|
fp_t inY = E_fields_inc[0][t].gety();
|
|
fp_t inZ = E_fields_inc[0][t].getz();
|
|
fout << t << "," << inX << "," << inY << "," << inZ << "\n";
|
|
}
|
|
fout.close();
|
|
std::cout << "[Info] Probe 0 time series written to probe0_timeseries.csv\n";
|
|
|
|
|
|
// Clear the original time domain fields to free memory
|
|
E_fields_sc.clear();
|
|
H_fields_sc.clear();
|
|
E_fields_inc.clear();
|
|
H_fields_inc.clear();
|
|
|
|
// --------------------------------------------------------------- //
|
|
// Normalization between scattered fields and incident fields //
|
|
// --------------------------------------------------------------- //
|
|
|
|
std::cout << "Compute frequency axis..." << std::endl;
|
|
int numTimeSamples = E_fields_inc[0].size();
|
|
std::vector<double> freqAxis = computeFreqAxis(numTimeSamples, dt);
|
|
|
|
std::cout << "Analyzing FFT Amplitude Band..." << std::endl;
|
|
auto band_results = analyzeFFTAmplitudeBand(E_fields_inc_fft[0], freqAxis);
|
|
|
|
// Allocate normalized output fields
|
|
FFTSeries3D E_fields_sc_fft_norm(num_nodes);
|
|
FFTSeries3D H_fields_sc_fft_norm(num_nodes);
|
|
|
|
// Initialize normalized fields to zero
|
|
for (auto& v : E_fields_sc_fft_norm)
|
|
v.resize(freqAxis.size());
|
|
for (auto& v : H_fields_sc_fft_norm)
|
|
v.resize(freqAxis.size());
|
|
|
|
if (!band_results.valid)
|
|
{
|
|
std::cerr << "[Warning] Invalid FFT band — skipping normalization." << std::endl;
|
|
}
|
|
else
|
|
{
|
|
std::cout << "Normalizing scattered E-fields within band..." << std::endl;
|
|
|
|
// Parallelization across nodes and frequency bins
|
|
#pragma omp parallel for
|
|
for (int n = 0; n < num_nodes; ++n)
|
|
{
|
|
for (int f = band_results.idx_min; f <= band_results.idx_max; ++f)
|
|
{
|
|
const cVtr& E_inc = E_fields_inc_fft[n][f];
|
|
const cVtr& E_scat = E_fields_sc_fft[n][f];
|
|
const cVtr& H_scat = H_fields_sc_fft[n][f];
|
|
|
|
double inc_mag = std::sqrt(
|
|
std::norm(E_inc.getx()) +
|
|
std::norm(E_inc.gety()) +
|
|
std::norm(E_inc.getz())
|
|
);
|
|
|
|
if (inc_mag > 1e-12)
|
|
{
|
|
cVtr E_norm = E_scat / inc_mag;
|
|
cVtr H_norm = H_scat / inc_mag;
|
|
E_fields_sc_fft_norm[n][f] = E_norm;
|
|
H_fields_sc_fft_norm[n][f] = H_norm;
|
|
}
|
|
else
|
|
{
|
|
cVtr zero;
|
|
zero.setcvtr(0, 0, 0);
|
|
E_fields_sc_fft_norm[n][f] = zero;
|
|
H_fields_sc_fft_norm[n][f] = zero;
|
|
}
|
|
}
|
|
}
|
|
}
|
|
|
|
E_fields_inc_fft.clear(); // Free memory for incident fields
|
|
H_fields_sc_fft.clear(); // Free memory for scattered fields
|
|
H_fields_sc_fft.clear(); // Free memory for scattered fields
|
|
|
|
|
|
// ---------------------------------- //
|
|
// Compute the tangential fields //
|
|
// ---------------------------------- //
|
|
|
|
std::cout << "Allocating curJ and curM..." << std::endl;
|
|
|
|
// Assuming number of frequency points
|
|
int numFreq = E_fields_sc_fft_norm[0].size();
|
|
|
|
// ✅ Step 1: Allocation of 3D current arrays
|
|
current3D curJ(numFreq);
|
|
current3D curM(numFreq);
|
|
|
|
for (int f = 0; f < numFreq; ++f)
|
|
{
|
|
curJ[f].resize(outputMesh.num_triangles);
|
|
curM[f].resize(outputMesh.num_triangles);
|
|
for (int t = 0; t < outputMesh.num_triangles; ++t)
|
|
{
|
|
curJ[f][t].resize(3); // 3 local nodes
|
|
curM[f][t].resize(3);
|
|
}
|
|
}
|
|
|
|
|
|
// ✅ Step 2: Compute the tangential fields
|
|
std::cout << "Computing tangential fields..." << std::endl;
|
|
#pragma omp parallel for
|
|
for (int t = 0; t < outputMesh.num_triangles; ++t)
|
|
{
|
|
std::vector<int> tri_nodes = outputMesh.getTriangle(t);
|
|
std::vector<double> normal_d = outputMesh.getNormal(t);
|
|
vtr n(normal_d[0], normal_d[1], normal_d[2]); // unit normal
|
|
|
|
for (int j = 0; j < 3; ++j)
|
|
{
|
|
int nodeIdx = tri_nodes[j];
|
|
|
|
for (int f = 0; f < numFreq; ++f)
|
|
{
|
|
const cVtr& E = E_fields_sc_fft_norm[nodeIdx][f];
|
|
const cVtr& H = H_fields_sc_fft_norm[nodeIdx][f];
|
|
|
|
// M = -n × E = E × n
|
|
cVtr M_vec = E.cross(n);
|
|
// J = n × H = -H × n
|
|
cVtr J_vec = H.cross(n) * (-1.0f);
|
|
|
|
curM[f][t][j] = M_vec;
|
|
curJ[f][t][j] = J_vec;
|
|
}
|
|
}
|
|
}
|
|
|
|
// clear up fft data
|
|
E_fields_sc_fft_norm.clear(); // Free memory for normalized fields
|
|
H_fields_sc_fft_norm.clear(); // Free memory for normalized fields
|
|
|
|
std::cout << "Tangential fields computed." << std::endl;
|
|
|
|
|
|
|
|
// --------------------------------- //
|
|
// Fast Near to Far Computation //
|
|
// --------------------------------- //
|
|
|
|
|
|
std::vector<double> freqList;
|
|
std::vector<double> rcsVList;
|
|
std::vector<double> rcsHList;
|
|
|
|
|
|
|
|
// Compute Far Field
|
|
for (int f = band_results.idx_min; f <= band_results.idx_max; ++f)
|
|
{
|
|
std::cout << "Processing frequency "
|
|
<< std::fixed << std::setprecision(3)
|
|
<< freqAxis[f] / 1e9 << " GHz" << std::endl;
|
|
|
|
double freq = freqAxis[f] / 1e6;
|
|
std::string FREQ_STR = std::to_string(int(freq));
|
|
std::string project_name(prjName);
|
|
std::string outfile = "ff_" + FREQ_STR;
|
|
|
|
|
|
// Run Near to Far field computation
|
|
auto [rcs_v, rcs_h] = process_RCS( curJ[f], curM[f], freq, FREQ_STR, project_name, outfile, Theta_sc, Phi_sc, writeMode, output_opt);
|
|
|
|
|
|
// Collect data
|
|
freqList.push_back(freq/1e3);
|
|
rcsVList.push_back(rcs_v);
|
|
rcsHList.push_back(rcs_h);
|
|
}
|
|
|
|
// Write RCS summary to CSV file
|
|
rapidcsv::Document doc("", rapidcsv::LabelParams(0, -1));
|
|
doc.InsertColumn(0, freqList, "Frequency_MHz");
|
|
doc.InsertColumn(1, rcsVList, "RCS_V");
|
|
doc.InsertColumn(2, rcsHList, "RCS_H");
|
|
doc.Save("rcs_summary.csv");
|
|
|
|
|
|
|
|
return 0;
|
|
}
|
|
|
|
|
|
|