13#define __FUNCT__ "SetupPostProcessSwarm"
26 PetscFunctionBeginUser;
28 char *pipeline_copy, *step_token, *step_saveptr;
29 PetscBool finalize_needed = PETSC_FALSE;
31 ierr = DMCreate(PETSC_COMM_WORLD, &user->
post_swarm); CHKERRQ(ierr);
32 ierr = DMSetType(user->
post_swarm, DMSWARM); CHKERRQ(ierr);
33 ierr = DMSetDimension(user->
post_swarm, 3); CHKERRQ(ierr);
34 ierr = DMSwarmSetType(user->
post_swarm, DMSWARM_BASIC); CHKERRQ(ierr);
37 ierr = DMSwarmSetCellDM(user->
post_swarm, user->
da); CHKERRQ(ierr);
50 step_token = strtok_r(pipeline_copy,
";", &step_saveptr);
53 if (strlen(step_token) == 0) { step_token = strtok_r(NULL,
";", &step_saveptr);
continue; }
55 char *keyword = strtok(step_token,
":");
56 char *args_str = strtok(NULL,
"");
58 PetscInt output_field_dimensions = 1;
60 if (strcasecmp(keyword,
"ComputeSpecificKE") == 0) {
61 if (!args_str) SETERRQ(PETSC_COMM_SELF, 1,
"Error (ComputeSpecificKE): Missing arguments.");
62 char *inputs_str = strtok(args_str,
">");
63 char *output_field = strtok(NULL,
">");
64 output_field_dimensions = 1;
65 if (!output_field) SETERRQ(PETSC_COMM_SELF, 1,
"Error (ComputeSpecificKE): Missing output field in 'in>out' syntax.");
70 finalize_needed = PETSC_TRUE;
77 step_token = strtok_r(NULL,
";", &step_saveptr);
80 ierr = PetscFree(pipeline_copy); CHKERRQ(ierr);
83 if (finalize_needed) {
85 ierr = DMSwarmFinalizeFieldRegister(user->
post_swarm); CHKERRQ(ierr);
91 PetscFunctionReturn(0);
96#define __FUNCT__ "EulerianDataProcessingPipeline"
113 char *pipeline_copy, *step_token, *step_saveptr;
115 PetscFunctionBeginUser;
122 PetscFunctionReturn(0);
128 ierr = PetscStrallocpy(pps->
process_pipeline, &pipeline_copy); CHKERRQ(ierr);
131 step_token = strtok_r(pipeline_copy,
";", &step_saveptr);
134 if (strlen(step_token) == 0) {
135 step_token = strtok_r(NULL,
";", &step_saveptr);
139 char *keyword = strtok(step_token,
":");
140 char *args_str = strtok(NULL,
"");
143 step_token = strtok_r(NULL,
";", &step_saveptr);
150 LOG_ALLOW(
GLOBAL,
LOG_INFO,
"Executing Transformation: '%s' on args: '%s'\n", keyword, args_str ? args_str :
"None");
153 if (strcasecmp(keyword,
"CellToNodeAverage") == 0) {
154 if (!args_str) SETERRQ(PETSC_COMM_SELF, 1,
"CellToNodeAverage requires arguments in 'in_field>out_field' format.");
155 char *in_field = strtok(args_str,
">");
156 char *out_field = strtok(NULL,
">");
157 if (!in_field || !out_field) SETERRQ(PETSC_COMM_SELF, 1,
"CellToNodeAverage requires 'in>out' syntax (e.g., P>P_nodal).");
158 if(strcmp(in_field,out_field)==0) SETERRQ(PETSC_COMM_SELF, 1,
"CellToNodeAverage input and output fields must be different.");
159 if(user->
simCtx->
np == 0 && (strcmp(out_field,
"Psi_nodal")==0 || strcmp(in_field,
"Psi_nodal")==0)){
160 LOG(
GLOBAL,
LOG_WARNING,
"CellToNodeAverage cannot process 'Psi_nodal' when no particles are present in the simulation.\n");
161 step_token = strtok_r(NULL,
";", &step_saveptr);
167 else if (strcasecmp(keyword,
"ComputeQCriterion") == 0) {
170 else if (strcasecmp(keyword,
"NormalizeRelativeField") == 0) {
171 if (!args_str) SETERRQ(PETSC_COMM_SELF, 1,
"NormalizePressure requires the pressure field name (e.g., 'P') as an argument.");
180 step_token = strtok_r(NULL,
";", &step_saveptr);
183 ierr = PetscFree(pipeline_copy); CHKERRQ(ierr);
186 PetscFunctionReturn(0);
191#define __FUNCT__ "WriteEulerianFile"
223 PetscFunctionBeginUser;
227 LOG_ALLOW(
GLOBAL,
LOG_DEBUG,
"No instantaneous fields requested for output at ti=%" PetscInt_FMT
". Skipping.\n", ti);
228 PetscFunctionReturn(0);
234 ierr = PetscMemzero(&meta,
sizeof(
VTKMetaData)); CHKERRQ(ierr);
239 LOG_ALLOW(
GLOBAL,
LOG_DEBUG,
"Using coords linearization order: fast=i mid=j slow=k (sizes: %" PetscInt_FMT
" x %" PetscInt_FMT
" x %" PetscInt_FMT
")\n",
240 meta.
mx, meta.
my, meta.
mz);
244 char *fields_copy, *field_name;
247 field_name = strtok(fields_copy,
",");
250 if (!*field_name) { field_name = strtok(NULL,
",");
continue; }
254 Vec field_vec = NULL;
255 PetscInt num_components = 0;
257 if (!strcasecmp(field_name,
"P_nodal")) {
258 field_vec = user->
P_nodal; num_components = 1;
259 }
else if (!strcasecmp(field_name,
"Ucat_nodal")) {
260 field_vec = user->
Ucat_nodal; num_components = 3;
261 }
else if (!strcasecmp(field_name,
"Qcrit")) {
262 field_vec = user->
Qcrit; num_components = 1;
263 }
else if (!strcasecmp(field_name,
"Psi_nodal")){
266 field_name = strtok(NULL,
",");
269 field_vec = user->
Psi_nodal; num_components = 1;
272 field_name = strtok(NULL,
",");
278 field_name = strtok(NULL,
",");
425 field_name = strtok(NULL,
",");
428 ierr = PetscFree(fields_copy); CHKERRQ(ierr);
494 snprintf(filename,
sizeof(filename),
"%s_%05" PetscInt_FMT
".vts", pps->
output_prefix, ti);
499 PetscFunctionReturn(0);
504#define __FUNCT__ "ParticleDataProcessingPipeline"
527 char *pipeline_copy, *step_token, *step_saveptr;
529 PetscFunctionBeginUser;
534 PetscInt n_global_source;
535 ierr = DMSwarmGetSize(user->
swarm, &n_global_source); CHKERRQ(ierr);
541 PetscFunctionReturn(0);
549 step_token = strtok_r(pipeline_copy,
";", &step_saveptr);
552 if (strlen(step_token) == 0) { step_token = strtok_r(NULL,
";", &step_saveptr);
continue; }
554 char *keyword = strtok(step_token,
":");
555 char *args_str = strtok(NULL,
"");
559 LOG_ALLOW(
GLOBAL,
LOG_INFO,
"Executing Particle Transformation: '%s' on args: '%s'\n", keyword, args_str ? args_str :
"None");
561 if (strcasecmp(keyword,
"ComputeSpecificKE") == 0) {
562 char *velocity_field = strtok(args_str,
">");
563 char *ske_field = strtok(NULL,
">");
572 step_token = strtok_r(NULL,
";", &step_saveptr);
574 ierr = PetscFree(pipeline_copy); CHKERRQ(ierr);
579 PetscFunctionReturn(0);
583#define __FUNCT__ "WriteParticleFile"
592 PetscInt n_total_particles_before_subsample;
594 PetscFunctionBeginUser;
599 PetscFunctionReturn(0);
602 ierr = DMSwarmGetSize(user->
swarm, &n_global); CHKERRQ(ierr);
605 PetscFunctionReturn(0);
608 ierr = PetscMemzero(&part_meta,
sizeof(
VTKMetaData)); CHKERRQ(ierr);
616 LOG_ALLOW(
GLOBAL,
LOG_INFO,
"--- Starting VTP Particle File Writing for ti = %" PetscInt_FMT
" (writing %" PetscInt_FMT
" of %" PetscInt_FMT
" particles) ---\n",
617 ti, part_meta.
npoints, n_total_particles_before_subsample);
637 PetscFunctionReturn(0);
641#define __FUNCT__ "main"
648 ierr = PetscInitialize(&argc, &argv, (
char *)0,
"Unified Post-Processing Tool"); CHKERRQ(ierr);
674 SETERRQ(PETSC_COMM_SELF,1,
"Particle output requested but np=0. Please set np>0 to enable particle output.");
708 if(simCtx->
rank == 0){
710 PetscReal currentTime = (PetscReal)ti*simCtx->
dt;
718 if (simCtx->
rank == 0) {
719 PetscInt endTime = pps->
endTime-1;
721 PetscReal endTimeValue = (PetscReal)pps->
endTime*simCtx->
dt;
723 PetscPrintf(PETSC_COMM_SELF,
"\n");
734 ierr = PetscFinalize();
PetscErrorCode ResizeSwarmGlobally(DM swarm, PetscInt N_target)
PetscErrorCode PreCheckAndResizeSwarm(UserCtx *user, PetscInt ti, const char *ext)
Checks particle count in the reference file and resizes the swarm if needed.
PetscErrorCode InitializeParticleSwarm(SimCtx *simCtx)
Perform particle swarm initialization, particle-grid interaction, and related operations.
PetscErrorCode RegisterSwarmField(DM swarm, const char *fieldName, PetscInt fieldDim, PetscDataType dtype)
Registers a swarm field without finalizing registration.
PetscErrorCode ReadSimulationFields(UserCtx *user, PetscInt ti)
Reads binary field data for velocity, pressure, and other required vectors.
void TrimWhitespace(char *str)
Helper function to trim leading/trailing whitespace from a string.
PetscErrorCode ReadAllSwarmFields(UserCtx *user, PetscInt ti)
Reads multiple fields (positions, velocity, CellID, and weight) into a DMSwarm.
PetscInt CreateVTKFileFromMetadata(const char *filename, const VTKMetaData *meta, MPI_Comm comm)
#define LOCAL
Logging scope definitions for controlling message output.
#define GLOBAL
Scope for global logging across all processes.
#define LOG_ALLOW(scope, level, fmt,...)
Logging macro that checks both the log level and whether the calling function is in the allowed-funct...
#define PROFILE_FUNCTION_END
Marks the end of a profiled code block.
PetscErrorCode ProfilingFinalize(SimCtx *simCtx)
the profiling excercise and build a profiling summary which is then printed to a log file.
#define LOG(scope, level, fmt,...)
Logging macro for PETSc-based applications with scope control.
void PrintProgressBar(PetscInt step, PetscInt startStep, PetscInt totalSteps, PetscReal currentTime)
Prints a progress bar to the console.
LogLevel get_log_level()
Retrieves the current logging level from the environment variable LOG_LEVEL.
@ LOG_ERROR
Critical errors that may halt the program.
@ LOG_INFO
Informational messages about program execution.
@ LOG_WARNING
Non-critical issues that warrant attention.
@ LOG_DEBUG
Detailed debugging information.
#define PROFILE_FUNCTION_BEGIN
Marks the beginning of a profiled code block (typically a function).
PetscErrorCode ComputeQCriterion(UserCtx *user)
Computes the Q-Criterion, a scalar value identifying vortex cores.
PetscErrorCode ComputeSpecificKE(UserCtx *user, const char *velocity_field, const char *ske_field)
Computes the specific kinetic energy (KE per unit mass) for each particle.
PetscErrorCode NormalizeRelativeField(UserCtx *user, const char *relative_field_name)
Normalizes a relative field by subtracting a reference value.
PetscErrorCode ComputeNodalAverage(UserCtx *user, const char *in_field_name, const char *out_field_name)
Computes node-centered data by averaging 8 surrounding cell-centered values, exactly replicating the ...
PetscErrorCode EulerianDataProcessingPipeline(UserCtx *user, PostProcessParams *pps)
Parses the processing pipeline string from the config and executes the requested kernels in sequence.
PetscErrorCode WriteEulerianFile(UserCtx *user, PostProcessParams *pps, PetscInt ti)
Writes Eulerian field data to a VTK structured grid file (.vts) for visualization.
int main(int argc, char **argv)
PetscErrorCode ParticleDataProcessingPipeline(UserCtx *user, PostProcessParams *pps)
Parses and executes the particle pipeline using a robust two-pass approach.
PetscErrorCode WriteParticleFile(UserCtx *user, PostProcessParams *pps, PetscInt ti)
Writes particle data to a VTP file using the Prepare-Write-Cleanup pattern.
PetscErrorCode SetupPostProcessSwarm(UserCtx *user, PostProcessParams *pps)
Creates a new, dedicated DMSwarm for post-processing tasks.
PetscErrorCode SetupDomainRankInfo(SimCtx *simCtx)
Sets up the full rank communication infrastructure, including neighbor ranks and bounding box exchang...
PetscErrorCode SetupGridAndSolvers(SimCtx *simCtx)
The main orchestrator for setting up all grid-related components.
PetscErrorCode SetupSimulationEnvironment(SimCtx *simCtx)
Verifies and prepares the complete I/O environment for a simulation run.
PetscErrorCode CreateSimulationContext(int argc, char **argv, SimCtx **p_simCtx)
Allocates and populates the master SimulationContext object.
PetscErrorCode SetupBoundaryConditions(SimCtx *simCtx)
(Orchestrator) Sets up all boundary conditions for the simulation.
#define MAX_POINT_DATA_FIELDS
Defines the maximum number of data fields for VTK point data.
char particle_output_prefix[256]
PetscInt num_point_data_fields
SimCtx * simCtx
Back-pointer to the master simulation context.
#define MAX_FILENAME_LENGTH
char output_fields_instantaneous[1024]
char particle_pipeline[1024]
char process_pipeline[1024]
char particle_fields[1024]
PetscBool outputParticles
VTKFieldInfo point_data_fields[20]
@ EXEC_MODE_POSTPROCESSOR
#define MAX_VTK_FIELD_NAME_LENGTH
Maximum length for VTK field names.
Holds all configuration parameters for a post-processing run.
The master context for the entire simulation.
User-defined context containing data specific to a single computational grid level.
Stores all necessary information for a single data array in a VTK file.
PetscErrorCode PrepareOutputEulerianFieldData(UserCtx *user, Vec field_vec, PetscInt num_components, PetscScalar **out_data)
Creates a C array of field data corresponding to a subsampled (legacy-style) grid.
PetscErrorCode PrepareOutputCoordinates(UserCtx *user, PetscScalar **out_coords, PetscInt *out_nx, PetscInt *out_ny, PetscInt *out_nz, PetscInt *out_npoints)
Creates a C array of coordinates corresponding to a subsampled (legacy-style) grid.
PetscErrorCode PrepareOutputParticleData(UserCtx *user, PostProcessParams *pps, VTKMetaData *meta, PetscInt *p_n_total)
Gathers, subsamples, and prepares all particle data for VTK output.