PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
Loading...
Searching...
No Matches
logging.c
Go to the documentation of this file.
1// logging.c
2#include "logging.h"
3
4/* Maximum temporary buffer size for converting numbers to strings */
5#define TMP_BUF_SIZE 128
6
7// --------------------- Static Variable for Log Level ---------------------
8
9/**
10 * @brief Static variable to cache the current logging level.
11 *
12 * Initialized to -1 to indicate that the log level has not been set yet.
13 */
15
16// --------------------- Static Variables for Allow-List -------------------
17
18/**
19 * @brief Global/static array of function names allowed to log.
20 */
21static char** gAllowedFunctions = NULL;
22
23/**
24 * @brief Number of entries in the gAllowedFunctions array.
25 */
26static int gNumAllowed = 0;
27
28// --------------------- Function Implementations ---------------------
29
30/**
31 * @brief Retrieves the current logging level from the environment variable `LOG_LEVEL`.
32 *
33 * The function checks the `LOG_LEVEL` environment variable and sets the logging level accordingly.
34 * Supported levels are "DEBUG", "INFO", "WARNING", and defaults to "ERROR" if not set or unrecognized.
35 * The log level is cached after the first call to avoid repeated environment variable checks.
36 *
37 * @return LogLevel The current logging level.
38 */
40 if (current_log_level == -1) { // Log level not set yet
41 const char *env = getenv("LOG_LEVEL");
42 if (!env) {
43 current_log_level = LOG_ERROR; // Default level
44 }
45 else if (strcmp(env, "DEBUG") == 0) {
47 }
48 else if (strcmp(env, "INFO") == 0) {
50 }
51 else if (strcmp(env, "WARNING") == 0) {
53 }
54 else if (strcmp(env, "PROFILE") == 0) { // <-- New profile level
56 }
57 else if (strcmp(env, "VERBOSE") == 0) {
59 }
60 else if (strcmp(env, "TRACE") == 0) {
62 }
63 else {
64 current_log_level = LOG_ERROR; // Default if unrecognized
65 }
66 }
67 return current_log_level;
68}
69
70/**
71 * @brief Prints the current logging level to the console.
72 *
73 * This function retrieves the log level using `get_log_level()` and prints
74 * the corresponding log level name. It helps verify the logging configuration
75 * at runtime.
76 *
77 * @note The log level is determined from the `LOG_LEVEL` environment variable.
78 * If `LOG_LEVEL` is not set, it defaults to `LOG_INFO`.
79 *
80 * @see get_log_level()
81 */
82PetscErrorCode print_log_level(void)
83{
84 PetscMPIInt rank;
85 PetscErrorCode ierr;
86 int level;
87 const char *level_name;
88
89 PetscFunctionBeginUser;
90 /* get MPI rank */
91 ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank); CHKERRMPI(ierr);
92
93 /* decide level name */
94 level = get_log_level();
95 level_name = (level == LOG_ERROR) ? "ERROR" :
96 (level == LOG_WARNING) ? "WARNING" :
97 (level == LOG_INFO) ? "INFO" :
98 (level == LOG_DEBUG) ? "DEBUG" :
99 (level == LOG_VERBOSE) ? "VERBOSE" :
100 (level == LOG_TRACE) ? "TRACE" :
101 (level == LOG_PROFILE) ? "PROFILE" : "UNKNOWN";
102
103 /* print it out */
104 ierr = PetscPrintf(PETSC_COMM_SELF,
105 "Current log level: %s (%d) | rank: %d\n",
106 level_name, level, (int)rank);
107 CHKERRMPI(ierr);
108
109 PetscFunctionReturn(PETSC_SUCCESS);
110}
111
112/**
113 * @brief Sets the global list of function names that are allowed to log.
114 *
115 * @param functionList An array of function name strings (e.g., {"foo", "bar"}).
116 * @param count The number of entries in the array.
117 *
118 * The existing allow-list is cleared and replaced by the new one.
119 * If you pass an empty list (count = 0), then no function will be allowed
120 * unless you change it later.
121 */
122void set_allowed_functions(const char** functionList, int count)
123{
124 // 1. Free any existing entries
125 if (gAllowedFunctions) {
126 for (int i = 0; i < gNumAllowed; ++i) {
127 free(gAllowedFunctions[i]); // each was strdup'ed
128 }
129 free(gAllowedFunctions);
130 gAllowedFunctions = NULL;
131 gNumAllowed = 0;
132 }
133
134 // 2. Allocate new array
135 if (count > 0) {
136 gAllowedFunctions = (char**)malloc(sizeof(char*) * count);
137 }
138
139 // 3. Copy the new entries
140 for (int i = 0; i < count; ++i) {
141 // strdup is a POSIX function. If not available, implement your own string copy.
142 gAllowedFunctions[i] = strdup(functionList[i]);
143 }
144 gNumAllowed = count;
145}
146
147/**
148 * @brief Checks if the given function name is in the allow-list.
149 *
150 * @param functionName The name of the function to check.
151 * @return PETSC_TRUE if functionName is allowed, otherwise PETSC_FALSE.
152 *
153 * If no functions are in the list, nothing is allowed by default.
154 * You can reverse this logic if you prefer to allow everything
155 * unless specified otherwise.
156 */
157PetscBool is_function_allowed(const char* functionName)
158{
159 /* no list ⇒ allow all */
160 if (gNumAllowed == 0) {
161 return PETSC_TRUE;
162 }
163
164 /* otherwise only the listed functions are allowed */
165 for (int i = 0; i < gNumAllowed; ++i) {
166 if (strcmp(gAllowedFunctions[i], functionName) == 0) {
167 return PETSC_TRUE;
168 }
169 }
170 return PETSC_FALSE;
171}
172
173/**
174 * @brief Prints the coordinates of a cell's vertices.
175 *
176 * This function iterates through the eight vertices of a given cell and prints their
177 * coordinates. It is primarily used for debugging purposes to verify the correctness
178 * of cell vertex assignments.
179 *
180 * @param[in] cell Pointer to a `Cell` structure representing the cell, containing its vertices.
181 * @param[in] rank MPI rank for identification (useful in parallel environments).
182 * @return PetscErrorCode Returns 0 to indicate successful execution. Non-zero on failure.
183 *
184 * @note
185 * - Ensure that the `cell` pointer is not `NULL` before calling this function..
186 */
187PetscErrorCode LOG_CELL_VERTICES(const Cell *cell, PetscMPIInt rank)
188{
189
190 // Validate input pointers
191 if (cell == NULL) {
192 LOG_ALLOW(LOCAL,LOG_ERROR, "'cell' is NULL.\n");
193 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "LOG_CELL_VERTICES - Input parameter 'cell' is NULL.");
194 }
195
196 LOG_ALLOW(LOCAL,LOG_VERBOSE, "Rank %d, Cell Vertices:\n", rank);
197 for(int i = 0; i < 8; i++){
198 LOG_ALLOW(LOCAL,LOG_VERBOSE, " Vertex[%d]: (%.2f, %.2f, %.2f)\n",
199 i, cell->vertices[i].x, cell->vertices[i].y, cell->vertices[i].z);
200 }
201
202 return 0; // Indicate successful execution
203}
204
205
206/**
207 * @brief Prints the signed distances to each face of the cell.
208 *
209 * This function iterates through the six signed distances from a point to each face of a given cell
210 * and prints their values. It is primarily used for debugging purposes to verify the correctness
211 * of distance calculations.
212 *
213 * @param[in] d An array of six `PetscReal` values representing the signed distances.
214 * The indices correspond to:
215 * - d[LEFT]: Left Face
216 * - d[RIGHT]: Right Face
217 * - d[BOTTOM]: Bottom Face
218 * - d[TOP]: Top Face
219 * - d[FRONT]: Front Face
220 * - d[BACK]: Back Face
221 *
222 * @return PetscErrorCode Returns 0 to indicate successful execution. Non-zero on failure.
223 *
224 * @note
225 * - Ensure that the `d` array is correctly populated with signed distances before calling this function.
226 */
227PetscErrorCode LOG_FACE_DISTANCES(PetscReal* d)
228{
229
230 // Validate input array
231 if (d == NULL) {
232 LOG_ALLOW(LOCAL,LOG_ERROR, " 'd' is NULL.\n");
233 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, " Input array 'd' is NULL.");
234 }
235
236 PetscPrintf(PETSC_COMM_SELF, " Face Distances:\n");
237 PetscPrintf(PETSC_COMM_SELF, " LEFT(%d): %.15f\n", LEFT, d[LEFT]);
238 PetscPrintf(PETSC_COMM_SELF, " RIGHT(%d): %.15f\n", RIGHT, d[RIGHT]);
239 PetscPrintf(PETSC_COMM_SELF, " BOTTOM(%d): %.15f\n", BOTTOM, d[BOTTOM]);
240 PetscPrintf(PETSC_COMM_SELF, " TOP(%d): %.15f\n", TOP, d[TOP]);
241 PetscPrintf(PETSC_COMM_SELF, " FRONT(%d): %.15f\n", FRONT, d[FRONT]);
242 PetscPrintf(PETSC_COMM_SELF, " BACK(%d): %.15f\n", BACK, d[BACK]);
243
244 return 0; // Indicate successful execution
245}
246
247/*
248 * Helper function: Converts an integer (of type int) to a string.
249 */
250static void IntToStr(int value, char *buf, size_t bufsize)
251{
252 snprintf(buf, bufsize, "%d", value);
253}
254
255/*
256 * Helper function: Converts a 64‐bit integer to a string.
257 */
258static void Int64ToStr(PetscInt64 value, char *buf, size_t bufsize)
259{
260 snprintf(buf, bufsize, "%ld", value);
261}
262
263/*
264 * Helper function: Converts three integers into a formatted string "(i, j, k)".
265 */
266static void CellToStr(const PetscInt *cell, char *buf, size_t bufsize)
267{
268 snprintf(buf, bufsize, "(%d, %d, %d)", cell[0], cell[1], cell[2]);
269}
270
271/*
272 * Helper function: Converts three PetscReal values into a formatted string "(x, y, z)".
273 */
274static void TripleRealToStr(const PetscReal *arr, char *buf, size_t bufsize)
275{
276 snprintf(buf, bufsize, "(%.4f, %.4f, %.4f)", arr[0], arr[1], arr[2]);
277}
278
279/*
280 * Helper function: Computes the maximum string length for each column (across all particles).
281 *
282 * The function examines every particle (from 0 to nParticles-1) and converts the value to a
283 * string using the helper functions above. The maximum length is stored in the pointers provided.
284 *
285 * @param nParticles Number of particles.
286 * @param ranks Array of particle MPI ranks.
287 * @param pids Array of particle IDs.
288 * @param cellIDs Array of cell IDs (stored consecutively, 3 per particle).
289 * @param positions Array of positions (3 per particle).
290 * @param velocities Array of velocities (3 per particle).
291 * @param weights Array of weights (3 per particle).
292 * @param wRank [out] Maximum width for Rank column.
293 * @param wPID [out] Maximum width for PID column.
294 * @param wCell [out] Maximum width for Cell column.
295 * @param wPos [out] Maximum width for Position column.
296 * @param wVel [out] Maximum width for Velocity column.
297 * @param wWt [out] Maximum width for Weights column.
298 */
299static PetscErrorCode ComputeMaxColumnWidths(PetscInt nParticles,
300 const PetscMPIInt *ranks,
301 const PetscInt64 *pids,
302 const PetscInt *cellIDs,
303 const PetscReal *positions,
304 const PetscReal *velocities,
305 const PetscReal *weights,
306 int *wRank, int *wPID, int *wCell,
307 int *wPos, int *wVel, int *wWt)
308{
309 char tmp[TMP_BUF_SIZE];
310
311 *wRank = strlen("Rank"); /* Start with the header label lengths */
312 *wPID = strlen("PID");
313 *wCell = strlen("Cell (i,j,k)");
314 *wPos = strlen("Position (x,y,z)");
315 *wVel = strlen("Velocity (x,y,z)");
316 *wWt = strlen("Weights (a1,a2,a3)");
317
318 for (PetscInt i = 0; i < nParticles; i++) {
319 /* Rank */
320 IntToStr(ranks[i], tmp, TMP_BUF_SIZE);
321 if ((int)strlen(tmp) > *wRank) *wRank = (int)strlen(tmp);
322
323 /* PID */
324 Int64ToStr(pids[i], tmp, TMP_BUF_SIZE);
325 if ((int)strlen(tmp) > *wPID) *wPID = (int)strlen(tmp);
326
327 /* Cell: use the three consecutive values */
328 CellToStr(&cellIDs[3 * i], tmp, TMP_BUF_SIZE);
329 if ((int)strlen(tmp) > *wCell) *wCell = (int)strlen(tmp);
330
331 /* Position */
332 TripleRealToStr(&positions[3 * i], tmp, TMP_BUF_SIZE);
333 if ((int)strlen(tmp) > *wPos) *wPos = (int)strlen(tmp);
334
335 /* Velocity */
336 TripleRealToStr(&velocities[3 * i], tmp, TMP_BUF_SIZE);
337 if ((int)strlen(tmp) > *wVel) *wVel = (int)strlen(tmp);
338
339 /* Weights */
340 TripleRealToStr(&weights[3 * i], tmp, TMP_BUF_SIZE);
341 if ((int)strlen(tmp) > *wWt) *wWt = (int)strlen(tmp);
342 }
343 return 0;
344}
345
346/*
347 * Helper function: Builds a format string for a table row.
348 *
349 * The format string will include proper width specifiers for each column.
350 * For example, it might create something like:
351 *
352 * "| %-6s | %-8s | %-20s | %-25s | %-25s | %-25s |\n"
353 *
354 * @param wRank Maximum width for the Rank column.
355 * @param wPID Maximum width for the PID column.
356 * @param wCell Maximum width for the Cell column.
357 * @param wPos Maximum width for the Position column.
358 * @param wVel Maximum width for the Velocity column.
359 * @param wWt Maximum width for the Weights column.
360 * @param fmtStr Buffer in which to build the format string.
361 * @param bufSize Size of fmtStr.
362 */
363static void BuildRowFormatString(PetscMPIInt wRank, PetscInt wPID, PetscInt wCell, PetscInt wPos, PetscInt wVel, PetscInt wWt, char *fmtStr, size_t bufSize)
364{
365 // Build a format string using snprintf.
366 // We assume that the Rank is an int (%d), PID is a 64-bit int (%ld)
367 // and the remaining columns are strings (which have been formatted already).
368 snprintf(fmtStr, bufSize,
369 "| %%-%dd | %%-%dd | %%-%ds | %%-%ds | %%-%ds | %%-%ds |\n",
370 wRank, wPID, wCell, wPos, wVel, wWt);
371}
372
373/*
374 * Helper function: Builds a header string for the table using column titles.
375 */
376static void BuildHeaderString(char *headerStr, size_t bufSize, PetscMPIInt wRank, PetscInt wPID, PetscInt wCell, PetscInt wPos, PetscInt wVel, PetscInt wWt)
377{
378 snprintf(headerStr, bufSize,
379 "| %-*s | %-*s | %-*s | %-*s | %-*s | %-*s |\n",
380 (int)wRank, "Rank",
381 (int)wPID, "PID",
382 (int)wCell, "Cell (i,j,k)",
383 (int)wPos, "Position (x,y,z)",
384 (int)wVel, "Velocity (x,y,z)",
385 (int)wWt, "Weights (a1,a2,a3)");
386}
387
388/**
389 * @brief Prints particle fields in a table that automatically adjusts its column widths.
390 *
391 * This function retrieves data from the particle swarm and prints a table where the
392 * width of each column is determined by the maximum width needed to display the data.
393 * Only every 'printInterval'-th particle is printed.
394 *
395 * @param[in] user Pointer to the UserCtx structure.
396 * @param[in] printInterval Only every printInterval‑th particle is printed.
397 *
398 * @return PetscErrorCode Returns 0 on success.
399 */
400PetscErrorCode LOG_PARTICLE_FIELDS(UserCtx* user, PetscInt printInterval)
401{
402 DM swarm = user->swarm;
403 PetscErrorCode ierr;
404 PetscInt localNumParticles;
405 PetscReal *positions = NULL;
406 PetscInt64 *particleIDs = NULL;
407 PetscMPIInt *particleRanks = NULL;
408 PetscInt *cellIDs = NULL;
409 PetscReal *weights = NULL;
410 PetscReal *velocities = NULL;
411 PetscMPIInt rank;
412
413 ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank); CHKERRQ(ierr);
414 LOG_ALLOW(LOCAL,LOG_INFO, "Rank %d is retrieving particle data.\n", rank);
415
416 ierr = DMSwarmGetLocalSize(swarm, &localNumParticles); CHKERRQ(ierr);
417 LOG_ALLOW(LOCAL,LOG_DEBUG,"Rank %d has %d particles.\n", rank, localNumParticles);
418
419 ierr = DMSwarmGetField(swarm, "position", NULL, NULL, (void**)&positions); CHKERRQ(ierr);
420 ierr = DMSwarmGetField(swarm, "DMSwarm_pid", NULL, NULL, (void**)&particleIDs); CHKERRQ(ierr);
421 ierr = DMSwarmGetField(swarm, "DMSwarm_rank", NULL, NULL, (void**)&particleRanks); CHKERRQ(ierr);
422 ierr = DMSwarmGetField(swarm, "DMSwarm_CellID", NULL, NULL, (void**)&cellIDs); CHKERRQ(ierr);
423 ierr = DMSwarmGetField(swarm, "weight", NULL, NULL, (void**)&weights); CHKERRQ(ierr);
424 ierr = DMSwarmGetField(swarm, "velocity", NULL, NULL, (void**)&velocities); CHKERRQ(ierr);
425
426 /* Compute maximum column widths. */
427 int wRank, wPID, wCell, wPos, wVel, wWt;
428 wRank = wPID = wCell = wPos = wVel = wWt = 0;
429 ierr = ComputeMaxColumnWidths(localNumParticles, particleRanks, particleIDs, cellIDs,
430 positions, velocities, weights,
431 &wRank, &wPID, &wCell, &wPos, &wVel, &wWt); CHKERRQ(ierr);
432
433 /* Build a header string and a row format string. */
434 char headerFmt[256];
435 char rowFmt[256];
436 BuildHeaderString(headerFmt, sizeof(headerFmt), wRank, wPID, wCell, wPos, wVel, wWt);
437 BuildRowFormatString(wRank, wPID, wCell, wPos, wVel, wWt, rowFmt, sizeof(rowFmt));
438
439 /* Print header (using synchronized printing for parallel output). */
440 ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, "--------------------------------------------------------------------------------------------------------------\n"); CHKERRQ(ierr);
441 ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s", headerFmt); CHKERRQ(ierr);
442 ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, "--------------------------------------------------------------------------------------------------------------\n"); CHKERRQ(ierr);
443
444 /* Loop over particles and print every printInterval-th row. */
445 char rowStr[256];
446 for (PetscInt i = 0; i < localNumParticles; i++) {
447 if (i % printInterval == 0) {
448 // ------- DEBUG
449 //char cellStr[TMP_BUF_SIZE], posStr[TMP_BUF_SIZE], velStr[TMP_BUF_SIZE], wtStr[TMP_BUF_SIZE];
450 //CellToStr(&cellIDs[3*i], cellStr, TMP_BUF_SIZE);
451 //TripleRealToStr(&positions[3*i], posStr, TMP_BUF_SIZE);
452 //TripleRealToStr(&velocities[3*i], velStr, TMP_BUF_SIZE);
453 // TripleRealToStr(&weights[3*i], wtStr, TMP_BUF_SIZE);
454
455 // if (rank == 0) { // Or whatever rank is Rank 0
456 //PetscPrintf(PETSC_COMM_SELF, "[Rank 0 DEBUG LPF] Particle %lld: PID=%lld, Rank=%d\n", (long long)i, (long long)particleIDs[i], particleRanks[i]);
457 //PetscPrintf(PETSC_COMM_SELF, "[Rank 0 DEBUG LPF] Raw Pos: (%.10e, %.10e, %.10e)\n", positions[3*i+0], positions[3*i+1], positions[3*i+2]);
458 //PetscPrintf(PETSC_COMM_SELF, "[Rank 0 DEBUG LPF] Str Pos: %s\n", posStr);
459 //PetscPrintf(PETSC_COMM_SELF, "[Rank 0 DEBUG LPF] Raw Vel: (%.10e, %.10e, %.10e)\n", velocities[3*i+0], velocities[3*i+1], velocities[3*i+2]);
460 // PetscPrintf(PETSC_COMM_SELF, "[Rank 0 DEBUG LPF] Str Vel: %s\n", velStr);
461 // Add similar for cell, weights
462 // PetscPrintf(PETSC_COMM_SELF, "[Rank 0 DEBUG LPF] About to build rowStr for particle %lld\n", (long long)i);
463 // fflush(stdout);
464 // }
465
466 // snprintf(rowStr, sizeof(rowStr), rowFmt,
467 // particleRanks[i],
468 // particleIDs[i],
469 // cellStr,
470 // posStr,
471 // velStr,
472 // wtStr);
473
474
475 // ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s", rowStr); CHKERRQ(ierr);
476
477 // ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s", rowStr); CHKERRQ(ierr);
478
479 // -------- DEBUG
480 /* Format the row by converting each field to a string first.
481 * We use temporary buffers and then build the row string.
482 */
483
484 char cellStr[TMP_BUF_SIZE], posStr[TMP_BUF_SIZE], velStr[TMP_BUF_SIZE], wtStr[TMP_BUF_SIZE];
485 CellToStr(&cellIDs[3*i], cellStr, TMP_BUF_SIZE);
486 TripleRealToStr(&positions[3*i], posStr, TMP_BUF_SIZE);
487 TripleRealToStr(&velocities[3*i], velStr, TMP_BUF_SIZE);
488 TripleRealToStr(&weights[3*i], wtStr, TMP_BUF_SIZE);
489
490 /* Build the row string. Note that for the integer fields we can use the row format string. */
491 snprintf(rowStr, sizeof(rowStr), rowFmt,
492 particleRanks[i],
493 particleIDs[i],
494 cellStr,
495 posStr,
496 velStr,
497 wtStr);
498 ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s", rowStr); CHKERRQ(ierr);
499 }
500 }
501
502
503 ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, "--------------------------------------------------------------------------------------------------------------\n"); CHKERRQ(ierr);
504 ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, "\n"); CHKERRQ(ierr);
505 ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT); CHKERRQ(ierr);
506
507 LOG_ALLOW_SYNC(GLOBAL,LOG_DEBUG,"Completed printing on Rank %d.\n", rank);
508
509 /* Restore fields */
510 ierr = DMSwarmRestoreField(swarm, "position", NULL, NULL, (void**)&positions); CHKERRQ(ierr);
511 ierr = DMSwarmRestoreField(swarm, "DMSwarm_pid", NULL, NULL, (void**)&particleIDs); CHKERRQ(ierr);
512 ierr = DMSwarmRestoreField(swarm, "DMSwarm_rank", NULL, NULL, (void**)&particleRanks); CHKERRQ(ierr);
513 ierr = DMSwarmRestoreField(swarm, "DMSwarm_CellID", NULL, NULL, (void**)&cellIDs); CHKERRQ(ierr);
514 ierr = DMSwarmRestoreField(swarm, "weight", NULL, NULL, (void**)&weights); CHKERRQ(ierr);
515 ierr = DMSwarmRestoreField(swarm, "velocity", NULL, NULL, (void**)&velocities); CHKERRQ(ierr);
516
517 LOG_ALLOW(LOCAL,LOG_DEBUG, "Restored all particle fields.\n");
518 return 0;
519}
520
521
522static void trim(char *s)
523{
524 if (!s) return;
525
526 /* ---- 1. strip leading blanks ----------------------------------- */
527 char *p = s;
528 while (*p && isspace((unsigned char)*p))
529 ++p;
530
531 if (p != s) /* move the trimmed text forward */
532 memmove(s, p, strlen(p) + 1); /* +1 to copy the final NUL */
533
534 /* ---- 2. strip trailing blanks ---------------------------------- */
535 size_t len = strlen(s);
536 while (len > 0 && isspace((unsigned char)s[len - 1]))
537 s[--len] = '\0';
538}
539
540/* ------------------------------------------------------------------------- */
541/**
542 * @brief Load function names from a text file.
543 *
544 * The file is expected to contain **one identifier per line**. Blank lines and
545 * lines whose first non‑blank character is a <tt>#</tt> are silently skipped so
546 * the file can include comments. Example:
547 *
548 * @code{.txt}
549 * # Allowed function list
550 * main
551 * InitializeSimulation
552 * InterpolateAllFieldsToSwarm # inline comments are OK, too
553 * @endcode
554 *
555 * The routine allocates memory as needed (growing an internal buffer with
556 * @c PetscRealloc()) and returns the resulting array and its length to the
557 * caller. Use FreeAllowedFunctions() to clean up when done.
558 *
559 * @param[in] filename Path of the configuration file to read.
560 * @param[out] funcsOut On success, points to a freshly‑allocated array of
561 * <tt>char*</tt> (size @p nOut).
562 * @param[out] nOut Number of valid entries in @p funcsOut.
563 *
564 * @return 0 on success, or a PETSc error code on failure (e.g. I/O error, OOM).
565 */
566PetscErrorCode LoadAllowedFunctionsFromFile(const char filename[],
567 char ***funcsOut,
568 PetscInt *nOut)
569{
570 FILE *fp = NULL;
571 char **funcs = NULL;
572 size_t cap = 16; /* initial capacity */
573 size_t n = 0; /* number of names */
574 char line[PETSC_MAX_PATH_LEN];
575 PetscErrorCode ierr;
576
577 PetscFunctionBegin;
578
579 /* ---------------------------------------------------------------------- */
580 /* 1. Open file */
581 fp = fopen(filename, "r");
582 if (!fp) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN,
583 "Cannot open %s", filename);
584
585 /* 2. Allocate initial pointer array */
586 ierr = PetscMalloc1(cap, &funcs); CHKERRQ(ierr);
587
588 /* 3. Read file line by line */
589 while (fgets(line, sizeof line, fp)) {
590 /* Strip everything after a comment character '#'. */
591 char *hash = strchr(line, '#');
592 if (hash) *hash = '\0';
593
594 trim(line); /* remove leading/trailing blanks */
595 if (!*line) continue; /* skip if empty */
596
597 /* Grow the array if necessary */
598 if (n == cap) {
599 cap *= 2;
600 ierr = PetscRealloc(cap * sizeof(*funcs), (void **)&funcs); CHKERRQ(ierr);
601 }
602
603 /* Deep‑copy the cleaned identifier */
604 ierr = PetscStrallocpy(line, &funcs[n++]); CHKERRQ(ierr);
605 }
606 fclose(fp);
607
608 /* 4. Return results to caller */
609 *funcsOut = funcs;
610 *nOut = (PetscInt)n;
611
612 PetscFunctionReturn(0);
613}
614
615/* ------------------------------------------------------------------------- */
616/**
617 * @brief Free an array previously returned by LoadAllowedFunctionsFromFile().
618 *
619 * @param[in,out] funcs Array of strings to release (may be @c NULL).
620 * @param[in] n Number of entries in @p funcs. Ignored if @p funcs is
621 * @c NULL.
622 *
623 * @return 0 on success or a PETSc error code.
624 */
625PetscErrorCode FreeAllowedFunctions(char **funcs, PetscInt n)
626{
627 PetscErrorCode ierr;
628 PetscFunctionBegin;
629 if (funcs) {
630 for (PetscInt i = 0; i < n; ++i) {
631 ierr = PetscFree(funcs[i]); CHKERRQ(ierr);
632 }
633 ierr = PetscFree(funcs); CHKERRQ(ierr);
634 }
635 PetscFunctionReturn(0);
636}
637
638/**
639 * @brief Helper function to convert BCFace enum to a string representation.
640 * @param[in] face The BCFace enum value.
641 * @return Pointer to a constant string representing the face.
642 */
643const char* BCFaceToString(BCFace face) {
644 switch (face) {
645 case BC_FACE_NEG_X: return "-Xi (I-Min)";
646 case BC_FACE_POS_X: return "+Xi (I-Max)";
647 case BC_FACE_NEG_Y: return "-Eta (J-Min)";
648 case BC_FACE_POS_Y: return "+Eta (J-Max)";
649 case BC_FACE_NEG_Z: return "-Zeta (K-Min)";
650 case BC_FACE_POS_Z: return "+Zeta (K-Max)";
651 default: return "Unknown Face";
652 }
653}
654
655/**
656 * @brief Helper function to convert FieldInitialization to a string representation.
657 * @param[in] PetscInt The FieldInitialization value.
658 * @return Pointer to a constant string representing the FieldInitialization.
659 */
660const char* FieldInitializationToString(PetscInt FieldInitialization)
661{
662 switch(FieldInitialization){
663 case 0: return "Zero";
664 case 1: return "Constant Normal velocity";
665 case 2: return "Poiseuille Normal velocity";
666 default: return "Unknown Field Initialization";
667 }
668}
669
670/**
671 * @brief Helper function to convert ParticleInitialization to a string representation.
672 * @param[in] ParticleInitializationType The ParticleInitialization value.
673 * @return Pointer to a constant string representing the ParticleInitialization.
674 */
676{
677 switch(ParticleInitialization){
678 case PARTICLE_INIT_SURFACE_RANDOM: return "Surface: Random";
679 case PARTICLE_INIT_VOLUME: return "Volume";
680 case PARTICLE_INIT_POINT_SOURCE: return "Point Source";
681 case PARTICLE_INIT_SURFACE_EDGES: return "Surface: At edges";
682 default: return "Unknown Particle Initialization";
683 }
684}
685
686/**
687 * @brief Helper function to convert LES Flag to a string representation.
688 * @param[in] LESModelType The LES flag value.
689 * @return Pointer to a constant string representing the LES Model Type.
690 */
691const char* LESModelToString(LESModelType LESFlag)
692{
693 switch(LESFlag){
694 case NO_LES_MODEL: return "No LES";
695 case CONSTANT_SMAGORINSKY: return "Constant Smagorinsky";
696 case DYNAMIC_SMAGORINSKY: return "Dynamic Smagorinsky";
697 default: return "Unknown LES Flag";
698 }
699}
700
701/**
702 * @brief Helper function to convert Momentum Solver flag to a string representation.
703 * @param[in] MomentumSolverType The Momentum Solver flag value.
704 * @return Pointer to a constant string representing the MomentumSolverType.
705 */
707{
708 switch(SolverFlag){
709 case MOMENTUM_SOLVER_EXPLICIT_RK: return "Explicit 4 stage Runge-Kutta ";
710 case MOMENTUM_SOLVER_DUALTIME_PICARD_RK4: return "Dual Time Stepping with Picard Iterations and 4 stage Runge-Kutta Smoothing";
711 case MOMENTUM_SOLVER_DUALTIME_NK_ANALYTIC_JACOBIAN: return "Dual Time Stepping with Newton-Krylov Iterations and Analytic Jacobian";
712 case MOMENTUM_SOLVER_DUALTIME_NK_ARNOLDI: return "Dual Time Stepping with Newton-Krylov Arnoldi Iteration Method.";
713 default: return "Unknown Momentum Solver Type";
714 }
715}
716
717/**
718 * @brief Helper function to convert BCType enum to a string representation.
719 * @param[in] type The BCType enum value.
720 * @return Pointer to a constant string representing the BC type.
721 */
722const char* BCTypeToString(BCType type) {
723 switch (type) {
724 // case DIRICHLET: return "DIRICHLET";
725 // case NEUMANN: return "NEUMANN";
726 case WALL: return "WALL";
727 case INLET: return "INLET";
728 case OUTLET: return "OUTLET";
729 case FARFIELD: return "FARFIELD";
730 case PERIODIC: return "PERIODIC";
731 case INTERFACE: return "INTERFACE";
732
733 // case CUSTOM: return "CUSTOM";
734 default: return "Unknown BC Type";
735 }
736}
737
738/**
739 * @brief Converts a BCHandlerType enum to its string representation.
740 *
741 * Provides a descriptive string for a specific boundary condition implementation strategy.
742 * This is crucial for logging the exact behavior configured for a face.
743 *
744 * @param handler_type The BCHandlerType enum value (e.g., BC_HANDLER_WALL_NOSLIP).
745 * @return A constant character string corresponding to the enum. Returns
746 * "UNKNOWN_HANDLER" if the enum value is not recognized.
747 */
748const char* BCHandlerTypeToString(BCHandlerType handler_type) {
749 switch (handler_type) {
750 // Wall & Symmetry Handlers
751 case BC_HANDLER_WALL_NOSLIP: return "noslip";
752 case BC_HANDLER_WALL_MOVING: return "moving";
753 case BC_HANDLER_SYMMETRY_PLANE: return "symmetry_plane";
754
755 // Inlet Handlers
756 case BC_HANDLER_INLET_CONSTANT_VELOCITY: return "constant_velocity";
757 case BC_HANDLER_INLET_PULSATILE_FLUX: return "pulsatile_flux";
758 case BC_HANDLER_INLET_PARABOLIC: return "parabolic";
759
760 // Outlet Handlers
761 case BC_HANDLER_OUTLET_CONSERVATION: return "conservation";
762 case BC_HANDLER_OUTLET_PRESSURE: return "pressure";
763
764 // Other Physical Handlers
765 case BC_HANDLER_FARFIELD_NONREFLECTING: return "nonreflecting";
766
767 // Multi-Block / Interface Handlers
768 case BC_HANDLER_PERIODIC_GEOMETRIC: return "geometric";
769 case BC_HANDLER_PERIODIC_DRIVEN_CONSTANT_FLUX: return "constant flux";
770 case BC_HANDLER_PERIODIC_DRIVEN_INITIAL_FLUX: return "initial flux";
771 case BC_HANDLER_INTERFACE_OVERSET: return "overset";
772
773 // Default case
775 default: return "UNKNOWN_HANDLER";
776 }
777}
778
779/**
780 * @brief Destroys the DualMonitorCtx.
781 *
782 * This function is passed to KSPMonitorSet to ensure the viewer is
783 * properly destroyed and the context memory is freed when the KSP is destroyed.
784 * @param Ctx a pointer to the context pointer to be destroyed
785 * @return PetscErrorCode
786 */
787PetscErrorCode DualMonitorDestroy(void **ctx)
788{
789 DualMonitorCtx *monctx = (DualMonitorCtx*)*ctx;
790 PetscErrorCode ierr;
791 PetscMPIInt rank;
792
793 PetscFunctionBeginUser;
794 ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank); CHKERRQ(ierr);
795 if(!rank && monctx->file_handle){
796 fclose(monctx->file_handle);
797 }
798
799 ierr = PetscFree(monctx); CHKERRQ(ierr);
800 *ctx = NULL;
801 PetscFunctionReturn(0);
802}
803
804/**
805 * @brief A custom KSP monitor that logs the true residual to a file and optionally to the console.
806 *
807 * This function replicates the behavior of KSPMonitorTrueResidualNorm by calculating
808 * the true residual norm ||b - Ax|| itself. It unconditionally logs to a file
809 * viewer and conditionally logs to the console based on a flag in the context.
810 *
811 * @param ksp The Krylov subspace context.
812 * @param it The current iteration number.
813 * @param rnorm The preconditioned residual norm (ignored, we compute our own).
814 * @param ctx A pointer to the DualMonitorCtx structure.
815 * @return PetscErrorCode 0 on success.
816 */
817#undef __FUNCT__
818#define __FUNCT__ "DualKSPMonitor"
819PetscErrorCode DualKSPMonitor(KSP ksp, PetscInt it, PetscReal rnorm, void *ctx)
820{
821 DualMonitorCtx *monctx = (DualMonitorCtx*)ctx;
822 PetscErrorCode ierr;
823 PetscReal trnorm, relnorm;
824 Vec r;
825 char norm_buf[256];
826 PetscMPIInt rank;
827
828 PetscFunctionBeginUser;
829 ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank); CHKERRQ(ierr);
830
831 // 1. Calculate the true residual norm.
832 ierr = KSPBuildResidual(ksp, NULL, NULL, &r); CHKERRQ(ierr);
833 ierr = VecNorm(r, NORM_2, &trnorm); CHKERRQ(ierr);
834
835 // 2. On the first iteration, compute and store the norm of the RHS vector `b`.
836 if (it == 0) {
837 Vec b;
838 ierr = KSPGetRhs(ksp, &b); CHKERRQ(ierr);
839 ierr = VecNorm(b, NORM_2, &monctx->bnorm); CHKERRQ(ierr);
840 }
841
842 if(!rank){
843 // 3. Compute the relative norm and format the output string.
844 if (monctx->bnorm > 1.e-15) {
845 relnorm = trnorm / monctx->bnorm;
846 sprintf(norm_buf, "ts: %-5d | block: %-2d | iter: %-3d | Unprecond Norm: %12.5e | True Norm: %12.5e | Rel Norm: %12.5e",(int)monctx->step, (int)monctx->block_id, (int)it, (double)rnorm, (double)trnorm, (double)relnorm);
847 } else {
848 sprintf(norm_buf,"ts: %-5d | block: %-2d | iter: %-3d | Unprecond Norm: %12.5e | True Norm: %12.5e",(int)monctx->step, (int)monctx->block_id, (int)it, (double)rnorm, (double)trnorm);
849 }
850
851 // 4. Log to the file viewer (unconditionally).
852 if(monctx->file_handle){
853 ierr = PetscFPrintf(PETSC_COMM_SELF,monctx->file_handle,"%s\n", norm_buf); CHKERRQ(ierr);
854 }
855 // 5. Log to the console (conditionally).
856 if (monctx->log_to_console) {
857 PetscFPrintf(PETSC_COMM_SELF,stdout, "%s\n", norm_buf); CHKERRQ(ierr);
858 }
859
860 } //rank
861
862 PetscFunctionReturn(0);
863}
864
865/**
866 * @brief Logs continuity metrics for a single block to a file.
867 *
868 * This function should be called for each block, once per timestep. It opens a
869 * central log file in append mode. To ensure the header is written only once,
870 * it checks if it is processing block 0 on the simulation's start step.
871 *
872 * @param user A pointer to the UserCtx for the specific block whose metrics
873 * are to be logged. The function accesses both global (SimCtx)
874 * and local (user->...) data.
875 * @return PetscErrorCode 0 on success.
876 */
877#undef __FUNCT__
878#define __FUNCT__ "LOG_CONTINUITY_METRICS"
879PetscErrorCode LOG_CONTINUITY_METRICS(UserCtx *user)
880{
881 PetscErrorCode ierr;
882 PetscMPIInt rank;
883 SimCtx *simCtx = user->simCtx; // Get the shared SimCtx
884 const PetscInt bi = user->_this; // Get this block's specific ID
885 const PetscInt ti = simCtx->step; // Get the current timestep
886
887 PetscFunctionBeginUser;
888 ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank); CHKERRQ(ierr);
889
890 // Only rank 0 performs file I/O.
891 if (!rank) {
892 FILE *f;
893 char filen[128];
894 sprintf(filen, "%s/Continuity_Metrics.log",simCtx->log_dir);
895
896 // Open the log file in append mode.
897 f = fopen(filen, "a");
898 if (!f) {
899 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open log file: %s", filen);
900 }
901
902 // Write a header only for the very first block (bi=0) on the very
903 // first timestep (ti=StartStep + 1). This ensures it's written only once.
904 if (ti == simCtx->StartStep + 1 && bi == 0) {
905 PetscFPrintf(PETSC_COMM_SELF, f, "%-10s | %-6s | %-18s | %-30s | %-18s | %-18s | %-18s | %-18s\n",
906 "Timestep", "Block", "Max Divergence", "Max Divergence Location ([k][j][i]=idx)", "Sum(RHS)","Total Flux In", "Total Flux Out", "Net Flux");
907 PetscFPrintf(PETSC_COMM_SELF, f, "------------------------------------------------------------------------------------------------------------------------------------------\n");
908 }
909
910 // Prepare the data strings and values for the current block.
911 PetscReal net_flux = simCtx->FluxInSum - simCtx->FluxOutSum;
912 char location_str[64];
913 sprintf(location_str, "([%d][%d][%d] = %d)", (int)simCtx->MaxDivz, (int)simCtx->MaxDivy, (int)simCtx->MaxDivx, (int)simCtx->MaxDivFlatArg);
914
915 // Write the formatted line for the current block.
916 PetscFPrintf(PETSC_COMM_SELF, f, "%-10d | %-6d | %-18.10e | %-39s | %-18.10e | %-18.10e | %-18.10e | %-18.10e\n",
917 (int)ti,
918 (int)bi,
919 (double)simCtx->MaxDiv,
920 location_str,
921 (double)simCtx->summationRHS,
922 (double)simCtx->FluxInSum,
923 (double)simCtx->FluxOutSum,
924 (double)net_flux);
925
926 fclose(f);
927 }
928
929 PetscFunctionReturn(0);
930}
931
932/**
933 * @brief A function that outputs the name of the current level in the ParticleLocation enum.
934 * @param level The ParticleLocation enum value.
935 * @return A constant character string corresponding to the enum. Returns
936 * "UNKNOWN_LEVEL" if the enum value is not recognized.
937 */
939{
940 switch (level) {
941 case NEEDS_LOCATION: return "NEEDS_LOCATION";
942 case ACTIVE_AND_LOCATED: return "ACTIVE_AND_LOCATED";
943 case MIGRATING_OUT: return "MIGRATING_OUT";
944 case LOST: return "LOST";
945 case UNINITIALIZED: return "UNINITIALIZED";
946 default: return "UNKNOWN_LEVEL";
947 }
948}
949
950///////// Profiling System /////////
951
952// Data structure to hold profiling info for one function
953typedef struct {
954 const char *name;
959 double start_time; // Timer for the current call
960 PetscBool always_log;
962
963// Global registry for all profiled functions
965static PetscInt g_profiler_count = 0;
966static PetscInt g_profiler_capacity = 0;
967
968// Internal helper to find a function in the registry or create it
969static PetscErrorCode _FindOrCreateEntry(const char *func_name, PetscInt *idx)
970{
971 PetscFunctionBeginUser;
972 // Search for existing entry
973 for (PetscInt i = 0; i < g_profiler_count; ++i) {
974 if (strcmp(g_profiler_registry[i].name, func_name) == 0) {
975 *idx = i;
976 PetscFunctionReturn(0);
977 }
978 }
979
980 // Not found, create a new entry
982 PetscInt new_capacity = g_profiler_capacity == 0 ? 16 : g_profiler_capacity * 2;
983 PetscErrorCode ierr = PetscRealloc(sizeof(ProfiledFunction) * new_capacity, &g_profiler_registry); CHKERRQ(ierr);
984 g_profiler_capacity = new_capacity;
985 }
986
987 *idx = g_profiler_count;
988 g_profiler_registry[*idx].name = func_name;
994 g_profiler_registry[*idx].always_log = PETSC_FALSE;
996
997 PetscFunctionReturn(0);
998}
999
1000// --- Public API Implementation ---
1001/**
1002 * @brief Initializes the custom profiling system using configuration from SimCtx.
1003 *
1004 * This function sets up the internal data structures for tracking function
1005 * performance. It reads the list of "critical functions" from the provided
1006 * SimCtx and marks them for per-step logging at LOG_INFO level.
1007 *
1008 * It should be called once at the beginning of the application, after
1009 * CreateSimulationContext() but before the main time loop.
1010 *
1011 * @param simCtx The master simulation context, which contains the list of
1012 * critical function names to always log.
1013 * @return PetscErrorCode
1014 */
1015PetscErrorCode ProfilingInitialize(SimCtx *simCtx)
1016{
1017 PetscFunctionBeginUser;
1018 if (!simCtx) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "SimCtx cannot be null for ProfilingInitialize");
1019
1020 // Iterate through the list of critical functions provided in SimCtx
1021 for (PetscInt i = 0; i < simCtx->nCriticalFuncs; ++i) {
1022 PetscInt idx;
1023 const char *func_name = simCtx->criticalFuncs[i];
1024 PetscErrorCode ierr = _FindOrCreateEntry(func_name, &idx); CHKERRQ(ierr);
1025 g_profiler_registry[idx].always_log = PETSC_TRUE;
1026
1027 LOG_ALLOW(GLOBAL, LOG_DEBUG, "Marked '%s' as a critical function for profiling.\n", func_name);
1028 }
1029 PetscFunctionReturn(0);
1030}
1031
1032void _ProfilingStart(const char *func_name)
1033{
1034 PetscInt idx;
1035 if (_FindOrCreateEntry(func_name, &idx) != 0) return; // Fail silently
1036 PetscTime(&g_profiler_registry[idx].start_time);
1037}
1038
1039void _ProfilingEnd(const char *func_name)
1040{
1041 double end_time;
1042 PetscTime(&end_time);
1043
1044 PetscInt idx;
1045 if (_FindOrCreateEntry(func_name, &idx) != 0) return; // Fail silently
1046
1047 double elapsed = end_time - g_profiler_registry[idx].start_time;
1048 g_profiler_registry[idx].total_time += elapsed;
1049 g_profiler_registry[idx].current_step_time += elapsed;
1052}
1053
1055{
1056 PetscFunctionBeginUser;
1057 for (PetscInt i = 0; i < g_profiler_count; ++i) {
1060 }
1061 PetscFunctionReturn(0);
1062}
1063
1064PetscErrorCode ProfilingLogTimestepSummary(PetscInt step)
1065{
1066 LogLevel log_level = get_log_level();
1067 PetscBool should_print = PETSC_FALSE;
1068
1069 PetscFunctionBeginUser;
1070
1071 // Decide if we should print anything at all
1072 if (log_level >= LOG_PROFILE) {
1073 for (PetscInt i = 0; i < g_profiler_count; ++i) {
1074 if (g_profiler_registry[i].current_step_call_count > 0) {
1075 if (log_level == LOG_PROFILE || g_profiler_registry[i].always_log) {
1076 should_print = PETSC_TRUE;
1077 break;
1078 }
1079 }
1080 }
1081 }
1082
1083 if (should_print) {
1084 PetscPrintf(PETSC_COMM_SELF, "[PROFILE] ----- Timestep %d Summary -----\n", step);
1085 for (PetscInt i = 0; i < g_profiler_count; ++i) {
1086 if (g_profiler_registry[i].current_step_call_count > 0) {
1087 if (log_level == LOG_PROFILE || g_profiler_registry[i].always_log) {
1088 PetscPrintf(PETSC_COMM_SELF, "[PROFILE] %-25s: %.6f s (%lld calls)\n",
1089 g_profiler_registry[i].name,
1090 g_profiler_registry[i].current_step_time,
1091 g_profiler_registry[i].current_step_call_count);
1092 }
1093 }
1094 }
1095 }
1096
1097 // Reset per-step counters for the next iteration
1098 for (PetscInt i = 0; i < g_profiler_count; ++i) {
1101 }
1102 PetscFunctionReturn(0);
1103}
1104
1105// Comparison function for qsort to sort by total_time in descending order
1106static int _CompareProfiledFunctions(const void *a, const void *b)
1107{
1108 const ProfiledFunction *func_a = (const ProfiledFunction *)a;
1109 const ProfiledFunction *func_b = (const ProfiledFunction *)b;
1110
1111 if (func_a->total_time < func_b->total_time) return 1;
1112 if (func_a->total_time > func_b->total_time) return -1;
1113 return 0;
1114}
1115
1116/**
1117 * @brief the profiling excercise and build a profiling summary which is then printed to a log file.
1118 *
1119 * @param simCtx The Simulation Context Structure that can contains all the data regarding the simulation.
1120 *
1121 * @return PetscErrorCode 0 on success.
1122 */
1123PetscErrorCode ProfilingFinalize(SimCtx *simCtx)
1124{
1125 PetscInt rank = simCtx->rank;
1126 PetscFunctionBeginUser;
1127 if (!rank) {
1128
1129 const char exec_mode_modifier[MAX_FILENAME_LENGTH];
1130 if(simCtx->exec_mode == EXEC_MODE_SOLVER) strcpy(exec_mode_modifier,"Solver");
1131 else if(simCtx->exec_mode == EXEC_MODE_POSTPROCESSOR) strcpy(exec_mode_modifier,"PostProcessor");
1132 //--- Step 0: Create a file viewer for log file
1133 FILE *f;
1134 char filen[128];
1135 sprintf(filen, "%s/ProfilingSummary_%s.log",simCtx->log_dir,exec_mode_modifier);
1136
1137 // Open the log file in write mode.
1138 f = fopen(filen,"w");
1139 if(!f){
1140 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot Open log file: %s",filen);
1141 }
1142
1143 // --- Step 1: Sort the data for readability ---
1145
1146 // --- Step 2: Dynamically determine the width for the function name column ---
1147 PetscInt max_name_len = strlen("Function"); // Start with the header's length
1148 for (PetscInt i = 0; i < g_profiler_count; ++i) {
1149 if (g_profiler_registry[i].total_call_count > 0) {
1150 PetscInt len = strlen(g_profiler_registry[i].name);
1151 if (len > max_name_len) {
1152 max_name_len = len;
1153 }
1154 }
1155 }
1156 // Add a little padding
1157 max_name_len += 2;
1158
1159 // --- Step 3: Define fixed widths for numeric columns for consistent alignment ---
1160 const int time_width = 18;
1161 const int count_width = 15;
1162 const int avg_width = 22;
1163
1164 // --- Step 4: Print the formatted table ---
1165 PetscFPrintf(PETSC_COMM_SELF, f, "=================================================================================================================\n");
1166 PetscFPrintf(PETSC_COMM_SELF, f, " FINAL PROFILING SUMMARY (Sorted by Total Time)\n");
1167 PetscFPrintf(PETSC_COMM_SELF, f, "=================================================================================================================\n");
1168
1169 // Header Row
1170 PetscFPrintf(PETSC_COMM_SELF, f, "%-*s | %-*s | %-*s | %-*s\n",
1171 max_name_len, "Function",
1172 time_width, "Total Time (s)",
1173 count_width, "Call Count",
1174 avg_width, "Avg. Time/Call (ms)");
1175
1176 // Separator Line (dynamically sized)
1177 for (int i = 0; i < max_name_len; i++) PetscFPrintf(PETSC_COMM_SELF, f, "-");
1178 PetscFPrintf(PETSC_COMM_SELF, f, "-|-");
1179 for (int i = 0; i < time_width; i++) PetscFPrintf(PETSC_COMM_SELF, f, "-");
1180 PetscFPrintf(PETSC_COMM_SELF, f, "-|-");
1181 for (int i = 0; i < count_width; i++) PetscFPrintf(PETSC_COMM_SELF, f, "-");
1182 PetscFPrintf(PETSC_COMM_SELF, f, "-|-");
1183 for (int i = 0; i < avg_width; i++) PetscFPrintf(PETSC_COMM_SELF, f, "-");
1184 PetscFPrintf(PETSC_COMM_SELF, f, "\n");
1185
1186 // Data Rows
1187 for (PetscInt i = 0; i < g_profiler_count; ++i) {
1188 if (g_profiler_registry[i].total_call_count > 0) {
1189 double avg_time_ms = (g_profiler_registry[i].total_time / g_profiler_registry[i].total_call_count) * 1000.0;
1190 PetscFPrintf(PETSC_COMM_SELF, f, "%-*s | %*.*f | %*lld | %*.*f\n",
1191 max_name_len, g_profiler_registry[i].name,
1192 time_width, 6, g_profiler_registry[i].total_time,
1193 count_width, g_profiler_registry[i].total_call_count,
1194 avg_width, 6, avg_time_ms);
1195 PetscFPrintf(PETSC_COMM_SELF, f, "------------------------------------------------------------------------------------------------------------------\n");
1196 }
1197 }
1198 PetscFPrintf(PETSC_COMM_SELF, f, "==================================================================================================================\n");
1199
1200 fclose(f);
1201 }
1202
1203 // --- Final Cleanup ---
1204 PetscFree(g_profiler_registry);
1205 g_profiler_registry = NULL;
1206 g_profiler_count = 0;
1208 PetscFunctionReturn(0);
1209}
1210
1211/*================================================================================*
1212 * PROGRESS BAR UTILITY *
1213 *================================================================================*/
1214
1215/**
1216 * @brief Prints a progress bar to the console.
1217 *
1218 * This function should only be called by the root process (rank 0). It uses
1219 * a carriage return `\r` to overwrite the same line in the terminal, creating
1220 * a dynamic progress bar.
1221 *
1222 * @param step The current step index from the loop (e.g., from 0 to N-1).
1223 * @param startStep The global starting step number of the simulation.
1224 * @param totalSteps The total number of steps to be run in this simulation instance.
1225 * @param currentTime The current simulation time to display.
1226 */
1227void PrintProgressBar(PetscInt step, PetscInt startStep, PetscInt totalSteps, PetscReal currentTime)
1228{
1229 if (totalSteps <= 0) return;
1230
1231 // --- Configuration ---
1232 const int barWidth = 50;
1233
1234 // --- Calculation ---
1235 // Calculate progress as a fraction from 0.0 to 1.0
1236 PetscReal progress = (PetscReal)(step - startStep + 1) / totalSteps;
1237 // Ensure progress doesn't exceed 1.0 due to floating point inaccuracies
1238 if (progress > 1.0) progress = 1.0;
1239
1240 int pos = (int)(barWidth * progress);
1241
1242 // --- Printing ---
1243 // Carriage return moves cursor to the beginning of the line
1244 PetscPrintf(PETSC_COMM_SELF, "\rProgress: [");
1245
1246 for (int i = 0; i < barWidth; ++i) {
1247 if (i < pos) {
1248 PetscPrintf(PETSC_COMM_SELF, "=");
1249 } else if (i == pos) {
1250 PetscPrintf(PETSC_COMM_SELF, ">");
1251 } else {
1252 PetscPrintf(PETSC_COMM_SELF, " ");
1253 }
1254 }
1255
1256 // Print percentage, step count, and current time
1257 PetscPrintf(PETSC_COMM_SELF, "] %3d%% (Step %ld/%ld, t=%.4f)",
1258 (int)(progress * 100.0),
1259 step + 1,
1260 startStep + totalSteps,
1261 currentTime);
1262
1263 // Flush the output buffer to ensure the bar is displayed immediately
1264 fflush(stdout);
1265}
1266
1267#undef __FUNCT__
1268#define __FUNCT__ "LOG_FIELD_MIN_MAX"
1269/**
1270 * @brief Computes and logs the local and global min/max values of a specified field,
1271 * respecting the solver's data layout architecture.
1272 *
1273 * This utility function inspects a PETSc Vec and calculates the minimum and maximum
1274 * values for its components, both locally and globally.
1275 *
1276 * It is "architecture-aware" and adjusts its iteration range to only include
1277 * physically meaningful data points:
1278 *
1279 * - **Cell-Centered Fields ("Ucat", "P"):** It uses the "Shifted Index Architecture,"
1280 * iterating from index 1 to N-1 to exclude the ghost/tool values at indices 0 and N.
1281 * - **Node-Centered Fields ("Coordinates"):** It iterates from index 0 to N-1, covering all
1282 * physical nodes.
1283 * - **Face-Centered Fields ("Ucont"):** It iterates from index 0 to N-1, covering all
1284 * physical faces.
1285 *
1286 * The results are printed to the standard output in a formatted, easy-to-read table.
1287 *
1288 * @param[in] user Pointer to the user-defined context.
1289 * @param[in] fieldName A string descriptor for the field being analyzed.
1290 *
1291 * @return PetscErrorCode Returns 0 on success, non-zero on failure.
1292 */
1293PetscErrorCode LOG_FIELD_MIN_MAX(UserCtx *user, const char *fieldName)
1294{
1295 PetscErrorCode ierr;
1296 PetscInt i, j, k;
1297 DMDALocalInfo info;
1298
1299 Vec fieldVec = NULL;
1300 DM dm = NULL;
1301 PetscInt dof;
1302 char data_layout[20];
1303
1304 PetscFunctionBeginUser;
1305
1306 // --- 1. Map string name to PETSc objects and determine data layout ---
1307 if (strcasecmp(fieldName, "Ucat") == 0) {
1308 fieldVec = user->Ucat; dm = user->fda; dof = 3; strcpy(data_layout, "Cell-Centered");
1309 } else if (strcasecmp(fieldName, "P") == 0) {
1310 fieldVec = user->P; dm = user->da; dof = 1; strcpy(data_layout, "Cell-Centered");
1311 } else if (strcasecmp(fieldName, "Diffusivity") == 0) {
1312 fieldVec = user->Diffusivity; dm = user->da; dof = 1; strcpy(data_layout, "Cell-Centered");
1313 } else if (strcasecmp(fieldName, "DiffusivityGradient") == 0) {
1314 fieldVec = user->DiffusivityGradient; dm = user->fda; dof = 3; strcpy(data_layout, "Cell-Centered");
1315 } else if (strcasecmp(fieldName, "Ucont") == 0) {
1316 fieldVec = user->lUcont; dm = user->fda; dof = 3; strcpy(data_layout, "Face-Centered");
1317 } else if (strcasecmp(fieldName, "Coordinates") == 0) {
1318 ierr = DMGetCoordinates(user->da, &fieldVec); CHKERRQ(ierr);
1319 dm = user->fda; dof = 3; strcpy(data_layout, "Node-Centered");
1320 } else if (strcasecmp(fieldName,"Psi") == 0) {
1321 fieldVec = user->Psi; dm = user->da; dof = 1; strcpy(data_layout, "Cell-Centered"); // Assuming Psi is cell-centered
1322 } else {
1323 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Field %s not recognized.", fieldName);
1324 }
1325
1326 if (!fieldVec) {
1327 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Vector for field '%s' is NULL.", fieldName);
1328 }
1329 if (!dm) {
1330 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "DM for field '%s' is NULL.", fieldName);
1331 }
1332
1333 ierr = DMDAGetLocalInfo(dm, &info); CHKERRQ(ierr);
1334
1335 // --- 2. Define Architecture-Aware Loop Bounds ---
1336 PetscInt i_start, i_end, j_start, j_end, k_start, k_end;
1337
1338 if (strcmp(data_layout, "Cell-Centered") == 0) {
1339 // For cell-centered data, the physical values are stored from index 1 to N-1.
1340 // We find the intersection of the rank's owned range [xs, xe) with the
1341 // physical data range [1, IM-1).
1342 i_start = PetscMax(info.xs, 1); i_end = PetscMin(info.xs + info.xm, user->IM);
1343 j_start = PetscMax(info.ys, 1); j_end = PetscMin(info.ys + info.ym, user->JM);
1344 k_start = PetscMax(info.zs, 1); k_end = PetscMin(info.zs + info.zm, user->KM);
1345 } else { // For Node- or Face-Centered data
1346 // The physical values are stored from index 0 to N-1.
1347 // We find the intersection of the rank's owned range [xs, xe) with the
1348 // physical data range [0, IM-1].
1349 i_start = PetscMax(info.xs, 0); i_end = PetscMin(info.xs + info.xm, user->IM);
1350 j_start = PetscMax(info.ys, 0); j_end = PetscMin(info.ys + info.ym, user->JM);
1351 k_start = PetscMax(info.zs, 0); k_end = PetscMin(info.zs + info.zm, user->KM);
1352 }
1353
1354 // --- 3. Barrier for clean, grouped output ---
1355 ierr = MPI_Barrier(PETSC_COMM_WORLD); CHKERRQ(ierr);
1356 if (user->simCtx->rank == 0) {
1357 PetscPrintf(PETSC_COMM_SELF, "\n--- Field Ranges: [%s] (Layout: %s) ---\n", fieldName, data_layout);
1358 }
1359
1360 // --- 4. Branch on DoF and perform calculation with correct bounds ---
1361 if (dof == 1) {
1362 PetscReal localMin = PETSC_MAX_REAL, localMax = PETSC_MIN_REAL;
1363 PetscReal globalMin, globalMax;
1364 const PetscScalar ***array;
1365
1366 ierr = DMDAVecGetArrayRead(dm, fieldVec, &array); CHKERRQ(ierr);
1367 for (k = k_start; k < k_end; k++) {
1368 for (j = j_start; j < j_end; j++) {
1369 for (i = i_start; i < i_end; i++) {
1370 localMin = PetscMin(localMin, array[k][j][i]);
1371 localMax = PetscMax(localMax, array[k][j][i]);
1372 }
1373 }
1374 }
1375 ierr = DMDAVecRestoreArrayRead(dm, fieldVec, &array); CHKERRQ(ierr);
1376
1377 ierr = MPI_Allreduce(&localMin, &globalMin, 1, MPIU_REAL, MPI_MIN, PETSC_COMM_WORLD); CHKERRQ(ierr);
1378 ierr = MPI_Allreduce(&localMax, &globalMax, 1, MPIU_REAL, MPI_MAX, PETSC_COMM_WORLD); CHKERRQ(ierr);
1379
1380 PetscSynchronizedPrintf(PETSC_COMM_WORLD, " [Rank %d] Local Range: [ %11.4e , %11.4e ]\n", user->simCtx->rank, localMin, localMax);
1381 ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT); CHKERRQ(ierr);
1382 if (user->simCtx->rank == 0) {
1383 PetscPrintf(PETSC_COMM_SELF, " Global Range: [ %11.4e , %11.4e ]\n", globalMin, globalMax);
1384 }
1385
1386 } else if (dof == 3) {
1387 Cmpnts localMin = {PETSC_MAX_REAL, PETSC_MAX_REAL, PETSC_MAX_REAL};
1388 Cmpnts localMax = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL};
1389 Cmpnts globalMin, globalMax;
1390 const Cmpnts ***array;
1391
1392 ierr = DMDAVecGetArrayRead(dm, fieldVec, &array); CHKERRQ(ierr);
1393 for (k = k_start; k < k_end; k++) {
1394 for (j = j_start; j < j_end; j++) {
1395 for (i = i_start; i < i_end; i++) {
1396 localMin.x = PetscMin(localMin.x, array[k][j][i].x);
1397 localMin.y = PetscMin(localMin.y, array[k][j][i].y);
1398 localMin.z = PetscMin(localMin.z, array[k][j][i].z);
1399 localMax.x = PetscMax(localMax.x, array[k][j][i].x);
1400 localMax.y = PetscMax(localMax.y, array[k][j][i].y);
1401 localMax.z = PetscMax(localMax.z, array[k][j][i].z);
1402 }
1403 }
1404 }
1405 ierr = DMDAVecRestoreArrayRead(dm, fieldVec, &array); CHKERRQ(ierr);
1406
1407 ierr = MPI_Allreduce(&localMin, &globalMin, 3, MPIU_REAL, MPI_MIN, PETSC_COMM_WORLD); CHKERRQ(ierr);
1408 ierr = MPI_Allreduce(&localMax, &globalMax, 3, MPIU_REAL, MPI_MAX, PETSC_COMM_WORLD); CHKERRQ(ierr);
1409
1410 ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, " [Rank %d] Local X-Range: [ %11.4e , %11.4e ]\n", user->simCtx->rank, localMin.x, localMax.x);
1411 ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, " [Rank %d] Local Y-Range: [ %11.4e , %11.4e ]\n", user->simCtx->rank, localMin.y, localMax.y);
1412 ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, " [Rank %d] Local Z-Range: [ %11.4e , %11.4e ]\n", user->simCtx->rank, localMin.z, localMax.z);
1413 ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT); CHKERRQ(ierr);
1414
1415 if (user->simCtx->rank == 0) {
1416 PetscPrintf(PETSC_COMM_SELF, " [Global] X-Range: [ %11.4e , %11.4e ]\n", globalMin.x, globalMax.x);
1417 PetscPrintf(PETSC_COMM_SELF, " [Global] Y-Range: [ %11.4e , %11.4e ]\n", globalMin.y, globalMax.y);
1418 PetscPrintf(PETSC_COMM_SELF, " [Global] Z-Range: [ %11.4e , %11.4e ]\n", globalMin.z, globalMax.z);
1419 }
1420
1421 } else {
1422 SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "LogFieldStatistics only supports fields with 1 or 3 components, but field '%s' has %D.", fieldName, dof);
1423 }
1424
1425 // --- 5. Final barrier for clean output ordering ---
1426 ierr = MPI_Barrier(PETSC_COMM_WORLD); CHKERRQ(ierr);
1427 if (user->simCtx->rank == 0) {
1428 PetscPrintf(PETSC_COMM_SELF, "--------------------------------------------\n\n");
1429 }
1430
1431 PetscFunctionReturn(0);
1432}
1433
1434#undef __FUNCT__
1435#define __FUNCT__ "LOG_FIELD_ANATOMY"
1436/**
1437 * @brief Logs the anatomy of a specified field at key boundary locations,
1438 * respecting the solver's specific grid and variable architecture.
1439 *
1440 * This intelligent diagnostic function inspects a PETSc Vec and prints its values
1441 * at critical boundary locations (-Xi/+Xi, -Eta/+Eta, -Zeta/+Zeta). It is "architecture-aware":
1442 *
1443 * - **Cell-Centered Fields ("Ucat", "P"):** Correctly applies the "Shifted Index Architecture,"
1444 * where the value for geometric `Cell i` is stored at array index `i+1`.
1445 * - **Face-Centered Fields ("Ucont", "Csi"):** Recognizes the staggered nature. For example, an
1446 * X-face-centered field is treated as node-centered in the I-direction but cell-centered
1447 * in the J and K directions. For vector fields like "Ucont", each component is handled
1448 * with its appropriate directional anatomy.
1449 * - **Node-Centered Fields ("Coordinates"):** Uses a direct index mapping, where the value for
1450 * `Node i` is stored at index `i`.
1451 *
1452 * The output is synchronized across MPI ranks and focuses on a slice through the center of
1453 * the domain to be concise.
1454 *
1455 * @param user A pointer to the UserCtx structure containing the DMs and Vecs.
1456 * @param field_name A string identifier for the field to log (e.g., "Ucat", "P", "Ucont", "Coordinates").
1457 * @param stage_name A string identifier for the current simulation stage (e.g., "After Advection").
1458 * @return PetscErrorCode Returns 0 on success, non-zero on failure.
1459 */
1460PetscErrorCode LOG_FIELD_ANATOMY(UserCtx *user, const char *field_name, const char *stage_name)
1461{
1462 PetscErrorCode ierr;
1463 DMDALocalInfo info;
1464 PetscMPIInt rank;
1465
1466 Vec vec_local = NULL;
1467 DM dm = NULL;
1468 PetscInt dof = 0;
1469 char data_layout[20];
1470 char dominant_dir = '\0'; // 'x', 'y', 'z' for face-centered, 'm' for mixed (Ucont)
1471
1472 PetscFunctionBeginUser;
1473 ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank); CHKERRQ(ierr);
1474
1475 // --- 1. Map string name to PETSc objects and determine data layout ---
1476 if (strcasecmp(field_name, "Ucat") == 0) {
1477 vec_local = user->lUcat; dm = user->fda; dof = 3; strcpy(data_layout, "Cell-Centered");
1478 } else if (strcasecmp(field_name, "P") == 0) {
1479 vec_local = user->lP; dm = user->da; dof = 1; strcpy(data_layout, "Cell-Centered");
1480 } else if (strcasecmp(field_name, "Diffusivity") == 0) {
1481 vec_local = user->lDiffusivity; dm = user->da; dof = 1; strcpy(data_layout, "Cell-Centered");
1482 } else if (strcasecmp(field_name, "DiffusivityGradient") == 0) {
1483 vec_local = user->lDiffusivityGradient; dm = user->fda; dof = 3; strcpy(data_layout, "Cell-Centered");
1484 } else if (strcasecmp(field_name, "Psi") == 0) {
1485 vec_local = user->lPsi; dm = user->da; dof = 1; strcpy(data_layout, "Cell-Centered");
1486 } else if (strcasecmp(field_name, "Center-Coordinates") == 0) {
1487 vec_local = user->lCent; dm = user->fda; dof = 3; strcpy(data_layout, "Cell-Centered");
1488 } else if (strcasecmp(field_name, "Ucont") == 0) {
1489 vec_local = user->lUcont; dm = user->fda; dof = 3; strcpy(data_layout, "Face-Centered"); dominant_dir = 'm'; // Mixed
1490 } else if (strcasecmp(field_name, "Csi") == 0 || strcasecmp(field_name, "X-Face-Centers") == 0) {
1491 vec_local = (strcasecmp(field_name, "Csi") == 0) ? user->lCsi : user->Centx;
1492 dm = user->fda; dof = 3; strcpy(data_layout, "Face-Centered"); dominant_dir = 'x';
1493 } else if (strcasecmp(field_name, "Eta") == 0 || strcasecmp(field_name, "Y-Face-Centers") == 0) {
1494 vec_local = (strcasecmp(field_name, "Eta") == 0) ? user->lEta : user->Centy;
1495 dm = user->fda; dof = 3; strcpy(data_layout, "Face-Centered"); dominant_dir = 'y';
1496 } else if (strcasecmp(field_name, "Zet") == 0 || strcasecmp(field_name, "Z-Face-Centers") == 0) {
1497 vec_local = (strcasecmp(field_name, "Zet") == 0) ? user->lZet : user->Centz;
1498 dm = user->fda; dof = 3; strcpy(data_layout, "Face-Centered"); dominant_dir = 'z';
1499 } else if (strcasecmp(field_name, "Coordinates") == 0) {
1500 ierr = DMGetCoordinatesLocal(user->da, &vec_local); CHKERRQ(ierr);
1501 dm = user->fda; dof = 3; strcpy(data_layout, "Node-Centered");
1502 } else if (strcasecmp(field_name, "CornerField")== 0){
1503 vec_local = user->lCellFieldAtCorner; strcpy(data_layout, "Node-Centered");
1504 PetscInt bs = 1;
1505 ierr = VecGetBlockSize(user->CellFieldAtCorner, &bs); CHKERRQ(ierr);
1506 dof = bs;
1507 if(dof == 1) dm = user->da;
1508 else dm = user->fda;
1509 } else {
1510 SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Unknown field name for LOG_FIELD_ANATOMY: %s", field_name);
1511 }
1512
1513 // --- 2. Get Grid Info and Array Pointers ---
1514 ierr = DMDAGetLocalInfo(dm, &info); CHKERRQ(ierr);
1515
1516 ierr = PetscBarrier(NULL);
1517 PetscPrintf(PETSC_COMM_WORLD, "\n--- Field Anatomy Log: [%s] | Stage: [%s] | Layout: [%s] ---\n", field_name, stage_name, data_layout);
1518
1519 // Global physical dimensions (number of cells)
1520 PetscInt im_phys = user->IM;
1521 PetscInt jm_phys = user->JM;
1522 PetscInt km_phys = user->KM;
1523
1524 // Slice through the center of the local domain
1525 PetscInt i_mid = (PetscInt)(info.xs + 0.5 * info.xm) - 1;
1526 PetscInt j_mid = (PetscInt)(info.ys + 0.5 * info.ym) - 1;
1527 PetscInt k_mid = (PetscInt)(info.zs + 0.5 * info.zm) - 1;
1528
1529 // --- 3. Print Boundary Information based on Data Layout ---
1530
1531 // ======================================================================
1532 // === CASE 1: Cell-Centered Fields (Ucat, P) - USES SHIFTED INDEX ===
1533 // ======================================================================
1534 if (strcmp(data_layout, "Cell-Centered") == 0) {
1535 const void *l_arr;
1536 ierr = DMDAVecGetArrayRead(dm, vec_local, (void*)&l_arr); CHKERRQ(ierr);
1537
1538
1539 // --- I-Direction Boundaries ---
1540 if (info.xs == 0) { // Rank on -Xi boundary
1541 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, I-DIR]: Idx %2d (Ghost for Cell[k][j][0]) = ", rank, 0);
1542 if(dof==1) PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.5f)\n", ((const PetscReal***)l_arr)[k_mid][j_mid][0]);
1543 else PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.5f, %.5f, %.5f)\n", ((const Cmpnts***)l_arr)[k_mid][j_mid][0].x, ((const Cmpnts***)l_arr)[k_mid][j_mid][0].y, ((const Cmpnts***)l_arr)[k_mid][j_mid][0].z);
1544
1545 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, I-DIR]: Idx %2d (Value for Cell[k][j][0]) = ", rank, 1);
1546 if(dof==1) PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.5f)\n", ((const PetscReal***)l_arr)[k_mid][j_mid][1]);
1547 else PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.5f, %.5f, %.5f)\n", ((const Cmpnts***)l_arr)[k_mid][j_mid][1].x, ((const Cmpnts***)l_arr)[k_mid][j_mid][1].y, ((const Cmpnts***)l_arr)[k_mid][j_mid][1].z);
1548 }
1549 if (info.xs + info.xm == info.mx) { // Rank on +Xi boundary
1550 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, I-DIR]: Idx %2d (Value for Cell[k][j][%d]) = ", rank, im_phys - 1, im_phys - 2);
1551 if(dof==1) PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.5f)\n", ((const PetscReal***)l_arr)[k_mid][j_mid][im_phys - 1]);
1552 else PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.5f, %.5f, %.5f)\n", ((const Cmpnts***)l_arr)[k_mid][j_mid][im_phys - 1].x, ((const Cmpnts***)l_arr)[k_mid][j_mid][im_phys - 1].y, ((const Cmpnts***)l_arr)[k_mid][j_mid][im_phys - 1].z);
1553
1554 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, I-DIR]: Idx %2d (Ghost for Cell[k][j][%d]) = ", rank, im_phys, im_phys - 2);
1555 if(dof==1) PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.5f)\n", ((const PetscReal***)l_arr)[k_mid][j_mid][im_phys]);
1556 else PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.5f, %.5f, %.5f)\n", ((const Cmpnts***)l_arr)[k_mid][j_mid][im_phys].x, ((const Cmpnts***)l_arr)[k_mid][j_mid][im_phys].y, ((const Cmpnts***)l_arr)[k_mid][j_mid][im_phys].z);
1557 }
1558
1559 // --- J-Direction Boundaries ---
1560 if (info.ys == 0) { // Rank on -Eta boundary
1561 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, J-DIR]: Jdx %2d (Ghost for Cell[k][0][i]) = ", rank, 0);
1562 if(dof==1) PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.5f)\n", ((const PetscReal***)l_arr)[k_mid][0][i_mid]);
1563 else PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.5f, %.5f, %.5f)\n", ((const Cmpnts***)l_arr)[k_mid][0][i_mid].x, ((const Cmpnts***)l_arr)[k_mid][0][i_mid].y, ((const Cmpnts***)l_arr)[k_mid][0][i_mid].z);
1564
1565 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, J-DIR]: Jdx %2d (Value for Cell[k][0][i]) = ", rank, 1);
1566 if(dof==1) PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.5f)\n", ((const PetscReal***)l_arr)[k_mid][1][i_mid]);
1567 else PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.5f, %.5f, %.5f)\n", ((const Cmpnts***)l_arr)[k_mid][1][i_mid].x, ((const Cmpnts***)l_arr)[k_mid][1][i_mid].y, ((const Cmpnts***)l_arr)[k_mid][1][i_mid].z);
1568 }
1569
1570 if (info.ys + info.ym == info.my) { // Rank on +Eta boundary
1571 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, J-DIR]: Jdx %2d (Value for Cell[k][%d][i]) = ", rank, jm_phys - 1, jm_phys - 2);
1572 if(dof==1) PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.5f)\n", ((const PetscReal***)l_arr)[k_mid][jm_phys - 1][i_mid]);
1573 else PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.5f, %.5f, %.5f)\n", ((const Cmpnts***)l_arr)[k_mid][jm_phys - 1][i_mid].x, ((const Cmpnts***)l_arr)[k_mid][jm_phys - 1][i_mid].y, ((const Cmpnts***)l_arr)[k_mid][jm_phys - 1][i_mid].z);
1574
1575 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, J-DIR]: Jdx %2d (Ghost for Cell[k][%d][i]) = ", rank, jm_phys, jm_phys - 2);
1576 if(dof==1) PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.5f)\n", ((const PetscReal***)l_arr)[k_mid][jm_phys][i_mid]);
1577 else PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.5f, %.5f, %.5f)\n", ((const Cmpnts***)l_arr)[k_mid][jm_phys][i_mid].x, ((const Cmpnts***)l_arr)[k_mid][jm_phys][i_mid].y, ((const Cmpnts***)l_arr)[k_mid][jm_phys][i_mid].z);
1578 }
1579
1580 // --- K-Direction Boundaries ---
1581 if (info.zs == 0) { // Rank on -Zeta boundary
1582 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, K-DIR]: Kdx %2d (Ghost for Cell[0][j][i]) = ", rank, 0);
1583 if(dof==1) PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.5f)\n", ((const PetscReal***)l_arr)[0][j_mid][i_mid]);
1584 else PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.5f, %.5f, %.5f)\n", ((const Cmpnts***)l_arr)[0][j_mid][i_mid].x, ((const Cmpnts***)l_arr)[0][j_mid][i_mid].y, ((const Cmpnts***)l_arr)[0][j_mid][i_mid].z);
1585 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, K-DIR]: Kdx %2d (Value for Cell[0][j][i]) = ", rank, 1);
1586 if(dof==1) PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.5f)\n", ((const PetscReal***)l_arr)[1][j_mid][i_mid]);
1587 else PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.5f, %.5f, %.5f)\n", ((const Cmpnts***)l_arr)[1][j_mid][i_mid].x, ((const Cmpnts***)l_arr)[1][j_mid][i_mid].y, ((const Cmpnts***)l_arr)[1][j_mid][i_mid].z);
1588 }
1589 if (info.zs + info.zm == info.mz) { // Rank on +Zeta boundary
1590 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, K-DIR]: Kdx %2d (Value for Cell[%d][j][i]) = ", rank, km_phys - 1, km_phys - 2);
1591 if(dof==1) PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.5f)\n", ((const PetscReal***)l_arr)[km_phys - 1][j_mid][i_mid]);
1592 else PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.5f, %.5f, %.5f)\n", ((const Cmpnts***)l_arr)[km_phys - 1][j_mid][i_mid].x, ((const Cmpnts***)l_arr)[km_phys - 1][j_mid][i_mid].y, ((const Cmpnts***)l_arr)[km_phys - 1][j_mid][i_mid].z);
1593 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, K-DIR]: Kdx %2d (Ghost for Cell[%d][j][i]) = ", rank, km_phys, km_phys - 2);
1594 if(dof==1) PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.5f)\n", ((const PetscReal***)l_arr)[km_phys][j_mid][i_mid]);
1595 else PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.5f, %.5f, %.5f)\n", ((const Cmpnts***)l_arr)[km_phys][j_mid][i_mid].x, ((const Cmpnts***)l_arr)[km_phys][j_mid][i_mid].y, ((const Cmpnts***)l_arr)[km_phys][j_mid][i_mid].z);
1596 }
1597 ierr = DMDAVecRestoreArrayRead(dm, vec_local, (void*)&l_arr); CHKERRQ(ierr);
1598 }
1599 // ======================================================================
1600 // === CASE 2: Face-Centered Fields - NUANCED DIRECTIONAL LOGIC ===
1601 // ======================================================================
1602 else if (strcmp(data_layout, "Face-Centered") == 0) {
1603 const Cmpnts ***l_arr;
1604 ierr = DMDAVecGetArrayRead(dm, vec_local, (void*)&l_arr); CHKERRQ(ierr);
1605
1606 // --- I-Direction Boundaries ---
1607 if (info.xs == 0) { // Rank on -Xi boundary
1608 if (dominant_dir == 'x') { // Node-like in I-dir
1609 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, I-DIR]: Idx %2d (First Phys. X-Face) = (%.5f, %.5f, %.5f)\n", rank, 0, l_arr[k_mid][j_mid][0].x, l_arr[k_mid][j_mid][0].y, l_arr[k_mid][j_mid][0].z);
1610 } else if (dominant_dir == 'y' || dominant_dir == 'z') { // Cell-like in I-dir
1611 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, I-DIR]: Idx %2d (Ghost for Cell[k][j][0]) = (%.5f, %.5f, %.5f)\n", rank, 0, l_arr[k_mid][j_mid][0].x, l_arr[k_mid][j_mid][0].y, l_arr[k_mid][j_mid][0].z);
1612 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, I-DIR]: Idx %2d (Value for Cell[k][j][0]) = (%.5f, %.5f, %.5f)\n", rank, 1, l_arr[k_mid][j_mid][1].x, l_arr[k_mid][j_mid][1].y, l_arr[k_mid][j_mid][1].z);
1613 } else if (dominant_dir == 'm') { // Ucont: Mixed
1614 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, I-DIR]: u-comp @ Idx %2d (1st X-Face) = %.5f\n", rank, 0, l_arr[k_mid][j_mid][0].x);
1615 }
1616 }
1617 if (info.xs + info.xm == info.mx) { // Rank on +Xi boundary
1618 if (dominant_dir == 'x') { // Node-like in I-dir
1619 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, I-DIR]: Idx %2d (Last Phys. X-Face) = (%.5f, %.5f, %.5f)\n", rank, im_phys - 1, l_arr[k_mid][j_mid][im_phys - 1].x, l_arr[k_mid][j_mid][im_phys-1].y, l_arr[k_mid][j_mid][im_phys - 1].z);
1620 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, I-DIR]: Idx %2d (Ghost Location) = (%.5f, %.5f, %.5f)\n", rank, im_phys, l_arr[k_mid][j_mid][im_phys].x, l_arr[k_mid][j_mid][im_phys].y, l_arr[k_mid][j_mid][im_phys].z);
1621 } else if (dominant_dir == 'y' || dominant_dir == 'z') { // Cell-like in I-dir
1622 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, I-DIR]: Idx %2d (Value for Cell[k][j][%d]) = (%.5f, %.5f, %.5f)\n", rank, im_phys - 1, im_phys - 2, l_arr[k_mid][j_mid][im_phys - 1].x, l_arr[k_mid][j_mid][im_phys - 1].y, l_arr[k_mid][j_mid][im_phys-1].z);
1623 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, I-DIR]: Idx %2d (Ghost for Cell[k][j][%d]) = (%.5f, %.5f, %.5f)\n", rank, im_phys, im_phys - 2, l_arr[k_mid][j_mid][im_phys].x, l_arr[k_mid][j_mid][im_phys].y, l_arr[k_mid][j_mid][im_phys].z);
1624 } else if (dominant_dir == 'm') { // Ucont: Mixed
1625 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, I-DIR]: u-comp @ Idx %2d (Last X-Face) = %.5f\n", rank, im_phys - 1, l_arr[k_mid][j_mid][im_phys - 1].x);
1626 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, I-DIR]: u-comp @ Idx %2d (Ghost Location) = %.5f\n", rank, im_phys, l_arr[k_mid][j_mid][im_phys].x);
1627 }
1628 }
1629
1630 // --- J-Direction Boundaries ---
1631 if (info.ys == 0) { // Rank on -Eta boundary
1632 if (dominant_dir == 'y') { // Node-like in J-dir
1633 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, J-DIR]: Jdx %2d (First Phys. Y-Face) = (%.5f, %.5f, %.5f)\n", rank, 0, l_arr[k_mid][0][i_mid].x, l_arr[k_mid][0][i_mid].y, l_arr[k_mid][0][i_mid].z);
1634 } else if (dominant_dir == 'x' || dominant_dir == 'z') { // Cell-like in J-dir
1635 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, J-DIR]: Jdx %2d (Ghost for Cell[k][0][i]) = (%.5f, %.5f, %.5f)\n", rank, 0, l_arr[k_mid][0][i_mid].x, l_arr[k_mid][0][i_mid].y, l_arr[k_mid][0][i_mid].z);
1636 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, J-DIR]: Jdx %2d (Value for Cell[k][0][i]) = (%.5f, %.5f, %.5f)\n", rank, 1, l_arr[k_mid][1][i_mid].x, l_arr[k_mid][1][i_mid].y, l_arr[k_mid][1][i_mid].z);
1637 } else if (dominant_dir == 'm') { // Ucont: Mixed
1638 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, J-DIR]: v-comp @ Jdx %2d (1st Y-Face) = %.5f\n", rank, 0, l_arr[k_mid][0][i_mid].y);
1639 }
1640 }
1641 if (info.ys + info.ym == info.my) { // Rank on +Eta boundary
1642 if (dominant_dir == 'y') { // Node-like in J-dir
1643 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, J-DIR]: Jdx %2d (Last Phys. Y-Face) = (%.5f, %.5f, %.5f)\n", rank, jm_phys - 1, l_arr[k_mid][jm_phys - 1][i_mid].x, l_arr[k_mid][jm_phys - 1][i_mid].y, l_arr[k_mid][jm_phys - 1][i_mid].z);
1644 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, J-DIR]: Jdx %2d (Ghost Location) = (%.5f, %.5f, %.5f)\n", rank, jm_phys, l_arr[k_mid][jm_phys][i_mid].x, l_arr[k_mid][jm_phys][i_mid].y, l_arr[k_mid][jm_phys][i_mid].z);
1645 } else if (dominant_dir == 'x' || dominant_dir == 'z') { // Cell-like in J-dir
1646 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, J-DIR]: Jdx %2d (Value for Cell[k][%d][i]) = (%.5f, %.5f, %.5f)\n", rank, jm_phys-1, jm_phys-2, l_arr[k_mid][jm_phys - 1][i_mid].x, l_arr[k_mid][jm_phys - 1][i_mid].y, l_arr[k_mid][jm_phys - 1][i_mid].z);
1647 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, J-DIR]: Jdx %2d (Ghost for Cell[k][%d][i]) = (%.5f, %.5f, %.5f)\n", rank, jm_phys, jm_phys-2, l_arr[k_mid][jm_phys][i_mid].x, l_arr[k_mid][jm_phys][i_mid].y, l_arr[k_mid][jm_phys][i_mid].z);
1648 } else if (dominant_dir == 'm') { // Ucont: Mixed
1649 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, J-DIR]: v-comp @ Jdx %2d (Last Y-Face) = %.5f\n", rank, jm_phys - 1, l_arr[k_mid][jm_phys - 1][i_mid].y);
1650 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, J-DIR]: v-comp @ Jdx %2d (Ghost Location) = %.5f\n", rank, jm_phys, l_arr[k_mid][jm_phys][i_mid].y);
1651 }
1652 }
1653
1654 // --- K-Direction Boundaries ---
1655 if (info.zs == 0) { // Rank on -Zeta boundary
1656 if (dominant_dir == 'z') { // Node-like in K-dir
1657 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, K-DIR]: Kdx %2d (First Phys. Z-Face) = (%.5f, %.5f, %.5f)\n", rank, 0, l_arr[0][j_mid][i_mid].x, l_arr[0][j_mid][i_mid].y, l_arr[0][j_mid][i_mid].z);
1658 } else if (dominant_dir == 'x' || dominant_dir == 'y') { // Cell-like in K-dir
1659 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, K-DIR]: Kdx %2d (Ghost for Cell[0][j][i]) = (%.5f, %.5f, %.5f)\n", rank, 0, l_arr[0][j_mid][i_mid].x, l_arr[0][j_mid][i_mid].y, l_arr[0][j_mid][i_mid].z);
1660 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, K-DIR]: Kdx %2d (Value for Cell[0][j][i]) = (%.5f, %.5f, %.5f)\n", rank, 1, l_arr[1][j_mid][i_mid].x, l_arr[1][j_mid][i_mid].y, l_arr[1][j_mid][i_mid].z);
1661 } else if (dominant_dir == 'm') { // Ucont: Mixed
1662 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, K-DIR]: w-comp @ Idx %2d (1st Z-Face) = %.5f\n", rank, 0, l_arr[0][j_mid][i_mid].z);
1663 }
1664 }
1665 if (info.zs + info.zm == info.mz) { // Rank on +Zeta boundary
1666 if (dominant_dir == 'z') { // Node-like in K-dir
1667 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, K-DIR]: Idx %2d (Last Phys. Z-Face) = (%.5f, %.5f, %.5f)\n", rank, km_phys - 1, l_arr[km_phys - 1][j_mid][i_mid].x, l_arr[km_phys - 1][j_mid][i_mid].y, l_arr[km_phys - 1][j_mid][i_mid].z);
1668 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, K-DIR]: Idx %2d (Ghost Location) = (%.5f, %.5f, %.5f)\n", rank, km_phys, l_arr[km_phys][j_mid][i_mid].x, l_arr[km_phys][j_mid][i_mid].y, l_arr[km_phys][j_mid][i_mid].z);
1669 } else if (dominant_dir == 'x' || dominant_dir == 'y') { // Cell-like in K-dir
1670 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, K-DIR]: Idx %2d (Value for Cell[%d][j][i]) = (%.5f, %.5f, %.5f)\n", rank, km_phys-1, km_phys-2, l_arr[km_phys-1][j_mid][i_mid].x, l_arr[km_phys-1][j_mid][i_mid].y, l_arr[km_phys - 1][j_mid][i_mid].z);
1671 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, K-DIR]: Idx %2d (Ghost for Cell[%d][j][i]) = (%.5f, %.5f, %.5f)\n", rank, km_phys, km_phys-2, l_arr[km_phys][j_mid][i_mid].x, l_arr[km_phys][j_mid][i_mid].y, l_arr[km_phys][j_mid][i_mid].z);
1672 } else if (dominant_dir == 'm') { // Ucont: Mixed
1673 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, K-DIR]: w-comp @ Idx %2d (Last Z-Face) = %.5f\n", rank, km_phys - 1, l_arr[km_phys - 1][j_mid][i_mid].z);
1674 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, K-DIR]: w-comp @ Idx %2d (Ghost Loc.) = %.5f\n", rank, km_phys, l_arr[km_phys][j_mid][i_mid].z);
1675
1676 }
1677 }
1678 ierr = DMDAVecRestoreArrayRead(dm, vec_local, (void*)&l_arr); CHKERRQ(ierr);
1679 }
1680 // ======================================================================
1681 // === CASE 3: Node-Centered Fields - USES DIRECT INDEX ===
1682 // ======================================================================
1683 else if (strcmp(data_layout, "Node-Centered") == 0) {
1684 const Cmpnts ***l_arr;
1685 ierr = DMDAVecGetArrayRead(dm, vec_local, (void*)&l_arr); CHKERRQ(ierr);
1686
1687 // --- I-Direction Boundaries ---
1688 if (info.xs == 0) {
1689 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, I-DIR]: Idx %2d (First Phys. Node) = (%.5f, %.5f, %.5f)\n", rank, 0, l_arr[k_mid][j_mid][0].x, l_arr[k_mid][j_mid][0].y, l_arr[k_mid][j_mid][0].z);
1690 }
1691 if (info.xs + info.xm == info.mx) {
1692 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, I-DIR]: Idx %2d (Last Phys. Node) = (%.5f, %.5f, %.5f)\n", rank, im_phys - 1, l_arr[k_mid][j_mid][im_phys - 1].x, l_arr[k_mid][j_mid][im_phys - 1].y, l_arr[k_mid][j_mid][im_phys - 1].z);
1693 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, I-DIR]: Idx %2d (Unused/Ghost Loc) = (%.5f, %.5f, %.5f)\n", rank, im_phys, l_arr[k_mid][j_mid][im_phys].x, l_arr[k_mid][j_mid][im_phys].y, l_arr[k_mid][j_mid][im_phys].z);
1694 }
1695 // --- J-Direction Boundaries ---
1696 if (info.ys == 0) {
1697 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, J-DIR]: Jdx %2d (First Phys. Node) = (%.5f, %.5f, %.5f)\n", rank, 0, l_arr[k_mid][0][i_mid].x, l_arr[k_mid][0][i_mid].y, l_arr[k_mid][0][i_mid].z);
1698 }
1699 if (info.ys + info.ym == info.my) {
1700 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, J-DIR]: Jdx %2d (Last Phys. Node) = (%.5f, %.5f, %.5f)\n", rank, jm_phys - 1, l_arr[k_mid][jm_phys - 1][i_mid].x, l_arr[k_mid][jm_phys - 1][i_mid].y, l_arr[k_mid][jm_phys - 1][i_mid].z);
1701 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, J-DIR]: Jdx %2d (Unused/Ghost Loc) = (%.5f, %.5f, %.5f)\n", rank, jm_phys, l_arr[k_mid][jm_phys][i_mid].x, l_arr[k_mid][jm_phys][i_mid].y, l_arr[k_mid][jm_phys][i_mid].z);
1702 }
1703 // --- K-Direction Boundaries ---
1704 if (info.zs == 0) {
1705 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, K-DIR]: Kdx %2d (First Phys. Node) = (%.5f, %.5f, %.5f)\n", rank, 0, l_arr[0][j_mid][i_mid].x, l_arr[0][j_mid][i_mid].y, l_arr[0][j_mid][i_mid].z);
1706 }
1707 if(info.zs + info.zm == info.mz) {
1708 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, K-DIR]: Kdx %2d (Last Phys. Node) = (%.5f, %.5f, %.5f)\n", rank, km_phys - 1, l_arr[km_phys - 1][j_mid][i_mid].x, l_arr[km_phys - 1][j_mid][i_mid].y, l_arr[km_phys - 1][j_mid][i_mid].z);
1709 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, K-DIR]: Kdx %2d (Unused/Ghost Loc) = (%.5f, %.5f, %.5f)\n", rank, km_phys, l_arr[km_phys][j_mid][i_mid].x, l_arr[km_phys][j_mid][i_mid].y, l_arr[km_phys][j_mid][i_mid].z);
1710 }
1711 ierr = DMDAVecRestoreArrayRead(dm, vec_local, (void*)&l_arr); CHKERRQ(ierr);
1712 }
1713 else {
1714 SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "LOG_FIELD_ANATOMY encountered an unknown data layout: %s", data_layout);
1715 }
1716
1717 ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT); CHKERRQ(ierr);
1718 ierr = PetscBarrier(NULL);
1719 PetscFunctionReturn(0);
1720}
1721
1722#undef __FUNCT__
1723#define __FUNCT__ "LOG_INTERPOLATION_ERROR"
1724/**
1725@brief Logs the interpolation error between the analytical and computed solutions.
1726*/
1728{
1729 SimCtx *simCtx = user->simCtx;
1730 PetscErrorCode ierr;
1731 DM swarm = user->swarm;
1732 Vec positionVec, analyticalvelocityVec, velocityVec, errorVec;
1733 PetscReal Interpolation_error = 0.0;
1734 PetscReal Maximum_Interpolation_error = 0.0;
1735 PetscReal AnalyticalSolution_magnitude = 0.0;
1736 PetscReal ErrorPercentage = 0.0;
1737
1738 LOG_ALLOW(GLOBAL, LOG_DEBUG, "Creating global vectors.\n");
1739 ierr = DMSwarmCreateGlobalVectorFromField(swarm, "position", &positionVec); CHKERRQ(ierr);
1740 ierr = DMSwarmCreateGlobalVectorFromField(swarm, "velocity", &velocityVec); CHKERRQ(ierr);
1741
1742 ierr = VecDuplicate(positionVec, &analyticalvelocityVec); CHKERRQ(ierr);
1743 ierr = VecCopy(positionVec, analyticalvelocityVec); CHKERRQ(ierr);
1744
1745 LOG_ALLOW(GLOBAL, LOG_DEBUG, "Computing analytical solution.\n");
1746 ierr = SetAnalyticalSolutionForParticles(analyticalvelocityVec, simCtx); CHKERRQ(ierr);
1747
1748 ierr = VecDuplicate(analyticalvelocityVec, &errorVec); CHKERRQ(ierr);
1749 ierr = VecCopy(analyticalvelocityVec, errorVec); CHKERRQ(ierr);
1750
1751 ierr = VecNorm(analyticalvelocityVec, NORM_2, &AnalyticalSolution_magnitude); CHKERRQ(ierr);
1752
1753 LOG_ALLOW(GLOBAL, LOG_DEBUG, "Computing error.\n");
1754 ierr = VecAXPY(errorVec, -1.0, velocityVec); CHKERRQ(ierr);
1755 ierr = VecNorm(errorVec, NORM_2, &Interpolation_error); CHKERRQ(ierr);
1756 ierr = VecNorm(errorVec,NORM_INFINITY,&Maximum_Interpolation_error); CHKERRQ(ierr);
1757
1758 ErrorPercentage = (AnalyticalSolution_magnitude > 0) ?
1759 (Interpolation_error / AnalyticalSolution_magnitude * 100.0) : 0.0;
1760
1761 LOG_ALLOW(GLOBAL, LOG_INFO, "Interpolation error (%%): %g\n", ErrorPercentage);
1762 PetscPrintf(PETSC_COMM_WORLD, "Interpolation error (%%): %g\n", ErrorPercentage);
1763 LOG_ALLOW(GLOBAL, LOG_INFO, "Maximum Interpolation error: %g\n", Maximum_Interpolation_error);
1764 PetscPrintf(PETSC_COMM_WORLD, "Maximum Interpolation error: %g\n", Maximum_Interpolation_error);
1765
1766
1767 ierr = VecDestroy(&analyticalvelocityVec); CHKERRQ(ierr);
1768 ierr = VecDestroy(&errorVec); CHKERRQ(ierr);
1769 ierr = DMSwarmDestroyGlobalVectorFromField(swarm, "position", &positionVec); CHKERRQ(ierr);
1770 ierr = DMSwarmDestroyGlobalVectorFromField(swarm, "velocity", &velocityVec); CHKERRQ(ierr);
1771
1772 return 0;
1773}
1774
1775#undef __FUNCT__
1776#define __FUNCT__ "CalculateAdvancedParticleMetrics"
1777/**
1778 * @brief Computes advanced particle statistics and stores them in SimCtx.
1779 *
1780 * This function calculates:
1781 * - Particle load imbalance across MPI ranks.
1782 * - The total number of grid cells occupied by at least one particle.
1783 *
1784 * It requires that CalculateParticleCountPerCell() has been called prior to its
1785 * execution. It uses collective MPI operations and must be called by all ranks.
1786 *
1787 * @param user Pointer to the UserCtx.
1788 * @return PetscErrorCode 0 on success.
1789 */
1791{
1792 PetscErrorCode ierr;
1793 SimCtx *simCtx = user->simCtx;
1794 PetscMPIInt size, rank;
1795
1796 PetscFunctionBeginUser;
1797 ierr = MPI_Comm_size(PETSC_COMM_WORLD, &size); CHKERRQ(ierr);
1798 ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank); CHKERRQ(ierr);
1799
1800 // --- 1. Particle Load Imbalance ---
1801 PetscInt nLocal, nGlobal, nLocalMax;
1802 ierr = DMSwarmGetLocalSize(user->swarm, &nLocal); CHKERRQ(ierr);
1803 ierr = DMSwarmGetSize(user->swarm, &nGlobal); CHKERRQ(ierr);
1804 ierr = MPI_Allreduce(&nLocal, &nLocalMax, 1, MPIU_INT, MPI_MAX, PETSC_COMM_WORLD); CHKERRQ(ierr);
1805
1806 PetscReal avg_per_rank = (size > 0) ? ((PetscReal)nGlobal / size) : 0.0;
1807 // Handle division by zero if there are no particles
1808 simCtx->particleLoadImbalance = (avg_per_rank > 1e-9) ? (nLocalMax / avg_per_rank) : 1.0;
1809
1810
1811 // --- 2. Number of Occupied Cells ---
1812 // This part requires access to the user->ParticleCount vector.
1813 PetscInt local_occupied_cells = 0;
1814 PetscInt global_occupied_cells;
1815 const PetscScalar *count_array;
1816 PetscInt vec_local_size;
1817
1818 ierr = VecGetLocalSize(user->ParticleCount, &vec_local_size); CHKERRQ(ierr);
1819 ierr = VecGetArrayRead(user->ParticleCount, &count_array); CHKERRQ(ierr);
1820
1821 for (PetscInt i = 0; i < vec_local_size; ++i) {
1822 if (count_array[i] > 0.5) { // Use 0.5 to be safe with floating point
1823 local_occupied_cells++;
1824 }
1825 }
1826 ierr = VecRestoreArrayRead(user->ParticleCount, &count_array); CHKERRQ(ierr);
1827
1828 ierr = MPI_Allreduce(&local_occupied_cells, &global_occupied_cells, 1, MPIU_INT, MPI_SUM, PETSC_COMM_WORLD); CHKERRQ(ierr);
1829 simCtx->occupiedCellCount = global_occupied_cells;
1830
1831 LOG_ALLOW_SYNC(GLOBAL, LOG_INFO, "[Rank %d] Advanced Metrics: Imbalance=%.2f, OccupiedCells=%d\n", rank, simCtx->particleLoadImbalance, simCtx->occupiedCellCount);
1832
1833 PetscFunctionReturn(0);
1834}
1835
1836#undef __FUNCT__
1837#define __FUNCT__ "LOG_PARTICLE_METRICS"
1838/**
1839 * @brief Logs particle swarm metrics, adapting its behavior based on a boolean flag in SimCtx.
1840 *
1841 * This function serves a dual purpose:
1842 * 1. If simCtx->isInitializationPhase is PETSC_TRUE, it logs settlement
1843 * diagnostics to "Initialization_Metrics.log", using the provided stageName.
1844 * 2. If simCtx->isInitializationPhase is PETSC_FALSE, it logs regular
1845 * timestep metrics to "Particle_Metrics.log".
1846 *
1847 * @param user A pointer to the UserCtx.
1848 * @param stageName A descriptive string for the initialization stage (ignored in timestep mode).
1849 * @return PetscErrorCode 0 on success.
1850 */
1851PetscErrorCode LOG_PARTICLE_METRICS(UserCtx *user, const char *stageName)
1852{
1853 PetscErrorCode ierr;
1854 PetscMPIInt rank, size;
1855 SimCtx *simCtx = user->simCtx;
1856
1857 PetscFunctionBeginUser;
1858 ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank); CHKERRQ(ierr);
1859 ierr = MPI_Comm_size(PETSC_COMM_WORLD, &size); CHKERRQ(ierr);
1860
1861 PetscInt totalParticles;
1862 ierr = DMSwarmGetSize(user->swarm, &totalParticles); CHKERRQ(ierr);
1863
1864 if (!rank) {
1865 FILE *f;
1866 char filen[PETSC_MAX_PATH_LEN];
1867 sprintf(filen, "%s/Particle_Metrics.log", simCtx->log_dir);
1868 f = fopen(filen, "a");
1869 if (!f) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open particle log file: %s", filen);
1870
1871 if (simCtx->step == simCtx->StartStep + 1) {
1872 PetscFPrintf(PETSC_COMM_SELF, f, "%-10s | %-12s | %-10s | %-10s | %-15s | %-10s | %-10s\n",
1873 "Timestep", "Total Ptls", "Lost", "Migrated", "Occupied Cells", "Imbalance", "Mig Passes");
1874 PetscFPrintf(PETSC_COMM_SELF, f, "----------------------------------------------------------------------------------------------------------\n");
1875 }
1876
1877 PetscFPrintf(PETSC_COMM_SELF, f, "%-10d | %-12d | %-10d | %-10d | %-15d | %-10.2f | %-10d\n",
1878 (int)simCtx->step, (int)totalParticles, (int)simCtx->particlesLostLastStep,
1879 (int)simCtx->particlesMigratedLastStep, (int)simCtx->occupiedCellCount,
1880 (double)simCtx->particleLoadImbalance, (int)simCtx->migrationPassesLastStep);
1881 fclose(f);
1882 }
1883 PetscFunctionReturn(0);
1884}
PetscErrorCode SetAnalyticalSolutionForParticles(Vec tempVec, SimCtx *simCtx)
Applies the analytical solution to particle velocity vector.
void set_allowed_functions(const char **functionList, int count)
Sets the global list of function names that are allowed to log.
Definition logging.c:122
PetscBool always_log
Definition logging.c:960
PetscErrorCode LOG_PARTICLE_METRICS(UserCtx *user, const char *stageName)
Logs particle swarm metrics, adapting its behavior based on a boolean flag in SimCtx.
Definition logging.c:1851
const char * BCHandlerTypeToString(BCHandlerType handler_type)
Converts a BCHandlerType enum to its string representation.
Definition logging.c:748
PetscBool is_function_allowed(const char *functionName)
Checks if the given function name is in the allow-list.
Definition logging.c:157
static PetscInt g_profiler_count
Definition logging.c:965
const char * FieldInitializationToString(PetscInt FieldInitialization)
Helper function to convert FieldInitialization to a string representation.
Definition logging.c:660
PetscErrorCode DualMonitorDestroy(void **ctx)
Destroys the DualMonitorCtx.
Definition logging.c:787
#define TMP_BUF_SIZE
Definition logging.c:5
static char ** gAllowedFunctions
Global/static array of function names allowed to log.
Definition logging.c:21
static LogLevel current_log_level
Static variable to cache the current logging level.
Definition logging.c:14
PetscErrorCode LOG_INTERPOLATION_ERROR(UserCtx *user)
Logs the interpolation error between the analytical and computed solutions.
Definition logging.c:1727
const char * BCFaceToString(BCFace face)
Helper function to convert BCFace enum to a string representation.
Definition logging.c:643
PetscErrorCode FreeAllowedFunctions(char **funcs, PetscInt n)
Free an array previously returned by LoadAllowedFunctionsFromFile().
Definition logging.c:625
PetscErrorCode print_log_level(void)
Prints the current logging level to the console.
Definition logging.c:82
long long total_call_count
Definition logging.c:957
static PetscInt g_profiler_capacity
Definition logging.c:966
PetscErrorCode ProfilingFinalize(SimCtx *simCtx)
the profiling excercise and build a profiling summary which is then printed to a log file.
Definition logging.c:1123
static void BuildRowFormatString(PetscMPIInt wRank, PetscInt wPID, PetscInt wCell, PetscInt wPos, PetscInt wVel, PetscInt wWt, char *fmtStr, size_t bufSize)
Definition logging.c:363
static void trim(char *s)
Definition logging.c:522
PetscErrorCode LoadAllowedFunctionsFromFile(const char filename[], char ***funcsOut, PetscInt *nOut)
Load function names from a text file.
Definition logging.c:566
PetscErrorCode LOG_FIELD_MIN_MAX(UserCtx *user, const char *fieldName)
Computes and logs the local and global min/max values of a specified field, respecting the solver's d...
Definition logging.c:1293
static int gNumAllowed
Number of entries in the gAllowedFunctions array.
Definition logging.c:26
double total_time
Definition logging.c:955
static void BuildHeaderString(char *headerStr, size_t bufSize, PetscMPIInt wRank, PetscInt wPID, PetscInt wCell, PetscInt wPos, PetscInt wVel, PetscInt wWt)
Definition logging.c:376
void PrintProgressBar(PetscInt step, PetscInt startStep, PetscInt totalSteps, PetscReal currentTime)
Prints a progress bar to the console.
Definition logging.c:1227
PetscErrorCode LOG_FIELD_ANATOMY(UserCtx *user, const char *field_name, const char *stage_name)
Logs the anatomy of a specified field at key boundary locations, respecting the solver's specific gri...
Definition logging.c:1460
static PetscErrorCode _FindOrCreateEntry(const char *func_name, PetscInt *idx)
Definition logging.c:969
static void CellToStr(const PetscInt *cell, char *buf, size_t bufsize)
Definition logging.c:266
LogLevel get_log_level()
Retrieves the current logging level from the environment variable LOG_LEVEL.
Definition logging.c:39
PetscErrorCode LOG_FACE_DISTANCES(PetscReal *d)
Prints the signed distances to each face of the cell.
Definition logging.c:227
PetscErrorCode LOG_PARTICLE_FIELDS(UserCtx *user, PetscInt printInterval)
Prints particle fields in a table that automatically adjusts its column widths.
Definition logging.c:400
void _ProfilingEnd(const char *func_name)
Definition logging.c:1039
static void TripleRealToStr(const PetscReal *arr, char *buf, size_t bufsize)
Definition logging.c:274
static void Int64ToStr(PetscInt64 value, char *buf, size_t bufsize)
Definition logging.c:258
const char * BCTypeToString(BCType type)
Helper function to convert BCType enum to a string representation.
Definition logging.c:722
PetscErrorCode CalculateAdvancedParticleMetrics(UserCtx *user)
Computes advanced particle statistics and stores them in SimCtx.
Definition logging.c:1790
const char * ParticleLocationStatusToString(ParticleLocationStatus level)
A function that outputs the name of the current level in the ParticleLocation enum.
Definition logging.c:938
static int _CompareProfiledFunctions(const void *a, const void *b)
Definition logging.c:1106
PetscErrorCode DualKSPMonitor(KSP ksp, PetscInt it, PetscReal rnorm, void *ctx)
A custom KSP monitor that logs to a file and optionally to the console.
Definition logging.c:819
PetscErrorCode LOG_CONTINUITY_METRICS(UserCtx *user)
Logs continuity metrics for a single block to a file.
Definition logging.c:879
double current_step_time
Definition logging.c:956
long long current_step_call_count
Definition logging.c:958
static ProfiledFunction * g_profiler_registry
Definition logging.c:964
PetscErrorCode ProfilingInitialize(SimCtx *simCtx)
Initializes the custom profiling system using configuration from SimCtx.
Definition logging.c:1015
const char * LESModelToString(LESModelType LESFlag)
Helper function to convert LES Flag to a string representation.
Definition logging.c:691
PetscErrorCode LOG_CELL_VERTICES(const Cell *cell, PetscMPIInt rank)
Prints the coordinates of a cell's vertices.
Definition logging.c:187
static void IntToStr(int value, char *buf, size_t bufsize)
Definition logging.c:250
PetscErrorCode ProfilingResetTimestepCounters(void)
Definition logging.c:1054
const char * MomentumSolverTypeToString(MomentumSolverType SolverFlag)
Helper function to convert Momentum Solver flag to a string representation.
Definition logging.c:706
static PetscErrorCode ComputeMaxColumnWidths(PetscInt nParticles, const PetscMPIInt *ranks, const PetscInt64 *pids, const PetscInt *cellIDs, const PetscReal *positions, const PetscReal *velocities, const PetscReal *weights, int *wRank, int *wPID, int *wCell, int *wPos, int *wVel, int *wWt)
Definition logging.c:299
double start_time
Definition logging.c:959
PetscErrorCode ProfilingLogTimestepSummary(PetscInt step)
Logs the performance summary for the current timestep and resets timers.
Definition logging.c:1064
const char * ParticleInitializationToString(ParticleInitializationType ParticleInitialization)
Helper function to convert ParticleInitialization to a string representation.
Definition logging.c:675
void _ProfilingStart(const char *func_name)
Definition logging.c:1032
const char * name
Definition logging.c:954
Logging utilities and macros for PETSc-based applications.
#define LOG_ALLOW_SYNC(scope, level, fmt,...)
----— DEBUG ---------------------------------------— #define LOG_ALLOW(scope, level,...
Definition logging.h:267
PetscBool log_to_console
Definition logging.h:58
#define LOCAL
Logging scope definitions for controlling message output.
Definition logging.h:45
#define GLOBAL
Scope for global logging across all processes.
Definition logging.h:46
#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:200
PetscReal bnorm
Definition logging.h:59
PetscInt step
Definition logging.h:60
LogLevel
Enumeration of logging levels.
Definition logging.h:27
@ LOG_ERROR
Critical errors that may halt the program.
Definition logging.h:28
@ LOG_PROFILE
Exclusive log level for performance timing and profiling.
Definition logging.h:30
@ LOG_TRACE
Very fine-grained tracing information for in-depth debugging.
Definition logging.h:33
@ LOG_INFO
Informational messages about program execution.
Definition logging.h:31
@ LOG_WARNING
Non-critical issues that warrant attention.
Definition logging.h:29
@ LOG_DEBUG
Detailed debugging information.
Definition logging.h:32
@ LOG_VERBOSE
Extremely detailed logs, typically for development use only.
Definition logging.h:34
FILE * file_handle
Definition logging.h:57
PetscInt block_id
Definition logging.h:61
Context for a dual-purpose KSP monitor.
Definition logging.h:56
LESModelType
Identifies the six logical faces of a structured computational block.
Definition variables.h:447
@ DYNAMIC_SMAGORINSKY
Definition variables.h:450
@ NO_LES_MODEL
Definition variables.h:448
@ CONSTANT_SMAGORINSKY
Definition variables.h:449
Vec lDiffusivityGradient
Definition variables.h:759
Vec lCent
Definition variables.h:776
BCType
Defines the general mathematical/physical Category of a boundary.
Definition variables.h:210
@ INLET
Definition variables.h:217
@ INTERFACE
Definition variables.h:212
@ FARFIELD
Definition variables.h:218
@ OUTLET
Definition variables.h:216
@ PERIODIC
Definition variables.h:219
@ WALL
Definition variables.h:213
PetscMPIInt rank
Definition variables.h:588
SimCtx * simCtx
Back-pointer to the master simulation context.
Definition variables.h:731
ParticleInitializationType
Enumerator to identify the particle initialization strategy.
Definition variables.h:465
@ PARTICLE_INIT_SURFACE_RANDOM
Random placement on the inlet face.
Definition variables.h:466
@ PARTICLE_INIT_SURFACE_EDGES
Deterministic placement at inlet face edges.
Definition variables.h:469
@ PARTICLE_INIT_POINT_SOURCE
All particles at a fixed (psrc_x,psrc_y,psrc_z) — for validation.
Definition variables.h:468
@ PARTICLE_INIT_VOLUME
Random volumetric distribution across the domain.
Definition variables.h:467
ParticleLocationStatus
Defines the state of a particle with respect to its location and migration status during the iterativ...
Definition variables.h:135
@ LOST
Definition variables.h:139
@ NEEDS_LOCATION
Definition variables.h:136
@ ACTIVE_AND_LOCATED
Definition variables.h:137
@ UNINITIALIZED
Definition variables.h:140
@ MIGRATING_OUT
Definition variables.h:138
PetscReal FluxOutSum
Definition variables.h:658
Vec Centz
Definition variables.h:777
PetscInt particlesLostLastStep
Definition variables.h:682
PetscInt KM
Definition variables.h:737
Vec lZet
Definition variables.h:776
BCHandlerType
Defines the specific computational "strategy" for a boundary handler.
Definition variables.h:230
@ BC_HANDLER_INLET_PULSATILE_FLUX
Definition variables.h:239
@ BC_HANDLER_PERIODIC_GEOMETRIC
Definition variables.h:243
@ BC_HANDLER_INLET_PARABOLIC
Definition variables.h:236
@ BC_HANDLER_INLET_CONSTANT_VELOCITY
Definition variables.h:235
@ BC_HANDLER_PERIODIC_DRIVEN_INITIAL_FLUX
Definition variables.h:246
@ BC_HANDLER_INTERFACE_OVERSET
Definition variables.h:244
@ BC_HANDLER_PERIODIC_DRIVEN_CONSTANT_FLUX
Definition variables.h:245
@ BC_HANDLER_WALL_MOVING
Definition variables.h:233
@ BC_HANDLER_WALL_NOSLIP
Definition variables.h:232
@ BC_HANDLER_OUTLET_CONSERVATION
Definition variables.h:241
@ BC_HANDLER_FARFIELD_NONREFLECTING
Definition variables.h:240
@ BC_HANDLER_OUTLET_PRESSURE
Definition variables.h:242
@ BC_HANDLER_SYMMETRY_PLANE
Definition variables.h:234
@ BC_HANDLER_UNDEFINED
Definition variables.h:231
Vec lCellFieldAtCorner
Definition variables.h:764
PetscInt _this
Definition variables.h:741
PetscInt occupiedCellCount
Definition variables.h:685
Vec lPsi
Definition variables.h:801
char ** criticalFuncs
Definition variables.h:710
Vec DiffusivityGradient
Definition variables.h:759
#define MAX_FILENAME_LENGTH
Definition variables.h:491
PetscInt StartStep
Definition variables.h:595
MomentumSolverType
Enumerator to identify the various momentumsolvers.
Definition variables.h:457
@ MOMENTUM_SOLVER_DUALTIME_NK_ARNOLDI
Definition variables.h:460
@ MOMENTUM_SOLVER_EXPLICIT_RK
Definition variables.h:458
@ MOMENTUM_SOLVER_DUALTIME_PICARD_RK4
Definition variables.h:459
@ MOMENTUM_SOLVER_DUALTIME_NK_ANALYTIC_JACOBIAN
Definition variables.h:461
PetscScalar x
Definition variables.h:101
PetscReal MaxDiv
Definition variables.h:705
Vec Centx
Definition variables.h:777
PetscInt MaxDivx
Definition variables.h:706
PetscInt MaxDivy
Definition variables.h:706
PetscInt MaxDivz
Definition variables.h:706
char log_dir[PETSC_MAX_PATH_LEN]
Definition variables.h:609
PetscInt MaxDivFlatArg
Definition variables.h:706
PetscReal FluxInSum
Definition variables.h:658
Vec lCsi
Definition variables.h:776
PetscScalar z
Definition variables.h:101
Vec Ucat
Definition variables.h:755
Vec ParticleCount
Definition variables.h:800
Vec CellFieldAtCorner
Definition variables.h:764
PetscInt JM
Definition variables.h:737
PetscInt particlesMigratedLastStep
Definition variables.h:684
Vec lUcont
Definition variables.h:755
PetscInt step
Definition variables.h:593
Vec Diffusivity
Definition variables.h:758
Vec lUcat
Definition variables.h:755
PetscInt migrationPassesLastStep
Definition variables.h:683
PetscScalar y
Definition variables.h:101
@ EXEC_MODE_SOLVER
Definition variables.h:558
@ EXEC_MODE_POSTPROCESSOR
Definition variables.h:559
PetscInt IM
Definition variables.h:737
@ TOP
Definition variables.h:145
@ FRONT
Definition variables.h:145
@ BOTTOM
Definition variables.h:145
@ BACK
Definition variables.h:145
@ LEFT
Definition variables.h:145
@ RIGHT
Definition variables.h:145
Vec lEta
Definition variables.h:776
Vec lDiffusivity
Definition variables.h:758
Vec Centy
Definition variables.h:777
ExecutionMode exec_mode
Definition variables.h:603
PetscInt nCriticalFuncs
Definition variables.h:711
PetscReal summationRHS
Definition variables.h:704
Cmpnts vertices[8]
Coordinates of the eight vertices of the cell.
Definition variables.h:161
Vec Psi
Definition variables.h:801
PetscReal particleLoadImbalance
Definition variables.h:686
BCFace
Identifies the six logical faces of a structured computational block.
Definition variables.h:203
@ BC_FACE_NEG_X
Definition variables.h:204
@ BC_FACE_POS_Z
Definition variables.h:206
@ BC_FACE_POS_Y
Definition variables.h:205
@ BC_FACE_NEG_Z
Definition variables.h:206
@ BC_FACE_POS_X
Definition variables.h:204
@ BC_FACE_NEG_Y
Definition variables.h:205
Defines the vertices of a single hexahedral grid cell.
Definition variables.h:160
A 3D point or vector with PetscScalar components.
Definition variables.h:100
The master context for the entire simulation.
Definition variables.h:585
User-defined context containing data specific to a single computational grid level.
Definition variables.h:728