PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
Loading...
Searching...
No Matches
postprocessor.c
Go to the documentation of this file.
1/**
2 * @file postprocessor.c
3 * @brief Phase 2 implementation of the post-processing tool.
4 *
5 * This phase introduces a dedicated configuration system and performs a
6 * single-step data load to verify the I/O and data structures.
7 */
8
9#include "postprocessor.h" // Use our new header
10
11
12#undef __FUNCT__
13#define __FUNCT__ "SetupPostProcessSwarm"
14/**
15 * @brief Creates a new, dedicated DMSwarm for post-processing tasks.
16 *
17 * This function is called once at startup. It creates an empty DMSwarm and
18 * associates it with the same grid DM as the primary swarm and registers all the required fields.
19 * @param user The UserCtx where user->post_swarm will be created.
20 * @param pps The PostProcessParams containing the particle_pipeline string for field registration.
21 * @return PetscErrorCode
22 */
24{
25 PetscErrorCode ierr;
26 PetscFunctionBeginUser;
28 char *pipeline_copy, *step_token, *step_saveptr;
29 PetscBool finalize_needed = PETSC_FALSE;
30
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);
35 // Associate it with the same grid as the solver's swarm
36 if (user->da) {
37 ierr = DMSwarmSetCellDM(user->post_swarm, user->da); CHKERRQ(ierr);
38 LOG_ALLOW(LOCAL,LOG_INFO,"Associated DMSwarm with Cell DM (user->da).\n");
39 } else {
40 // If user->da is essential for your simulation logic with particles, this should be a fatal error.
41 LOG_ALLOW(GLOBAL, LOG_WARNING, "user->da (Cell DM for Swarm) is NULL. Cell-based swarm operations might fail.\n");
42 // SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONGSTATE, "CreateParticleSwarm - user->da (Cell DM) is NULL but required.");
43 }
44
45 LOG_ALLOW(GLOBAL, LOG_INFO, "Created dedicated DMSwarm for post-processing.\n");
46
47 LOG_ALLOW(GLOBAL, LOG_DEBUG, " --- Setting up Post-Processing Pipeline fields-- \n");
48
49 ierr = PetscStrallocpy(pps->particle_pipeline, &pipeline_copy); CHKERRQ(ierr);
50 step_token = strtok_r(pipeline_copy, ";", &step_saveptr);
51 while (step_token) {
52 TrimWhitespace(step_token);
53 if (strlen(step_token) == 0) { step_token = strtok_r(NULL, ";", &step_saveptr); continue; }
54
55 char *keyword = strtok(step_token, ":");
56 char *args_str = strtok(NULL, "");
57 TrimWhitespace(keyword);
58 PetscInt output_field_dimensions = 1; // Default to scalar output fields
59
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; // SKE is scalar
65 if (!output_field) SETERRQ(PETSC_COMM_SELF, 1, "Error (ComputeSpecificKE): Missing output field in 'in>out' syntax.");
66 TrimWhitespace(output_field);
67 // Register the output field
68 ierr = RegisterSwarmField(user->post_swarm, output_field, output_field_dimensions,PETSC_REAL); CHKERRQ(ierr);
69 LOG_ALLOW(GLOBAL, LOG_INFO, "Registered particle field '%s' (ComputeSpecificKE).\n", output_field);
70 finalize_needed = PETSC_TRUE;
71 } else {
72 LOG_ALLOW(GLOBAL, LOG_WARNING, "Warning: Unknown particle transformation keyword '%s'. Skipping.\n", keyword);
73 }
74
75 // Add other 'else if' blocks here for other kernels that create output fields
76
77 step_token = strtok_r(NULL, ";", &step_saveptr);
78 } // while step_token
79
80 ierr = PetscFree(pipeline_copy); CHKERRQ(ierr);
81
82 // --- FINALIZE STEP ---
83 if (finalize_needed) {
84 LOG_ALLOW(GLOBAL, LOG_INFO, "Finalizing new field registrations for particle pipeline.\n");
85 ierr = DMSwarmFinalizeFieldRegister(user->post_swarm); CHKERRQ(ierr);
86 }
87
88 LOG_ALLOW(GLOBAL, LOG_INFO, "Post-Processing DMSwarm setup complete.\n");
89
91 PetscFunctionReturn(0);
92}
93
94
95#undef __FUNCT__
96#define __FUNCT__ "EulerianDataProcessingPipeline"
97/**
98 * @brief Parses the processing pipeline string from the config and executes the requested kernels in sequence.
99 *
100 * This function uses a general-purpose parser to handle a syntax of the form:
101 * "Keyword1:in1>out1; Keyword2:in1,in2>out2; Keyword3:arg1;"
102 *
103 * It tokenizes the pipeline string and dispatches to the appropriate kernel function
104 * from processing_kernels.c with the specified field name arguments.
105 *
106 * @param user The UserCtx containing the data to be transformed.
107 * @param pps The PostProcessParams struct containing the pipeline string.
108 * @return PetscErrorCode
109 */
111{
112 PetscErrorCode ierr;
113 char *pipeline_copy, *step_token, *step_saveptr;
114
115 PetscFunctionBeginUser;
117 LOG_ALLOW(GLOBAL, LOG_INFO, "--- Starting Data Transformation Pipeline ---\n");
118
119 // Do nothing if the pipeline string is empty
120 if (pps->process_pipeline[0] == '\0') {
121 LOG_ALLOW(GLOBAL, LOG_INFO, "Processing pipeline is empty. No transformations will be run.\n");
122 PetscFunctionReturn(0);
123 }
124
125 LOG_ALLOW(GLOBAL, LOG_DEBUG, "Pipeline string: [%s]\n", pps->process_pipeline);
126
127 // Make a writable copy for strtok_r, as it modifies the string
128 ierr = PetscStrallocpy(pps->process_pipeline, &pipeline_copy); CHKERRQ(ierr);
129
130 // --- Outer Loop: Tokenize by Semicolon (;) to get each processing step ---
131 step_token = strtok_r(pipeline_copy, ";", &step_saveptr);
132 while (step_token) {
133 TrimWhitespace(step_token);
134 if (strlen(step_token) == 0) {
135 step_token = strtok_r(NULL, ";", &step_saveptr);
136 continue;
137 }
138
139 char *keyword = strtok(step_token, ":");
140 char *args_str = strtok(NULL, ""); // Get the rest of the string as arguments
141
142 if (!keyword) { // Should not happen with TrimWhitespace, but is a safe check
143 step_token = strtok_r(NULL, ";", &step_saveptr);
144 continue;
145 }
146
147 TrimWhitespace(keyword);
148 if (args_str) TrimWhitespace(args_str);
149
150 LOG_ALLOW(GLOBAL, LOG_INFO, "Executing Transformation: '%s' on args: '%s'\n", keyword, args_str ? args_str : "None");
151
152 // --- DISPATCHER: Route to the correct kernel based on the keyword ---
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);
162 continue;
163 }
164 TrimWhitespace(in_field); TrimWhitespace(out_field);
165 ierr = ComputeNodalAverage(user, in_field, out_field); CHKERRQ(ierr);
166 }
167 else if (strcasecmp(keyword, "ComputeQCriterion") == 0) {
168 ierr = ComputeQCriterion(user); CHKERRQ(ierr);
169 }
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.");
172 ierr = NormalizeRelativeField(user, args_str); CHKERRQ(ierr);
173 }
174 // *** Add new kernels here in the future using 'else if' ***
175 // else if (strcasecmp(keyword, "ComputeVorticity") == 0) { ... }
176 else {
177 LOG_ALLOW(GLOBAL, LOG_WARNING, "Unknown transformation keyword '%s'. Skipping.\n", keyword);
178 }
179
180 step_token = strtok_r(NULL, ";", &step_saveptr);
181 }
182
183 ierr = PetscFree(pipeline_copy); CHKERRQ(ierr);
184 LOG_ALLOW(GLOBAL, LOG_INFO, "--- Data Transformation Pipeline Complete ---\n");
186 PetscFunctionReturn(0);
187}
188
189
190#undef __FUNCT__
191#define __FUNCT__ "WriteEulerianFile"
192/**
193 * @brief Writes Eulerian field data to a VTK structured grid file (.vts) for visualization.
194 *
195 * This function exports instantaneous Eulerian field data (such as pressure, velocity,
196 * Q-criterion, and particle concentration) to a VTK structured grid file for the specified
197 * time step. The output includes subsampled interior grid coordinates and point data fields
198 * as specified in the PostProcessParams configuration.
199 *
200 * The function performs the following steps:
201 * 1. Initializes VTK metadata structure
202 * 2. Prepares output coordinates (subsampled interior grid)
203 * 3. Gathers and prepares requested field data on rank 0
204 * 4. Writes the VTK structured grid file with the naming convention: {prefix}_{ti:05d}.vts
205 *
206 * Only fields listed in pps->output_fields_instantaneous are written. If the field list is
207 * empty, the function returns immediately without creating a file.
208 *
209 * @note Fields containing particle data (e.g., Psi_nodal) are automatically skipped when
210 * no particles are present in the simulation (simCtx->np == 0).
211 *
212 * @param user The UserCtx containing simulation data and field vectors to be output.
213 * @param pps The PostProcessParams struct containing output configuration (field list, prefix).
214 * @param ti The time index/step number used for file naming.
215 * @return PetscErrorCode indicating success or failure.
216 */
217PetscErrorCode WriteEulerianFile(UserCtx* user, PostProcessParams* pps, PetscInt ti)
218{
219 PetscErrorCode ierr;
220 VTKMetaData meta;
221 char filename[MAX_FILENAME_LENGTH];
222
223 PetscFunctionBeginUser;
225
226 if (pps->output_fields_instantaneous[0] == '\0') {
227 LOG_ALLOW(GLOBAL, LOG_DEBUG, "No instantaneous fields requested for output at ti=%" PetscInt_FMT ". Skipping.\n", ti);
228 PetscFunctionReturn(0);
229 }
230
231 LOG_ALLOW(GLOBAL, LOG_INFO, "--- Starting VTK File Writing for ti = %" PetscInt_FMT " ---\n", ti);
232
233 /* 1) Metadata init */
234 ierr = PetscMemzero(&meta, sizeof(VTKMetaData)); CHKERRQ(ierr);
236
237 /* 2) Coordinates (subsampled interior) */
238 ierr = PrepareOutputCoordinates(user, &meta.coords, &meta.mx, &meta.my, &meta.mz, &meta.npoints); 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);
241
242 /* 3) Fields (rank 0) */
243 if (user->simCtx->rank == 0) {
244 char *fields_copy, *field_name;
245 ierr = PetscStrallocpy(pps->output_fields_instantaneous, &fields_copy); CHKERRQ(ierr);
246
247 field_name = strtok(fields_copy, ",");
248 while (field_name) {
249 TrimWhitespace(field_name);
250 if (!*field_name) { field_name = strtok(NULL, ","); continue; }
251
252 LOG_ALLOW(LOCAL, LOG_DEBUG, "Preparing field '%s' for output.\n", field_name);
253
254 Vec field_vec = NULL;
255 PetscInt num_components = 0;
256
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")){
264 if(user->simCtx->np==0){
265 LOG_ALLOW(LOCAL, LOG_WARNING, "Field 'Psi_nodal' requested but no particles are present. Skipping.\n");
266 field_name = strtok(NULL, ",");
267 continue;
268 }
269 field_vec = user->Psi_nodal; num_components = 1;
270 } else {
271 LOG_ALLOW(LOCAL, LOG_WARNING, "Field '%s' not recognized. Skipping.\n", field_name);
272 field_name = strtok(NULL, ",");
273 continue;
274 }
275
277 LOG_ALLOW(LOCAL, LOG_WARNING, "MAX_POINT_DATA_FIELDS reached. Cannot add '%s'.\n", field_name);
278 field_name = strtok(NULL, ",");
279 continue;
280 }
281
282 // --- Add field to metadata ---
283 VTKFieldInfo* current_field = &meta.point_data_fields[meta.num_point_data_fields];
284 strncpy(current_field->name, field_name, MAX_VTK_FIELD_NAME_LENGTH-1);
285 current_field->name[MAX_VTK_FIELD_NAME_LENGTH-1] = '\0';
286 current_field->num_components = num_components;
287
288 /* Build interior AoS from NATURAL gathered Vec */
289 ierr = PrepareOutputEulerianFieldData(user, field_vec, num_components, &current_field->data); CHKERRQ(ierr);
290
291 /*
292 // *** DEBUG: Dump Ucat_nodal details and add scalar companions Ux, Uy, Uz for easier visualization ***
293
294 // If this is Ucat_nodal, dump a few tuples and add scalar companions
295 if (!strcasecmp(field_name, "Ucat_nodal")) {
296 const PetscInt npts = meta.npoints;
297 const PetscScalar *a = (const PetscScalar*)current_field->data;
298
299 LOG_ALLOW(GLOBAL, LOG_INFO, "DBG Ucat_nodal: ptr=%p npoints=%" PetscInt_FMT " num_components=%" PetscInt_FMT "\n",
300 (void*)a, npts, current_field->num_components);
301
302 if (a && current_field->num_components == 3 && npts > 0) {
303 const PetscInt nshow = (npts < 5) ? npts : 5;
304 LOG_ALLOW(GLOBAL, LOG_INFO, "DBG Ucat_nodal: showing first %d of %" PetscInt_FMT " tuples (AoS x,y,z):\n",
305 (int)nshow, npts);
306 for (PetscInt t = 0; t < nshow; ++t) {
307 LOG_ALLOW(GLOBAL, LOG_INFO, " Ucat_nodal[%3" PetscInt_FMT "] = (%g, %g, %g)\n",
308 t, (double)a[3*t+0], (double)a[3*t+1], (double)a[3*t+2]);
309 }
310 if (npts > 10) {
311 PetscInt mid = npts / 2;
312 LOG_ALLOW(GLOBAL, LOG_INFO, " Ucat_nodal[mid=%" PetscInt_FMT "] = (%g, %g, %g)\n",
313 mid, (double)a[3*mid+0], (double)a[3*mid+1], (double)a[3*mid+2]);
314 }
315
316 // Add scalar companions from the AoS we just created
317 PetscScalar *Ux=NULL,*Uy=NULL,*Uz=NULL;
318 ierr = PetscMalloc1(meta.npoints, &Ux); CHKERRQ(ierr);
319 ierr = PetscMalloc1(meta.npoints, &Uy); CHKERRQ(ierr);
320 ierr = PetscMalloc1(meta.npoints, &Uz); CHKERRQ(ierr);
321 for (PetscInt i = 0; i < meta.npoints; ++i) {
322 Ux[i] = a[3*i+0]; Uy[i] = a[3*i+1]; Uz[i] = a[3*i+2];
323 }
324
325 if (meta.num_point_data_fields + 3 <= MAX_POINT_DATA_FIELDS) {
326 VTKFieldInfo *fx = &meta.point_data_fields[++meta.num_point_data_fields];
327 strncpy(fx->name, "Ux_debug", MAX_VTK_FIELD_NAME_LENGTH-1);
328 fx->name[MAX_VTK_FIELD_NAME_LENGTH-1] = '\0';
329 fx->num_components = 1; fx->data = Ux;
330
331 VTKFieldInfo *fy = &meta.point_data_fields[++meta.num_point_data_fields];
332 strncpy(fy->name, "Uy_debug", MAX_VTK_FIELD_NAME_LENGTH-1);
333 fy->name[MAX_VTK_FIELD_NAME_LENGTH-1] = '\0';
334 fy->num_components = 1; fy->data = Uy;
335
336 VTKFieldInfo *fz = &meta.point_data_fields[++meta.num_point_data_fields];
337 strncpy(fz->name, "Uz_debug", MAX_VTK_FIELD_NAME_LENGTH-1);
338 fz->name[MAX_VTK_FIELD_NAME_LENGTH-1] = '\0';
339 fz->num_components = 1; fz->data = Uz;
340
341 LOG_ALLOW(GLOBAL, LOG_INFO, "DBG: Added scalar companions Ux_debug, Uy_debug, Uz_debug.\n");
342 } else {
343 LOG_ALLOW(GLOBAL, LOG_WARNING, "DBG: Not enough slots to add Ux/Uy/Uz debug fields.\n");
344 PetscFree(Ux); PetscFree(Uy); PetscFree(Uz);
345 }
346
347 // Mid-plane CSV + AoS vs NATURAL compare (component X)
348
349 // Gather NATURAL again (small cost, but isolated and clear)
350 PetscInt Ng = 0;
351 double *nat_d = NULL;
352 DM dmU = NULL;
353 DMDALocalInfo infU;
354 ierr = VecGetDM(field_vec, &dmU); CHKERRQ(ierr);
355 ierr = DMDAGetLocalInfo(dmU, &infU); CHKERRQ(ierr);
356
357 const PetscInt M=infU.mx, N=infU.my, P=infU.mz;
358 const PetscInt mx = meta.mx, my = meta.my, mz = meta.mz;
359 const PetscInt iInnerMid = mx/2; // interior index [0..mx-1]
360 const PetscInt iGlob = iInnerMid;
361
362 ierr = VecToArrayOnRank0(field_vec, &Ng, &nat_d); CHKERRQ(ierr);
363
364 if (nat_d) {
365 const PetscScalar *nar = (const PetscScalar*)nat_d;
366 const char *base = pps->output_prefix;
367 char fn[512], fnc[512];
368 snprintf(fn, sizeof(fn), "%s_%05" PetscInt_FMT "_iMid.csv", base, ti);
369 snprintf(fnc, sizeof(fnc), "%s_%05" PetscInt_FMT "_iMid_compare.csv", base, ti);
370
371 FILE *fp = fopen(fn, "w");
372 FILE *fpc = fopen(fnc, "w");
373 if (fp) fprintf(fp, "jInner,kInner,Ux,Uy,Uz\n");
374 if (fpc) fprintf(fpc, "jInner,kInner,Ux_AoS,Ux_NAT,abs_diff\n");
375
376 double maxAbsDiff = 0.0, sumAbs = 0.0;
377 PetscInt count = 0;
378
379 for (PetscInt kInner = 0; kInner < mz; ++kInner) {
380 const PetscInt k = kInner;
381 for (PetscInt jInner = 0; jInner < my; ++jInner) {
382 const PetscInt j = jInner;
383
384 // AoS tuple index
385 const PetscInt t = iInnerMid + mx * (jInner + my * kInner);
386 const PetscScalar ux = a[3*t+0], uy = a[3*t+1], uz = a[3*t+2];
387
388 if (fp) fprintf(fp, "%d,%d,%.15e,%.15e,%.15e\n",
389 (int)jInner,(int)kInner,(double)ux,(double)uy,(double)uz);
390
391 // NATURAL base for (iGlob,j,k)
392 const PetscInt baseNat = 3 * (((k)*N + j)*M + iGlob);
393 const PetscScalar uxN = nar[baseNat + 0];
394
395 const double diff = fabs((double)ux - (double)uxN);
396 if (diff > maxAbsDiff) maxAbsDiff = diff;
397 sumAbs += diff; ++count;
398
399 if (fpc) fprintf(fpc, "%d,%d,%.15e,%.15e,%.15e\n",
400 (int)jInner,(int)kInner,(double)ux,(double)uxN,diff);
401 } // for jInner
402 } // for kInner
403 if (fp) fclose(fp);
404 if (fpc) fclose(fpc);
405
406 if (count > 0) {
407 const double meanAbs = sumAbs / (double)count;
408 LOG_ALLOW(GLOBAL, LOG_INFO,
409 "PETSc-Vec vs AoS (i-mid, Ux): max|Δ|=%.6e, mean|Δ|=%.6e -> CSV: %s\n",
410 maxAbsDiff, meanAbs, fnc);
411 LOG_ALLOW(GLOBAL, LOG_INFO, "Wrote i-mid plane CSV: %s\n", fn);
412 } // if count>0
413
414 ierr = PetscFree(nat_d); CHKERRQ(ierr);
415
416
417 } // if nat_d
418
419 } // if a && num_components==3 && npts>0
420 } // if Ucat_nodal
421 // --- END DEBUG BLOCK (Ucat_nodal) ---
422 */
423
424 meta.num_point_data_fields++; /* count the main field we just added */
425 field_name = strtok(NULL, ",");
426 }
427
428 ierr = PetscFree(fields_copy); CHKERRQ(ierr);
429
430 // --- DEBUG: Add sanity fields i_idx, j_idx, k_idx and x_pos, y_pos, z_pos ---
431 // These are the logical indices and physical coordinates of each point in the subsampled grid.
432 // They can be used to verify the grid structure and orientation in visualization tools.
433 // They are added as scalar fields with names "i_idx", "j_idx", "k_idx" and "x_pos", "y_pos", "z_pos".
434 // Note: these are only added if there is room in the MAX_POINT_DATA_FIELDS limit.
435 // They are allocated and owned here, and will be freed below.
436 // They are in the same linearization order as meta.coords (AoS x,y,z by point).
437 /*
438 // Append sanity fields i/j/k indices and coordinates
439
440 // Build i/j/k and x/y/z (length = npoints); these match the same linearization as coords
441 const PetscInt n = meta.npoints;
442 PetscScalar *i_idx=NULL,*j_idx=NULL,*k_idx=NULL,*x_pos=NULL,*y_pos=NULL,*z_pos=NULL;
443
444 ierr = PetscMalloc1(n, &i_idx); CHKERRQ(ierr);
445 ierr = PetscMalloc1(n, &j_idx); CHKERRQ(ierr);
446 ierr = PetscMalloc1(n, &k_idx); CHKERRQ(ierr);
447 ierr = PetscMalloc1(n, &x_pos); CHKERRQ(ierr);
448 ierr = PetscMalloc1(n, &y_pos); CHKERRQ(ierr);
449 ierr = PetscMalloc1(n, &z_pos); CHKERRQ(ierr);
450
451 // coords is length 3*n, AoS: (x,y,z) by point
452 const PetscScalar *c = (const PetscScalar*)meta.coords;
453 for (PetscInt k = 0; k < meta.mz; ++k) {
454 for (PetscInt j = 0; j < meta.my; ++j) {
455 for (PetscInt i = 0; i < meta.mx; ++i) {
456 const PetscInt t = i + meta.mx * (j + meta.my * k);
457 i_idx[t] = (PetscScalar)i;
458 j_idx[t] = (PetscScalar)j;
459 k_idx[t] = (PetscScalar)k;
460 x_pos[t] = c[3*t+0];
461 y_pos[t] = c[3*t+1];
462 z_pos[t] = c[3*t+2];
463 }
464 }
465 }
466
467 const char *nf[6] = {"i_idx","j_idx","k_idx","x_pos","y_pos","z_pos"};
468 PetscScalar *arrs[6] = {i_idx,j_idx,k_idx,x_pos,y_pos,z_pos};
469 for (int s=0; s<6; ++s) {
470 if (meta.num_point_data_fields < MAX_POINT_DATA_FIELDS) {
471 VTKFieldInfo *f = &meta.point_data_fields[meta.num_point_data_fields++];
472 strncpy(f->name, nf[s], MAX_VTK_FIELD_NAME_LENGTH-1);
473 f->name[MAX_VTK_FIELD_NAME_LENGTH-1] = '\0';
474 f->num_components = 1;
475 f->data = arrs[s];
476 } else {
477 LOG_ALLOW(GLOBAL, LOG_WARNING, "Sanity field '%s' dropped: MAX_POINT_DATA_FIELDS reached.\n", nf[s]);
478 PetscFree(arrs[s]);
479 }
480 }
481 LOG_ALLOW(GLOBAL, LOG_INFO, "DBG: Added sanity fields i_idx/j_idx/k_idx and x_pos/y_pos/z_pos.\n");
482 */
483 // --- END DEBUG BLOCK (sanity fields) ---
484
485 /* Field summary */
486 LOG_ALLOW(GLOBAL, LOG_INFO, "PointData fields to write: %d\n", (int)meta.num_point_data_fields);
487 for (PetscInt ii=0; ii<meta.num_point_data_fields; ++ii) {
488 LOG_ALLOW(GLOBAL, LOG_INFO, " # %2" PetscInt_FMT " Field Name = %s Components = %d\n",
489 ii, meta.point_data_fields[ii].name, (int)meta.point_data_fields[ii].num_components);
490 }
491 } // if rank 0
492
493 /* 4) Write the VTS */
494 snprintf(filename, sizeof(filename), "%s_%05" PetscInt_FMT ".vts", pps->output_prefix, ti);
495 ierr = CreateVTKFileFromMetadata(filename, &meta, PETSC_COMM_WORLD); CHKERRQ(ierr);
496
497 LOG_ALLOW(GLOBAL, LOG_INFO, "--- Eulerian File Writing for ti = %" PetscInt_FMT " Complete ---\n", ti);
499 PetscFunctionReturn(0);
500}
501
502
503#undef __FUNCT__
504#define __FUNCT__ "ParticleDataProcessingPipeline"
505/**
506 * @brief Parses and executes the particle pipeline using a robust two-pass approach.
507 *
508 * This function ensures correctness and efficiency by separating field registration
509 * from kernel execution.
510 *
511 * PASS 1 (Registration): The pipeline string is parsed to identify all new fields
512 * that will be created. These fields are registered with the DMSwarm.
513 *
514 * Finalize: After Pass 1, DMSwarmFinalizeFieldRegister is called exactly once if
515 * any new fields were added, preparing the swarm's memory layout.
516 *
517 * PASS 2 (Execution): The pipeline string is parsed again, and this time the
518 * actual compute kernels are executed, filling the now-valid fields.
519 *
520 * @param user The UserCtx containing the DMSwarm.
521 * @param pps The PostProcessParams struct containing the particle_pipeline string.
522 * @return PetscErrorCode
523 */
525{
526 PetscErrorCode ierr;
527 char *pipeline_copy, *step_token, *step_saveptr;
528
529 PetscFunctionBeginUser;
530
532
533 // --- Timestep Setup: Synchronize post_swarm size ---
534 PetscInt n_global_source;
535 ierr = DMSwarmGetSize(user->swarm, &n_global_source); CHKERRQ(ierr);
536
537 // Resize post_swarm to match source swarm
538 ierr = ResizeSwarmGlobally(user->post_swarm, n_global_source); CHKERRQ(ierr);
539
540 if (pps->particle_pipeline[0] == '\0') {
541 PetscFunctionReturn(0);
542 }
543
544 LOG_ALLOW(GLOBAL, LOG_INFO, "--- Starting Particle Data Transformation Pipeline ---\n");
545 LOG_ALLOW(GLOBAL, LOG_DEBUG, "Particle Pipeline string: [%s]\n", pps->particle_pipeline);
546
547 LOG_ALLOW(GLOBAL, LOG_DEBUG, "Executing compute kernels...\n");
548 ierr = PetscStrallocpy(pps->particle_pipeline, &pipeline_copy); CHKERRQ(ierr);
549 step_token = strtok_r(pipeline_copy, ";", &step_saveptr);
550 while (step_token) {
551 TrimWhitespace(step_token);
552 if (strlen(step_token) == 0) { step_token = strtok_r(NULL, ";", &step_saveptr); continue; }
553
554 char *keyword = strtok(step_token, ":");
555 char *args_str = strtok(NULL, "");
556 TrimWhitespace(keyword);
557 if (args_str) TrimWhitespace(args_str);
558
559 LOG_ALLOW(GLOBAL, LOG_INFO, "Executing Particle Transformation: '%s' on args: '%s'\n", keyword, args_str ? args_str : "None");
560
561 if (strcasecmp(keyword, "ComputeSpecificKE") == 0) {
562 char *velocity_field = strtok(args_str, ">");
563 char *ske_field = strtok(NULL, ">");
564 TrimWhitespace(velocity_field); TrimWhitespace(ske_field);
565
566 ierr = ComputeSpecificKE(user, velocity_field, ske_field); CHKERRQ(ierr);
567 }
568 else {
569 LOG_ALLOW(GLOBAL, LOG_WARNING, "Unknown particle transformation keyword '%s'. Skipping.\n", keyword);
570 }
571
572 step_token = strtok_r(NULL, ";", &step_saveptr);
573 }
574 ierr = PetscFree(pipeline_copy); CHKERRQ(ierr);
575
576 LOG_ALLOW(GLOBAL, LOG_INFO, "--- Particle Data Transformation Pipeline Complete ---\n");
577
579 PetscFunctionReturn(0);
580}
581
582#undef __FUNCT__
583#define __FUNCT__ "WriteParticleFile"
584/**
585 * @brief Writes particle data to a VTP file using the Prepare-Write-Cleanup pattern.
586 */
587PetscErrorCode WriteParticleFile(UserCtx* user, PostProcessParams* pps, PetscInt ti)
588{
589 PetscErrorCode ierr;
590 VTKMetaData part_meta;
591 char filename[MAX_FILENAME_LENGTH];
592 PetscInt n_total_particles_before_subsample;
593
594 PetscFunctionBeginUser;
596
597 // These checks can be done on all ranks
598 if (!pps->outputParticles || pps->particle_fields[0] == '\0') {
599 PetscFunctionReturn(0);
600 }
601 PetscInt n_global;
602 ierr = DMSwarmGetSize(user->swarm, &n_global); CHKERRQ(ierr);
603 if (n_global == 0) {
604 LOG_ALLOW(GLOBAL, LOG_DEBUG, "Swarm is empty for ti=%" PetscInt_FMT ". Skipping particle file write.\n", ti);
605 PetscFunctionReturn(0);
606 }
607
608 ierr = PetscMemzero(&part_meta, sizeof(VTKMetaData)); CHKERRQ(ierr);
609
610 // --- 1. PREPARE (Collective Call) ---
611 ierr = PrepareOutputParticleData(user, pps, &part_meta, &n_total_particles_before_subsample); CHKERRQ(ierr);
612
613 // --- 2. WRITE and CLEANUP (Rank 0 only) ---
614 if (user->simCtx->rank == 0) {
615 if (part_meta.npoints > 0) {
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);
618
619 /* Field summary */
620 LOG_ALLOW(GLOBAL, LOG_INFO, "Particle Data fields to write: %d\n", (int)part_meta.num_point_data_fields);
621 for (PetscInt ii=0; ii<part_meta.num_point_data_fields; ++ii) {
622 LOG_ALLOW(GLOBAL, LOG_INFO, " # %2" PetscInt_FMT " Field Name = %s Components = %d\n",
623 ii, part_meta.point_data_fields[ii].name, (int)part_meta.point_data_fields[ii].num_components);
624 }
625
626 snprintf(filename, sizeof(filename), "%s_%05" PetscInt_FMT ".vtp", pps->particle_output_prefix, ti);
627 ierr = CreateVTKFileFromMetadata(filename, &part_meta, PETSC_COMM_WORLD); CHKERRQ(ierr);
628
629 } else {
630 LOG_ALLOW(GLOBAL, LOG_DEBUG, "No particles to write at ti=%" PetscInt_FMT " after subsampling. Skipping.\n", ti);
631 }
632
633 }
634
635 LOG_ALLOW(GLOBAL, LOG_INFO, "--- Particle File Writing for ti = %" PetscInt_FMT " Complete ---\n", ti);
637 PetscFunctionReturn(0);
638}
639
640#undef __FUNCT__
641#define __FUNCT__ "main"
642int main(int argc, char **argv)
643{
644 PetscErrorCode ierr;
645 SimCtx *simCtx = NULL;
646
647 // === I. INITIALIZE PETSC & MPI ===========================================
648 ierr = PetscInitialize(&argc, &argv, (char *)0, "Unified Post-Processing Tool"); CHKERRQ(ierr);
649
650 // === II. CONFIGURE SIMULATION & POST-PROCESSING CONTEXTS =================
651 ierr = CreateSimulationContext(argc, argv, &simCtx); CHKERRQ(ierr);
652 // === IIB. SET EXECUTION MODE (SOLVER vs POST-PROCESSOR) =====
654 // == IIC. CONFIGURE SIMULATION ENVIRONMENT & DIRECTORIES =====
655 ierr = SetupSimulationEnvironment(simCtx); CHKERRQ(ierr);
656 // === III. SETUP GRID & DATA STRUCTURES ===================================
657 ierr = SetupGridAndSolvers(simCtx); CHKERRQ(ierr);
658 // === IV. SETUP DOMAIN DECOMPOSITION INFORMATION =========================
659 ierr = SetupDomainRankInfo(simCtx); CHKERRQ(ierr);
660 // === V. SETUP BOUNDARY CONDITIONS ====================================
661 ierr = SetupBoundaryConditions(simCtx); CHKERRQ(ierr);
662 // === VI. SETUP USER CONTEXT & DATA STRUCTURES ============================
663 // Get the finest-level user context, as this is where we'll load data
664 UserCtx *user = simCtx->usermg.mgctx[simCtx->usermg.mglevels-1].user;
665 PostProcessParams *pps = simCtx->pps;
666
667 // === VI. PARTICLE INITIALIZATION (if needed) ============================
668 if(pps->outputParticles) {
669 if(simCtx->np > 0){
670 ierr = InitializeParticleSwarm(simCtx); CHKERRQ(ierr);
671 // Create a post-processing specific DMSwarm
672 ierr = SetupPostProcessSwarm(user,pps); CHKERRQ(ierr);
673 }else{
674 SETERRQ(PETSC_COMM_SELF,1,"Particle output requested but np=0. Please set np>0 to enable particle output.");
675 }
676 }
677
678 LOG_ALLOW(GLOBAL, LOG_INFO, "=============================================================\n");
679
680
681 // === VII. MAIN POST-PROCESSING LOOP ======================================
682 for (PetscInt ti = pps->startTime; ti <= pps->endTime; ti += pps->timeStep) {
683 LOG_ALLOW(GLOBAL, LOG_INFO, "--- Processing Time Step %" PetscInt_FMT " ---\n", ti);
684
685 // 1. Load Data (UpdateLocalGhosts is called inside the kernels)
686 ierr = ReadSimulationFields(user, ti); CHKERRQ(ierr);
687
688 // 2. Transform Data
689 ierr = EulerianDataProcessingPipeline(user, pps); CHKERRQ(ierr);
690
691 // 3. Write Output
692 ierr = WriteEulerianFile(user, pps, ti); CHKERRQ(ierr);
693
694 if(pps->outputParticles) {
695 // 1. Resize swarm based on particle count in this timestep's file
696 ierr = PreCheckAndResizeSwarm(user, ti, pps->particleExt); CHKERRQ(ierr);
697
698 // 2. Load particle data into the correctly sized swarm
699 ierr = ReadAllSwarmFields(user, ti); CHKERRQ(ierr);
700
701 // 3. Transform particle data
702 ierr = ParticleDataProcessingPipeline(user, pps); CHKERRQ(ierr);
703
704 // 4. Write particle output
705 ierr = WriteParticleFile(user, pps, ti); CHKERRQ(ierr);
706 }
707
708 if(simCtx->rank == 0){
709 PetscInt StepsToRun = pps->endTime - pps->startTime;
710 PetscReal currentTime = (PetscReal)ti*simCtx->dt;
711 PrintProgressBar(ti-1,pps->startTime,StepsToRun,currentTime);
712 if(get_log_level()>LOG_ERROR)PetscPrintf(PETSC_COMM_SELF,"\n");
713 }
714 }
715
716 // After the loop, print the 100% complete bar on rank 0 and add a newline
717 // to ensure subsequent terminal output starts on a fresh line.
718 if (simCtx->rank == 0) {
719 PetscInt endTime = pps->endTime-1; // needs to be verified.
720 PetscInt StepsToRun = pps->endTime - pps->startTime;
721 PetscReal endTimeValue = (PetscReal)pps->endTime*simCtx->dt;
722 PrintProgressBar(endTime, pps->startTime, StepsToRun, endTimeValue);
723 PetscPrintf(PETSC_COMM_SELF, "\n");
724 fflush(stdout);
725 }
726
727 LOG_ALLOW(GLOBAL, LOG_INFO, "=============================================================\n");
728 LOG_ALLOW(GLOBAL, LOG_INFO, "Post-processing finished successfully.\n");
729
730
731 // === VIII. FINALIZE =========================================================
732 // ierr = FinalizeSimulation(simCtx); CHKERRQ(ierr);
733 ierr = ProfilingFinalize(simCtx); CHKERRQ(ierr);
734 ierr = PetscFinalize();
735 return ierr;
736}
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.
Definition io.c:1024
void TrimWhitespace(char *str)
Helper function to trim leading/trailing whitespace from a string.
Definition io.c:36
PetscErrorCode ReadAllSwarmFields(UserCtx *user, PetscInt ti)
Reads multiple fields (positions, velocity, CellID, and weight) into a DMSwarm.
Definition io.c:1327
PetscInt CreateVTKFileFromMetadata(const char *filename, const VTKMetaData *meta, MPI_Comm comm)
Definition vtk_io.c:128
#define LOCAL
Logging scope definitions for controlling message output.
Definition logging.h:46
#define GLOBAL
Scope for global logging across all processes.
Definition logging.h:47
#define LOG_ALLOW(scope, level, fmt,...)
Logging macro that checks both the log level and whether the calling function is in the allowed-funct...
Definition logging.h:201
#define PROFILE_FUNCTION_END
Marks the end of a profiled code block.
Definition logging.h:740
PetscErrorCode ProfilingFinalize(SimCtx *simCtx)
the profiling excercise and build a profiling summary which is then printed to a log file.
Definition logging.c:1106
#define LOG(scope, level, fmt,...)
Logging macro for PETSc-based applications with scope control.
Definition logging.h:85
void PrintProgressBar(PetscInt step, PetscInt startStep, PetscInt totalSteps, PetscReal currentTime)
Prints a progress bar to the console.
Definition logging.c:1207
LogLevel get_log_level()
Retrieves the current logging level from the environment variable LOG_LEVEL.
Definition logging.c:39
@ LOG_ERROR
Critical errors that may halt the program.
Definition logging.h:29
@ LOG_INFO
Informational messages about program execution.
Definition logging.h:32
@ LOG_WARNING
Non-critical issues that warrant attention.
Definition logging.h:30
@ LOG_DEBUG
Detailed debugging information.
Definition logging.h:33
#define PROFILE_FUNCTION_BEGIN
Marks the beginning of a profiled code block (typically a function).
Definition logging.h:731
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...
Definition setup.c:1981
PetscErrorCode SetupGridAndSolvers(SimCtx *simCtx)
The main orchestrator for setting up all grid-related components.
Definition setup.c:857
PetscErrorCode SetupSimulationEnvironment(SimCtx *simCtx)
Verifies and prepares the complete I/O environment for a simulation run.
Definition setup.c:558
PetscErrorCode CreateSimulationContext(int argc, char **argv, SimCtx **p_simCtx)
Allocates and populates the master SimulationContext object.
Definition setup.c:44
PetscErrorCode SetupBoundaryConditions(SimCtx *simCtx)
(Orchestrator) Sets up all boundary conditions for the simulation.
Definition setup.c:1354
Vec P_nodal
Definition variables.h:734
UserCtx * user
Definition variables.h:427
#define MAX_POINT_DATA_FIELDS
Defines the maximum number of data fields for VTK point data.
Definition variables.h:446
PetscInt npoints
Definition variables.h:499
PetscInt num_components
Definition variables.h:486
char particle_output_prefix[256]
Definition variables.h:469
PetscMPIInt rank
Definition variables.h:541
PetscInt num_point_data_fields
Definition variables.h:502
SimCtx * simCtx
Back-pointer to the master simulation context.
Definition variables.h:664
char output_prefix[256]
Definition variables.h:466
UserMG usermg
Definition variables.h:631
DM post_swarm
Definition variables.h:733
PetscReal dt
Definition variables.h:552
Vec Ucat_nodal
Definition variables.h:735
PetscInt timeStep
Definition variables.h:459
PetscInt np
Definition variables.h:616
#define MAX_FILENAME_LENGTH
Definition variables.h:444
Vec Qcrit
Definition variables.h:736
VTKFileType fileType
Definition variables.h:497
char output_fields_instantaneous[1024]
Definition variables.h:464
PetscScalar * data
Definition variables.h:487
char particle_pipeline[1024]
Definition variables.h:467
char name[64]
Definition variables.h:485
PetscScalar * coords
Definition variables.h:500
PetscInt mglevels
Definition variables.h:434
char process_pipeline[1024]
Definition variables.h:463
char particle_fields[1024]
Definition variables.h:468
PetscBool outputParticles
Definition variables.h:460
VTKFieldInfo point_data_fields[20]
Definition variables.h:501
Vec Psi_nodal
Definition variables.h:737
PetscInt mz
Definition variables.h:498
PostProcessParams * pps
Definition variables.h:648
@ EXEC_MODE_POSTPROCESSOR
Definition variables.h:512
char particleExt[8]
Definition variables.h:474
MGCtx * mgctx
Definition variables.h:437
ExecutionMode exec_mode
Definition variables.h:556
PetscInt startTime
Definition variables.h:457
PetscInt my
Definition variables.h:498
PetscInt mx
Definition variables.h:498
#define MAX_VTK_FIELD_NAME_LENGTH
Maximum length for VTK field names.
Definition variables.h:447
@ VTK_STRUCTURED
Definition variables.h:492
Holds all configuration parameters for a post-processing run.
Definition variables.h:452
The master context for the entire simulation.
Definition variables.h:538
User-defined context containing data specific to a single computational grid level.
Definition variables.h:661
Stores all necessary information for a single data array in a VTK file.
Definition variables.h:484
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.
Definition vtk_io.c:287
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.
Definition vtk_io.c:223
PetscErrorCode PrepareOutputParticleData(UserCtx *user, PostProcessParams *pps, VTKMetaData *meta, PetscInt *p_n_total)
Gathers, subsamples, and prepares all particle data for VTK output.
Definition vtk_io.c:421