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] PetscInt The ParticleInitialization value.
673 * @return Pointer to a constant string representing the FieldInitialization.
674 */
675const char* ParticleInitializationToString(PetscInt ParticleInitialization)
676{
677 switch(ParticleInitialization){
678 case 0: return "Surface: Random";
679 case 1: return "Volume";
680 case 3: return "Surface: At edges";
681 default: return "Unknown Particle Initialization";
682 }
683}
684
685/**
686 * @brief Helper function to convert LES Flag to a string representation.
687 * @param[in] PetscInt The LES flag value.
688 * @return Pointer to a constant string representing the FieldInitialization.
689 */
690const char* LESFlagToString(PetscInt LESFlag)
691{
692 switch(LESFlag){
693 case 0: return "No LES";
694 case 1: return "Constant Smagorinsky";
695 case 2: return "Dynamic Smagorinsky";
696 default: return "Unknown LES Flag";
697 }
698}
699
700/**
701 * @brief Helper function to convert BCType enum to a string representation.
702 * @param[in] type The BCType enum value.
703 * @return Pointer to a constant string representing the BC type.
704 */
705const char* BCTypeToString(BCType type) {
706 switch (type) {
707 // case DIRICHLET: return "DIRICHLET";
708 // case NEUMANN: return "NEUMANN";
709 case WALL: return "WALL";
710 case INLET: return "INLET";
711 case OUTLET: return "OUTLET";
712 case FARFIELD: return "FARFIELD";
713 case PERIODIC: return "PERIODIC";
714 case INTERFACE: return "INTERFACE";
715 case NOGRAD: return "NOGRAD";
716
717 // case CUSTOM: return "CUSTOM";
718 default: return "Unknown BC Type";
719 }
720}
721
722/**
723 * @brief Converts a BCHandlerType enum to its string representation.
724 *
725 * Provides a descriptive string for a specific boundary condition implementation strategy.
726 * This is crucial for logging the exact behavior configured for a face.
727 *
728 * @param handler_type The BCHandlerType enum value (e.g., BC_HANDLER_WALL_NOSLIP).
729 * @return A constant character string corresponding to the enum. Returns
730 * "UNKNOWN_HANDLER" if the enum value is not recognized.
731 */
732const char* BCHandlerTypeToString(BCHandlerType handler_type) {
733 switch (handler_type) {
734 // Wall & Symmetry Handlers
735 case BC_HANDLER_WALL_NOSLIP: return "noslip";
736 case BC_HANDLER_WALL_MOVING: return "moving";
737 case BC_HANDLER_SYMMETRY_PLANE: return "symmetry_plane";
738
739 // Inlet Handlers
740 case BC_HANDLER_INLET_CONSTANT_VELOCITY: return "constant_velocity";
741 case BC_HANDLER_INLET_PULSATILE_FLUX: return "pulsatile_flux";
742 case BC_HANDLER_INLET_PARABOLIC: return "parabolic";
743
744 // Outlet Handlers
745 case BC_HANDLER_OUTLET_CONSERVATION: return "conservation";
746 case BC_HANDLER_OUTLET_PRESSURE: return "pressure";
747
748 // Other Physical Handlers
749 case BC_HANDLER_FARFIELD_NONREFLECTING: return "nonreflecting";
750 case BC_HANDLER_NOGRAD_COPY_GHOST: return "no_gradient";
751
752 // Multi-Block / Interface Handlers
753 case BC_HANDLER_PERIODIC: return "periodic";
754 case BC_HANDLER_INTERFACE_OVERSET: return "overset";
755
756 // Default case
758 default: return "UNKNOWN_HANDLER";
759 }
760}
761
762/**
763 * @brief Destroys the DualMonitorCtx.
764 *
765 * This function is passed to KSPMonitorSet to ensure the viewer is
766 * properly destroyed and the context memory is freed when the KSP is destroyed.
767 * @param Ctx a pointer to the context pointer to be destroyed
768 * @return PetscErrorCode
769 */
770PetscErrorCode DualMonitorDestroy(void **ctx)
771{
772 DualMonitorCtx *monctx = (DualMonitorCtx*)*ctx;
773 PetscErrorCode ierr;
774 PetscMPIInt rank;
775
776 PetscFunctionBeginUser;
777 ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank); CHKERRQ(ierr);
778 if(!rank && monctx->file_handle){
779 fclose(monctx->file_handle);
780 }
781
782 ierr = PetscFree(monctx); CHKERRQ(ierr);
783 *ctx = NULL;
784 PetscFunctionReturn(0);
785}
786
787/**
788 * @brief A custom KSP monitor that logs the true residual to a file and optionally to the console.
789 *
790 * This function replicates the behavior of KSPMonitorTrueResidualNorm by calculating
791 * the true residual norm ||b - Ax|| itself. It unconditionally logs to a file
792 * viewer and conditionally logs to the console based on a flag in the context.
793 *
794 * @param ksp The Krylov subspace context.
795 * @param it The current iteration number.
796 * @param rnorm The preconditioned residual norm (ignored, we compute our own).
797 * @param ctx A pointer to the DualMonitorCtx structure.
798 * @return PetscErrorCode 0 on success.
799 */
800#undef __FUNCT__
801#define __FUNCT__ "DualKSPMonitor"
802PetscErrorCode DualKSPMonitor(KSP ksp, PetscInt it, PetscReal rnorm, void *ctx)
803{
804 DualMonitorCtx *monctx = (DualMonitorCtx*)ctx;
805 PetscErrorCode ierr;
806 PetscReal trnorm, relnorm;
807 Vec r;
808 char norm_buf[256];
809 PetscMPIInt rank;
810
811 PetscFunctionBeginUser;
812 ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank); CHKERRQ(ierr);
813
814 // 1. Calculate the true residual norm.
815 ierr = KSPBuildResidual(ksp, NULL, NULL, &r); CHKERRQ(ierr);
816 ierr = VecNorm(r, NORM_2, &trnorm); CHKERRQ(ierr);
817
818 // 2. On the first iteration, compute and store the norm of the RHS vector `b`.
819 if (it == 0) {
820 Vec b;
821 ierr = KSPGetRhs(ksp, &b); CHKERRQ(ierr);
822 ierr = VecNorm(b, NORM_2, &monctx->bnorm); CHKERRQ(ierr);
823 }
824
825 if(!rank){
826 // 3. Compute the relative norm and format the output string.
827 if (monctx->bnorm > 1.e-15) {
828 relnorm = trnorm / monctx->bnorm;
829 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);
830 } else {
831 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);
832 }
833
834 // 4. Log to the file viewer (unconditionally).
835 if(monctx->file_handle){
836 ierr = PetscFPrintf(PETSC_COMM_SELF,monctx->file_handle,"%s\n", norm_buf); CHKERRQ(ierr);
837 }
838 // 5. Log to the console (conditionally).
839 if (monctx->log_to_console) {
840 PetscFPrintf(PETSC_COMM_SELF,stdout, "%s\n", norm_buf); CHKERRQ(ierr);
841 }
842
843 } //rank
844
845 PetscFunctionReturn(0);
846}
847
848/**
849 * @brief Logs continuity metrics for a single block to a file.
850 *
851 * This function should be called for each block, once per timestep. It opens a
852 * central log file in append mode. To ensure the header is written only once,
853 * it checks if it is processing block 0 on the simulation's start step.
854 *
855 * @param user A pointer to the UserCtx for the specific block whose metrics
856 * are to be logged. The function accesses both global (SimCtx)
857 * and local (user->...) data.
858 * @return PetscErrorCode 0 on success.
859 */
860#undef __FUNCT__
861#define __FUNCT__ "LOG_CONTINUITY_METRICS"
862PetscErrorCode LOG_CONTINUITY_METRICS(UserCtx *user)
863{
864 PetscErrorCode ierr;
865 PetscMPIInt rank;
866 SimCtx *simCtx = user->simCtx; // Get the shared SimCtx
867 const PetscInt bi = user->_this; // Get this block's specific ID
868 const PetscInt ti = simCtx->step; // Get the current timestep
869
870 PetscFunctionBeginUser;
871 ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank); CHKERRQ(ierr);
872
873 // Only rank 0 performs file I/O.
874 if (!rank) {
875 FILE *f;
876 char filen[128];
877 sprintf(filen, "%s/Continuity_Metrics.log",simCtx->log_dir);
878
879 // Open the log file in append mode.
880 f = fopen(filen, "a");
881 if (!f) {
882 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open log file: %s", filen);
883 }
884
885 // Write a header only for the very first block (bi=0) on the very
886 // first timestep (ti=StartStep + 1). This ensures it's written only once.
887 if (ti == simCtx->StartStep + 1 && bi == 0) {
888 PetscFPrintf(PETSC_COMM_SELF, f, "%-10s | %-6s | %-18s | %-30s | %-18s | %-18s | %-18s | %-18s\n",
889 "Timestep", "Block", "Max Divergence", "Max Divergence Location ([k][j][i]=idx)", "Sum(RHS)","Total Flux In", "Total Flux Out", "Net Flux");
890 PetscFPrintf(PETSC_COMM_SELF, f, "------------------------------------------------------------------------------------------------------------------------------------------\n");
891 }
892
893 // Prepare the data strings and values for the current block.
894 PetscReal net_flux = simCtx->FluxInSum - simCtx->FluxOutSum;
895 char location_str[64];
896 sprintf(location_str, "([%d][%d][%d] = %d)", (int)simCtx->MaxDivz, (int)simCtx->MaxDivy, (int)simCtx->MaxDivx, (int)simCtx->MaxDivFlatArg);
897
898 // Write the formatted line for the current block.
899 PetscFPrintf(PETSC_COMM_SELF, f, "%-10d | %-6d | %-18.10e | %-39s | %-18.10e | %-18.10e | %-18.10e | %-18.10e\n",
900 (int)ti,
901 (int)bi,
902 (double)simCtx->MaxDiv,
903 location_str,
904 (double)simCtx->summationRHS,
905 (double)simCtx->FluxInSum,
906 (double)simCtx->FluxOutSum,
907 (double)net_flux);
908
909 fclose(f);
910 }
911
912 PetscFunctionReturn(0);
913}
914
915/**
916 * @brief A function that outputs the name of the current level in the ParticleLocation enum.
917 * @param level The ParticleLocation enum value.
918 * @return A constant character string corresponding to the enum. Returns
919 * "UNKNOWN_LEVEL" if the enum value is not recognized.
920 */
922{
923 switch (level) {
924 case NEEDS_LOCATION: return "NEEDS_LOCATION";
925 case ACTIVE_AND_LOCATED: return "ACTIVE_AND_LOCATED";
926 case MIGRATING_OUT: return "MIGRATING_OUT";
927 case LOST: return "LOST";
928 case UNINITIALIZED: return "UNINITIALIZED";
929 default: return "UNKNOWN_LEVEL";
930 }
931}
932
933///////// Profiling System /////////
934
935// Data structure to hold profiling info for one function
936typedef struct {
937 const char *name;
942 double start_time; // Timer for the current call
943 PetscBool always_log;
945
946// Global registry for all profiled functions
948static PetscInt g_profiler_count = 0;
949static PetscInt g_profiler_capacity = 0;
950
951// Internal helper to find a function in the registry or create it
952static PetscErrorCode _FindOrCreateEntry(const char *func_name, PetscInt *idx)
953{
954 PetscFunctionBeginUser;
955 // Search for existing entry
956 for (PetscInt i = 0; i < g_profiler_count; ++i) {
957 if (strcmp(g_profiler_registry[i].name, func_name) == 0) {
958 *idx = i;
959 PetscFunctionReturn(0);
960 }
961 }
962
963 // Not found, create a new entry
965 PetscInt new_capacity = g_profiler_capacity == 0 ? 16 : g_profiler_capacity * 2;
966 PetscErrorCode ierr = PetscRealloc(sizeof(ProfiledFunction) * new_capacity, &g_profiler_registry); CHKERRQ(ierr);
967 g_profiler_capacity = new_capacity;
968 }
969
970 *idx = g_profiler_count;
971 g_profiler_registry[*idx].name = func_name;
977 g_profiler_registry[*idx].always_log = PETSC_FALSE;
979
980 PetscFunctionReturn(0);
981}
982
983// --- Public API Implementation ---
984/**
985 * @brief Initializes the custom profiling system using configuration from SimCtx.
986 *
987 * This function sets up the internal data structures for tracking function
988 * performance. It reads the list of "critical functions" from the provided
989 * SimCtx and marks them for per-step logging at LOG_INFO level.
990 *
991 * It should be called once at the beginning of the application, after
992 * CreateSimulationContext() but before the main time loop.
993 *
994 * @param simCtx The master simulation context, which contains the list of
995 * critical function names to always log.
996 * @return PetscErrorCode
997 */
998PetscErrorCode ProfilingInitialize(SimCtx *simCtx)
999{
1000 PetscFunctionBeginUser;
1001 if (!simCtx) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "SimCtx cannot be null for ProfilingInitialize");
1002
1003 // Iterate through the list of critical functions provided in SimCtx
1004 for (PetscInt i = 0; i < simCtx->nCriticalFuncs; ++i) {
1005 PetscInt idx;
1006 const char *func_name = simCtx->criticalFuncs[i];
1007 PetscErrorCode ierr = _FindOrCreateEntry(func_name, &idx); CHKERRQ(ierr);
1008 g_profiler_registry[idx].always_log = PETSC_TRUE;
1009
1010 LOG_ALLOW(GLOBAL, LOG_DEBUG, "Marked '%s' as a critical function for profiling.\n", func_name);
1011 }
1012 PetscFunctionReturn(0);
1013}
1014
1015void _ProfilingStart(const char *func_name)
1016{
1017 PetscInt idx;
1018 if (_FindOrCreateEntry(func_name, &idx) != 0) return; // Fail silently
1019 PetscTime(&g_profiler_registry[idx].start_time);
1020}
1021
1022void _ProfilingEnd(const char *func_name)
1023{
1024 double end_time;
1025 PetscTime(&end_time);
1026
1027 PetscInt idx;
1028 if (_FindOrCreateEntry(func_name, &idx) != 0) return; // Fail silently
1029
1030 double elapsed = end_time - g_profiler_registry[idx].start_time;
1031 g_profiler_registry[idx].total_time += elapsed;
1032 g_profiler_registry[idx].current_step_time += elapsed;
1035}
1036
1038{
1039 PetscFunctionBeginUser;
1040 for (PetscInt i = 0; i < g_profiler_count; ++i) {
1043 }
1044 PetscFunctionReturn(0);
1045}
1046
1047PetscErrorCode ProfilingLogTimestepSummary(PetscInt step)
1048{
1049 LogLevel log_level = get_log_level();
1050 PetscBool should_print = PETSC_FALSE;
1051
1052 PetscFunctionBeginUser;
1053
1054 // Decide if we should print anything at all
1055 if (log_level >= LOG_PROFILE) {
1056 for (PetscInt i = 0; i < g_profiler_count; ++i) {
1057 if (g_profiler_registry[i].current_step_call_count > 0) {
1058 if (log_level == LOG_PROFILE || g_profiler_registry[i].always_log) {
1059 should_print = PETSC_TRUE;
1060 break;
1061 }
1062 }
1063 }
1064 }
1065
1066 if (should_print) {
1067 PetscPrintf(PETSC_COMM_SELF, "[PROFILE] ----- Timestep %d Summary -----\n", step);
1068 for (PetscInt i = 0; i < g_profiler_count; ++i) {
1069 if (g_profiler_registry[i].current_step_call_count > 0) {
1070 if (log_level == LOG_PROFILE || g_profiler_registry[i].always_log) {
1071 PetscPrintf(PETSC_COMM_SELF, "[PROFILE] %-25s: %.6f s (%lld calls)\n",
1072 g_profiler_registry[i].name,
1073 g_profiler_registry[i].current_step_time,
1074 g_profiler_registry[i].current_step_call_count);
1075 }
1076 }
1077 }
1078 }
1079
1080 // Reset per-step counters for the next iteration
1081 for (PetscInt i = 0; i < g_profiler_count; ++i) {
1084 }
1085 PetscFunctionReturn(0);
1086}
1087
1088// Comparison function for qsort to sort by total_time in descending order
1089static int _CompareProfiledFunctions(const void *a, const void *b)
1090{
1091 const ProfiledFunction *func_a = (const ProfiledFunction *)a;
1092 const ProfiledFunction *func_b = (const ProfiledFunction *)b;
1093
1094 if (func_a->total_time < func_b->total_time) return 1;
1095 if (func_a->total_time > func_b->total_time) return -1;
1096 return 0;
1097}
1098
1099/**
1100 * @brief the profiling excercise and build a profiling summary which is then printed to a log file.
1101 *
1102 * @param simCtx The Simulation Context Structure that can contains all the data regarding the simulation.
1103 *
1104 * @return PetscErrorCode 0 on success.
1105 */
1106PetscErrorCode ProfilingFinalize(SimCtx *simCtx)
1107{
1108 PetscInt rank = simCtx->rank;
1109 PetscFunctionBeginUser;
1110 if (!rank) {
1111
1112 //--- Step 0: Create a file viewer for log file
1113 FILE *f;
1114 char filen[128];
1115 sprintf(filen, "%s/ProfilingSummary.log",simCtx->log_dir);
1116
1117 // Open the log file in write mode.
1118 f = fopen(filen,"w");
1119 if(!f){
1120 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot Open log file: %s",filen);
1121 }
1122
1123 // --- Step 1: Sort the data for readability ---
1125
1126 // --- Step 2: Dynamically determine the width for the function name column ---
1127 PetscInt max_name_len = strlen("Function"); // Start with the header's length
1128 for (PetscInt i = 0; i < g_profiler_count; ++i) {
1129 if (g_profiler_registry[i].total_call_count > 0) {
1130 PetscInt len = strlen(g_profiler_registry[i].name);
1131 if (len > max_name_len) {
1132 max_name_len = len;
1133 }
1134 }
1135 }
1136 // Add a little padding
1137 max_name_len += 2;
1138
1139 // --- Step 3: Define fixed widths for numeric columns for consistent alignment ---
1140 const int time_width = 18;
1141 const int count_width = 15;
1142 const int avg_width = 22;
1143
1144 // --- Step 4: Print the formatted table ---
1145 PetscFPrintf(PETSC_COMM_SELF, f, "=================================================================================================================\n");
1146 PetscFPrintf(PETSC_COMM_SELF, f, " FINAL PROFILING SUMMARY (Sorted by Total Time)\n");
1147 PetscFPrintf(PETSC_COMM_SELF, f, "=================================================================================================================\n");
1148
1149 // Header Row
1150 PetscFPrintf(PETSC_COMM_SELF, f, "%-*s | %-*s | %-*s | %-*s\n",
1151 max_name_len, "Function",
1152 time_width, "Total Time (s)",
1153 count_width, "Call Count",
1154 avg_width, "Avg. Time/Call (ms)");
1155
1156 // Separator Line (dynamically sized)
1157 for (int i = 0; i < max_name_len; i++) PetscFPrintf(PETSC_COMM_SELF, f, "-");
1158 PetscFPrintf(PETSC_COMM_SELF, f, "-|-");
1159 for (int i = 0; i < time_width; i++) PetscFPrintf(PETSC_COMM_SELF, f, "-");
1160 PetscFPrintf(PETSC_COMM_SELF, f, "-|-");
1161 for (int i = 0; i < count_width; i++) PetscFPrintf(PETSC_COMM_SELF, f, "-");
1162 PetscFPrintf(PETSC_COMM_SELF, f, "-|-");
1163 for (int i = 0; i < avg_width; i++) PetscFPrintf(PETSC_COMM_SELF, f, "-");
1164 PetscFPrintf(PETSC_COMM_SELF, f, "\n");
1165
1166 // Data Rows
1167 for (PetscInt i = 0; i < g_profiler_count; ++i) {
1168 if (g_profiler_registry[i].total_call_count > 0) {
1169 double avg_time_ms = (g_profiler_registry[i].total_time / g_profiler_registry[i].total_call_count) * 1000.0;
1170 PetscFPrintf(PETSC_COMM_SELF, f, "%-*s | %*.*f | %*lld | %*.*f\n",
1171 max_name_len, g_profiler_registry[i].name,
1172 time_width, 6, g_profiler_registry[i].total_time,
1173 count_width, g_profiler_registry[i].total_call_count,
1174 avg_width, 6, avg_time_ms);
1175 PetscFPrintf(PETSC_COMM_SELF, f, "------------------------------------------------------------------------------------------------------------------\n");
1176 }
1177 }
1178 PetscFPrintf(PETSC_COMM_SELF, f, "==================================================================================================================\n");
1179
1180 fclose(f);
1181 }
1182
1183 // --- Final Cleanup ---
1184 PetscFree(g_profiler_registry);
1185 g_profiler_registry = NULL;
1186 g_profiler_count = 0;
1188 PetscFunctionReturn(0);
1189}
1190
1191/*================================================================================*
1192 * PROGRESS BAR UTILITY *
1193 *================================================================================*/
1194
1195/**
1196 * @brief Prints a progress bar to the console.
1197 *
1198 * This function should only be called by the root process (rank 0). It uses
1199 * a carriage return `\r` to overwrite the same line in the terminal, creating
1200 * a dynamic progress bar.
1201 *
1202 * @param step The current step index from the loop (e.g., from 0 to N-1).
1203 * @param startStep The global starting step number of the simulation.
1204 * @param totalSteps The total number of steps to be run in this simulation instance.
1205 * @param currentTime The current simulation time to display.
1206 */
1207void PrintProgressBar(PetscInt step, PetscInt startStep, PetscInt totalSteps, PetscReal currentTime)
1208{
1209 if (totalSteps <= 0) return;
1210
1211 // --- Configuration ---
1212 const int barWidth = 50;
1213
1214 // --- Calculation ---
1215 // Calculate progress as a fraction from 0.0 to 1.0
1216 PetscReal progress = (PetscReal)(step - startStep + 1) / totalSteps;
1217 // Ensure progress doesn't exceed 1.0 due to floating point inaccuracies
1218 if (progress > 1.0) progress = 1.0;
1219
1220 int pos = (int)(barWidth * progress);
1221
1222 // --- Printing ---
1223 // Carriage return moves cursor to the beginning of the line
1224 PetscPrintf(PETSC_COMM_SELF, "\rProgress: [");
1225
1226 for (int i = 0; i < barWidth; ++i) {
1227 if (i < pos) {
1228 PetscPrintf(PETSC_COMM_SELF, "=");
1229 } else if (i == pos) {
1230 PetscPrintf(PETSC_COMM_SELF, ">");
1231 } else {
1232 PetscPrintf(PETSC_COMM_SELF, " ");
1233 }
1234 }
1235
1236 // Print percentage, step count, and current time
1237 PetscPrintf(PETSC_COMM_SELF, "] %3d%% (Step %ld/%ld, t=%.4f)",
1238 (int)(progress * 100.0),
1239 step + 1,
1240 startStep + totalSteps,
1241 currentTime);
1242
1243 // Flush the output buffer to ensure the bar is displayed immediately
1244 fflush(stdout);
1245}
1246
1247#undef __FUNCT__
1248#define __FUNCT__ "LOG_FIELD_MIN_MAX"
1249/**
1250 * @brief Computes and logs the local and global min/max values of a specified field,
1251 * respecting the solver's data layout architecture.
1252 *
1253 * This utility function inspects a PETSc Vec and calculates the minimum and maximum
1254 * values for its components, both locally and globally.
1255 *
1256 * It is "architecture-aware" and adjusts its iteration range to only include
1257 * physically meaningful data points:
1258 *
1259 * - **Cell-Centered Fields ("Ucat", "P"):** It uses the "Shifted Index Architecture,"
1260 * iterating from index 1 to N-1 to exclude the ghost/tool values at indices 0 and N.
1261 * - **Node-Centered Fields ("Coordinates"):** It iterates from index 0 to N-1, covering all
1262 * physical nodes.
1263 * - **Face-Centered Fields ("Ucont"):** It iterates from index 0 to N-1, covering all
1264 * physical faces.
1265 *
1266 * The results are printed to the standard output in a formatted, easy-to-read table.
1267 *
1268 * @param[in] user Pointer to the user-defined context.
1269 * @param[in] fieldName A string descriptor for the field being analyzed.
1270 *
1271 * @return PetscErrorCode Returns 0 on success, non-zero on failure.
1272 */
1273PetscErrorCode LOG_FIELD_MIN_MAX(UserCtx *user, const char *fieldName)
1274{
1275 PetscErrorCode ierr;
1276 PetscInt i, j, k;
1277 DMDALocalInfo info;
1278
1279 Vec fieldVec = NULL;
1280 DM dm = NULL;
1281 PetscInt dof;
1282 char data_layout[20];
1283
1284 PetscFunctionBeginUser;
1285
1286 // --- 1. Map string name to PETSc objects and determine data layout ---
1287 if (strcasecmp(fieldName, "Ucat") == 0) {
1288 fieldVec = user->Ucat; dm = user->fda; dof = 3; strcpy(data_layout, "Cell-Centered");
1289 } else if (strcasecmp(fieldName, "P") == 0) {
1290 fieldVec = user->P; dm = user->da; dof = 1; strcpy(data_layout, "Cell-Centered");
1291 } else if (strcasecmp(fieldName, "Ucont") == 0) {
1292 fieldVec = user->lUcont; dm = user->fda; dof = 3; strcpy(data_layout, "Face-Centered");
1293 } else if (strcasecmp(fieldName, "Coordinates") == 0) {
1294 ierr = DMGetCoordinates(user->da, &fieldVec); CHKERRQ(ierr);
1295 dm = user->fda; dof = 3; strcpy(data_layout, "Node-Centered");
1296 } else if (strcasecmp(fieldName,"Psi") == 0) {
1297 fieldVec = user->Psi; dm = user->da; dof = 1; strcpy(data_layout, "Cell-Centered"); // Assuming Psi is cell-centered
1298 } else {
1299 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Field %s not recognized.", fieldName);
1300 }
1301
1302 if (!fieldVec) {
1303 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Vector for field '%s' is NULL.", fieldName);
1304 }
1305 if (!dm) {
1306 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "DM for field '%s' is NULL.", fieldName);
1307 }
1308
1309 ierr = DMDAGetLocalInfo(dm, &info); CHKERRQ(ierr);
1310
1311 // --- 2. Define Architecture-Aware Loop Bounds ---
1312 PetscInt i_start, i_end, j_start, j_end, k_start, k_end;
1313
1314 if (strcmp(data_layout, "Cell-Centered") == 0) {
1315 // For cell-centered data, the physical values are stored from index 1 to N-1.
1316 // We find the intersection of the rank's owned range [xs, xe) with the
1317 // physical data range [1, IM-1).
1318 i_start = PetscMax(info.xs, 1); i_end = PetscMin(info.xs + info.xm, user->IM);
1319 j_start = PetscMax(info.ys, 1); j_end = PetscMin(info.ys + info.ym, user->JM);
1320 k_start = PetscMax(info.zs, 1); k_end = PetscMin(info.zs + info.zm, user->KM);
1321 } else { // For Node- or Face-Centered data
1322 // The physical values are stored from index 0 to N-1.
1323 // We find the intersection of the rank's owned range [xs, xe) with the
1324 // physical data range [0, IM-1].
1325 i_start = PetscMax(info.xs, 0); i_end = PetscMin(info.xs + info.xm, user->IM);
1326 j_start = PetscMax(info.ys, 0); j_end = PetscMin(info.ys + info.ym, user->JM);
1327 k_start = PetscMax(info.zs, 0); k_end = PetscMin(info.zs + info.zm, user->KM);
1328 }
1329
1330 // --- 3. Barrier for clean, grouped output ---
1331 ierr = MPI_Barrier(PETSC_COMM_WORLD); CHKERRQ(ierr);
1332 if (user->simCtx->rank == 0) {
1333 PetscPrintf(PETSC_COMM_SELF, "\n--- Field Ranges: [%s] (Layout: %s) ---\n", fieldName, data_layout);
1334 }
1335
1336 // --- 4. Branch on DoF and perform calculation with correct bounds ---
1337 if (dof == 1) {
1338 PetscReal localMin = PETSC_MAX_REAL, localMax = PETSC_MIN_REAL;
1339 PetscReal globalMin, globalMax;
1340 const PetscScalar ***array;
1341
1342 ierr = DMDAVecGetArrayRead(dm, fieldVec, &array); CHKERRQ(ierr);
1343 for (k = k_start; k < k_end; k++) {
1344 for (j = j_start; j < j_end; j++) {
1345 for (i = i_start; i < i_end; i++) {
1346 localMin = PetscMin(localMin, array[k][j][i]);
1347 localMax = PetscMax(localMax, array[k][j][i]);
1348 }
1349 }
1350 }
1351 ierr = DMDAVecRestoreArrayRead(dm, fieldVec, &array); CHKERRQ(ierr);
1352
1353 ierr = MPI_Allreduce(&localMin, &globalMin, 1, MPIU_REAL, MPI_MIN, PETSC_COMM_WORLD); CHKERRQ(ierr);
1354 ierr = MPI_Allreduce(&localMax, &globalMax, 1, MPIU_REAL, MPI_MAX, PETSC_COMM_WORLD); CHKERRQ(ierr);
1355
1356 PetscSynchronizedPrintf(PETSC_COMM_WORLD, " [Rank %d] Local Range: [ %11.4e , %11.4e ]\n", user->simCtx->rank, localMin, localMax);
1357 ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT); CHKERRQ(ierr);
1358 if (user->simCtx->rank == 0) {
1359 PetscPrintf(PETSC_COMM_SELF, " Global Range: [ %11.4e , %11.4e ]\n", globalMin, globalMax);
1360 }
1361
1362 } else if (dof == 3) {
1363 Cmpnts localMin = {PETSC_MAX_REAL, PETSC_MAX_REAL, PETSC_MAX_REAL};
1364 Cmpnts localMax = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL};
1365 Cmpnts globalMin, globalMax;
1366 const Cmpnts ***array;
1367
1368 ierr = DMDAVecGetArrayRead(dm, fieldVec, &array); CHKERRQ(ierr);
1369 for (k = k_start; k < k_end; k++) {
1370 for (j = j_start; j < j_end; j++) {
1371 for (i = i_start; i < i_end; i++) {
1372 localMin.x = PetscMin(localMin.x, array[k][j][i].x);
1373 localMin.y = PetscMin(localMin.y, array[k][j][i].y);
1374 localMin.z = PetscMin(localMin.z, array[k][j][i].z);
1375 localMax.x = PetscMax(localMax.x, array[k][j][i].x);
1376 localMax.y = PetscMax(localMax.y, array[k][j][i].y);
1377 localMax.z = PetscMax(localMax.z, array[k][j][i].z);
1378 }
1379 }
1380 }
1381 ierr = DMDAVecRestoreArrayRead(dm, fieldVec, &array); CHKERRQ(ierr);
1382
1383 ierr = MPI_Allreduce(&localMin, &globalMin, 3, MPIU_REAL, MPI_MIN, PETSC_COMM_WORLD); CHKERRQ(ierr);
1384 ierr = MPI_Allreduce(&localMax, &globalMax, 3, MPIU_REAL, MPI_MAX, PETSC_COMM_WORLD); CHKERRQ(ierr);
1385
1386 ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, " [Rank %d] Local X-Range: [ %11.4e , %11.4e ]\n", user->simCtx->rank, localMin.x, localMax.x);
1387 ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, " [Rank %d] Local Y-Range: [ %11.4e , %11.4e ]\n", user->simCtx->rank, localMin.y, localMax.y);
1388 ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, " [Rank %d] Local Z-Range: [ %11.4e , %11.4e ]\n", user->simCtx->rank, localMin.z, localMax.z);
1389 ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT); CHKERRQ(ierr);
1390
1391 if (user->simCtx->rank == 0) {
1392 PetscPrintf(PETSC_COMM_SELF, " [Global] X-Range: [ %11.4e , %11.4e ]\n", globalMin.x, globalMax.x);
1393 PetscPrintf(PETSC_COMM_SELF, " [Global] Y-Range: [ %11.4e , %11.4e ]\n", globalMin.y, globalMax.y);
1394 PetscPrintf(PETSC_COMM_SELF, " [Global] Z-Range: [ %11.4e , %11.4e ]\n", globalMin.z, globalMax.z);
1395 }
1396
1397 } else {
1398 SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "LogFieldStatistics only supports fields with 1 or 3 components, but field '%s' has %D.", fieldName, dof);
1399 }
1400
1401 // --- 5. Final barrier for clean output ordering ---
1402 ierr = MPI_Barrier(PETSC_COMM_WORLD); CHKERRQ(ierr);
1403 if (user->simCtx->rank == 0) {
1404 PetscPrintf(PETSC_COMM_SELF, "--------------------------------------------\n\n");
1405 }
1406
1407 PetscFunctionReturn(0);
1408}
1409
1410#undef __FUNCT__
1411#define __FUNCT__ "LOG_FIELD_ANATOMY"
1412/**
1413 * @brief Logs the anatomy of a specified field at key boundary locations,
1414 * respecting the solver's specific grid and variable architecture.
1415 *
1416 * This intelligent diagnostic function inspects a PETSc Vec and prints its values
1417 * at critical boundary locations (-Xi/+Xi, -Eta/+Eta, -Zeta/+Zeta). It is "architecture-aware":
1418 *
1419 * - **Cell-Centered Fields ("Ucat", "P"):** It correctly applies the "Shifted Index Architecture,"
1420 * where the value for geometric `Cell i` is stored at array index `i+1`. It labels
1421 * the output to clearly distinguish between true physical values and ghost values.
1422 * - **Face-Centered Fields ("Ucont"):** It uses a direct index mapping, where the value for
1423 * the face at `Node i` is stored at index `i`.
1424 * - **Node-Centered Fields ("Coordinates"):** It uses a direct index mapping, where the value for
1425 * `Node i` is stored at index `i`.
1426 *
1427 * The output is synchronized across MPI ranks to ensure readability and focuses on a
1428 * slice through the center of the domain to be concise.
1429 *
1430 * @param user A pointer to the UserCtx structure containing the DMs and Vecs.
1431 * @param field_name A string identifier for the field to log (e.g., "Ucat", "P", "Ucont", "Coordinates").
1432 * @param stage_name A string identifier for the current simulation stage (e.g., "After Advection").
1433 * @return PetscErrorCode Returns 0 on success, non-zero on failure.
1434 */
1435PetscErrorCode LOG_FIELD_ANATOMY(UserCtx *user, const char *field_name, const char *stage_name)
1436{
1437 PetscErrorCode ierr;
1438 DMDALocalInfo info;
1439 PetscMPIInt rank;
1440
1441 Vec vec_local = NULL;
1442 DM dm = NULL;
1443 PetscInt dof = 0;
1444 char data_layout[20]; // To store "Cell-Centered", "Face-Centered", etc.
1445
1446 PetscFunctionBeginUser;
1447 ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank); CHKERRQ(ierr);
1448
1449 // --- 1. Map string name to PETSc objects and determine data layout ---
1450 if (strcasecmp(field_name, "Ucat") == 0) {
1451 vec_local = user->lUcat; dm = user->fda; dof = 3; strcpy(data_layout, "Cell-Centered");
1452 } else if (strcasecmp(field_name, "P") == 0) {
1453 vec_local = user->lP; dm = user->da; dof = 1; strcpy(data_layout, "Cell-Centered");
1454 } else if (strcasecmp(field_name, "Psi") == 0) {
1455 vec_local = user->lPsi; dm = user->da; dof = 1; strcpy(data_layout, "Cell-Centered");
1456 } else if (strcasecmp(field_name, "Ucont") == 0) {
1457 vec_local = user->lUcont; dm = user->fda; dof = 3; strcpy(data_layout, "Face-Centered");
1458 } else if (strcasecmp(field_name, "Coordinates") == 0) {
1459 ierr = DMGetCoordinatesLocal(user->da, &vec_local); CHKERRQ(ierr);
1460 dm = user->fda; dof = 3; strcpy(data_layout, "Node-Centered");
1461 } else if (strcasecmp(field_name, "CornerField")== 0){
1462 vec_local = user->lCellFieldAtCorner; dm = user->fda; dof = 3; strcpy(data_layout, "Node-Centered");
1463 } else {
1464 SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Unknown field name for LOG_FIELD_ANATOMY: %s", field_name);
1465 }
1466
1467 // --- 2. Get Grid Info and Array Pointers ---
1468 ierr = DMDAGetLocalInfo(dm, &info); CHKERRQ(ierr);
1469
1470 // Synchronize for clean output
1471 ierr = PetscBarrier(NULL);
1472 PetscPrintf(PETSC_COMM_WORLD, "\n--- Field Anatomy Log: [%s] | Stage: [%s] | Layout: [%s] ---\n", field_name, stage_name, data_layout);
1473
1474 // We will check a slice at the center of the other two dimensions (in the rank's local domain)
1475
1476 PetscInt im_phys = user->IM; // Physical size in I-direction
1477 PetscInt jm_phys = user->JM; // Physical size in J-direction
1478 PetscInt km_phys = user->KM; // Physical size in K-direction
1479
1480 PetscInt i_mid = (PetscInt)(info.xs + 0.5*info.xm);
1481 PetscInt j_mid = (PetscInt)(info.ys + 0.5*info.ym);
1482 PetscInt k_mid = (PetscInt)(info.zs + 0.5*info.zm);
1483
1484 // --- 3. Print Boundary Information based on Data Layout ---
1485
1486 // ======================================================================
1487 // === CASE 1: Cell-Centered Fields (Ucat, P) - USES SHIFTED INDEX ===
1488 // ======================================================================
1489 if (strcmp(data_layout, "Cell-Centered") == 0) {
1490 const void *l_arr; // Use void pointer for generic array access
1491 ierr = DMDAVecGetArrayRead(dm, vec_local, (void*)&l_arr); CHKERRQ(ierr);
1492
1493 // --- I-Direction Boundaries ---
1494 if (info.xs == 0) { // Rank on -Xi boundary
1495 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, I-DIR]: Idx %2d (Ghost for Cell 0) = ", rank, 0);
1496 if(dof==1) PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.3e)\n", ((const PetscReal***)l_arr)[k_mid][j_mid][0]);
1497 else PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.3e, %.3e, %.3e)\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);
1498
1499 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, I-DIR]: Idx %2d (Value for Cell 0) = ", rank, 1);
1500 if(dof==1) PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.3e)\n", ((const PetscReal***)l_arr)[k_mid][j_mid][1]);
1501 else PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.3e, %.3e, %.3e)\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);
1502 }
1503 if (info.xs + info.xm == info.mx) { // Rank on +Xi boundary
1504 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, I-DIR]: Idx %2d (Value for Cell %d) = ", rank, im_phys - 1, im_phys - 2);
1505 if(dof==1) PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.3e)\n", ((const PetscReal***)l_arr)[k_mid][j_mid][im_phys - 1]);
1506 else PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.3e, %.3e, %.3e)\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);
1507
1508 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, I-DIR]: Idx %2d (Ghost for Cell %d) = ", rank, im_phys, im_phys - 2);
1509 if(dof==1) PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.3e)\n", ((const PetscReal***)l_arr)[k_mid][j_mid][im_phys]);
1510 else PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.3e, %.3e, %.3e)\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);
1511 }
1512
1513 // --- J-Direction Boundaries ---
1514 if (info.ys == 0) { // Rank on -Eta boundary
1515 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, J-DIR]: Idx %2d (Ghost for Cell 0) = ", rank, 0);
1516 if(dof==1) PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.3e)\n", ((const PetscReal***)l_arr)[k_mid][0][i_mid]);
1517 else PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.3e, %.3e, %.3e)\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);
1518
1519 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, J-DIR]: Idx %2d (Value for Cell 0) = ", rank, 1);
1520 if(dof==1) PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.3e)\n", ((const PetscReal***)l_arr)[k_mid][1][i_mid]);
1521 else PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.3e, %.3e, %.3e)\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);
1522 }
1523
1524 if (info.ys + info.ym == info.my) { // Rank on +Eta boundary
1525 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, J-DIR]: Idx %2d (Value for Cell %d) = ", rank, jm_phys - 1, jm_phys - 2);
1526 if(dof==1) PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.3e)\n", ((const PetscReal***)l_arr)[k_mid][jm_phys - 1][i_mid]);
1527 else PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.3e, %.3e, %.3e)\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);
1528
1529 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, J-DIR]: Idx %2d (Ghost for Cell %d) = ", rank, jm_phys, jm_phys - 2);
1530 if(dof==1) PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.3e)\n", ((const PetscReal***)l_arr)[k_mid][jm_phys][i_mid]);
1531 else PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.3e, %.3e, %.3e)\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);
1532 }
1533
1534 // --- K-Direction Boundaries ---
1535 if (info.zs == 0) { // Rank on -Zeta boundary
1536 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, K-DIR]: Idx %2d (Ghost for Cell 0) = ", rank, 0);
1537 if(dof==1) PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.3e)\n", ((const PetscReal***)l_arr)[0][j_mid][i_mid]);
1538 else PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.3e, %.3e, %.3e)\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);
1539 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, K-DIR]: Idx %2d (Value for Cell 0) = ", rank, 1);
1540 if(dof==1) PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.3e)\n", ((const PetscReal***)l_arr)[1][j_mid][i_mid]);
1541 else PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.3e, %.3e, %.3e)\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);
1542 }
1543 if (info.zs + info.zm == info.mz) { // Rank on +Zeta boundary
1544 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, K-DIR]: Idx %2d (Value for Cell %d) = ", rank, km_phys - 1, km_phys - 2);
1545 if(dof==1) PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.3e)\n", ((const PetscReal***)l_arr)[km_phys - 1][j_mid][i_mid]);
1546 else PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.3e, %.3e, %.3e)\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);
1547 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, K-DIR]: Idx %2d (Ghost for Cell %d) = ", rank, km_phys, km_phys - 2);
1548 if(dof==1) PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.3e)\n", ((const PetscReal***)l_arr)[km_phys][j_mid][i_mid]);
1549 else PetscSynchronizedPrintf(PETSC_COMM_WORLD, "(%.3e, %.3e, %.3e)\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);
1550 }
1551
1552 ierr = DMDAVecRestoreArrayRead(dm, vec_local, (void*)&l_arr); CHKERRQ(ierr);
1553 }
1554 // ======================================================================
1555 // === CASE 2: Face-Centered & Node-Centered Fields - USES DIRECT INDEX ===
1556 // ======================================================================
1557 else if (strcmp(data_layout, "Face-Centered") == 0 || strcmp(data_layout, "Node-Centered") == 0) {
1558 const Cmpnts ***l_arr; // Both Ucont and Coordinates are Cmpnts
1559 ierr = DMDAVecGetArrayRead(dm, vec_local, (void*)&l_arr); CHKERRQ(ierr);
1560
1561 // --- I-Direction Boundaries ---
1562 if (info.xs == 0) { // Rank on -Xi boundary
1563 if(strcmp(data_layout, "Face-Centered")==0){
1564 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, I-DIR]: Idx %2d (First Physical Face) = (%.3e)\n", rank, 0,
1565 l_arr[k_mid][j_mid][0].x);
1566 }else if(strcmp(data_layout, "Node-Centered")==0){// Node-Centered
1567 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, I-DIR]: Idx %2d (First Physical Face/Node) = (%.3e, %.3e, %.3e)\n", rank, 0,
1568 l_arr[k_mid][j_mid][0].x, l_arr[k_mid][j_mid][0].y, l_arr[k_mid][j_mid][0].z);
1569 }
1570 }
1571 if (info.xs + info.xm == info.mx) { // Rank on +Xi boundary
1572 if(strcmp(data_layout, "Face-Centered")==0){
1573 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, I-DIR]: Idx %2d (Last Physical Face) = (%.3e)\n", rank, im_phys - 1,
1574 l_arr[k_mid][j_mid][im_phys - 1].x);
1575 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, I-DIR]: Idx %2d (Unused/Ghost Location) = (%.3e)\n", rank, im_phys,
1576 l_arr[k_mid][j_mid][im_phys].x);
1577 }else if(strcmp(data_layout, "Node-Centered")==0){// Node-Centered
1578 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, I-DIR]: Idx %2d (Last Physical Face/Node) = (%.3e, %.3e, %.3e)\n", rank, im_phys - 1,
1579 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);
1580 // Also show the value at the unused memory location for verification
1581 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, I-DIR]: Idx %2d (Unused/Ghost Location) = (%.3e, %.3e, %.3e)\n", rank, im_phys,
1582 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);
1583 }
1584 }
1585
1586 // --- J-Direction
1587
1588 if (info.ys == 0) { // Rank on -Eta boundary
1589 if(strcmp(data_layout, "Face-Centered")==0){
1590 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, J-DIR]: Idx %2d (First Physical Face) = (%.3e)\n", rank, 0,
1591 l_arr[k_mid][0][i_mid].y);
1592 }else if(strcmp(data_layout, "Node-Centered")==0){// Node-Centered
1593 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, J-DIR]: Idx %2d (First Physical Face/Node) = (%.3e, %.3e, %.3e)\n", rank, 0,
1594 l_arr[k_mid][0][i_mid].x, l_arr[k_mid][0][i_mid].y, l_arr[k_mid][0][i_mid].z);
1595 }
1596 }
1597 if (info.ys + info.ym == info.my) { // Rank on +Eta boundary
1598 if(strcmp(data_layout, "Face-Centered")==0){
1599 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, J-DIR]: Idx %2d (Last Physical Face) = (%.3e)\n", rank, jm_phys - 1,
1600 l_arr[k_mid][jm_phys - 1][i_mid].y);
1601 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, J-DIR]: Idx %2d (Unused/Ghost Location) = (%.3e)\n", rank, jm_phys,
1602 l_arr[k_mid][jm_phys][i_mid].y);
1603 }else if(strcmp(data_layout, "Node-Centered")==0){// Node-Centered
1604 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, J-DIR]: Idx %2d (Last Physical Face/Node) = (%.3e, %.3e, %.3e)\n", rank, jm_phys - 1,
1605 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);
1606 // Also show the value at the unused memory location for verification
1607 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, J-DIR]: Idx %2d (Unused/Ghost Location) = (%.3e, %.3e, %.3e)\n", rank, jm_phys,
1608 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);
1609 }
1610 }
1611 // --- K-Direction ---
1612 if (info.zs == 0) { // Rank on -Zeta boundary
1613 if(strcmp(data_layout, "Face-Centered")==0){
1614 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, K-DIR]: Idx %2d (First Physical Face) = (%.3e)\n", rank, 0,
1615 l_arr[0][j_mid][i_mid].z);
1616 }else if(strcmp(data_layout, "Node-Centered")==0){// Node-Centered
1617 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, K-DIR]: Idx %2d (First Physical Face/Node) = (%.3e, %.3e, %.3e)\n", rank, 0,
1618 l_arr[0][j_mid][i_mid].x, l_arr[0][j_mid][i_mid].y, l_arr[0][j_mid][i_mid].z);
1619 }
1620 }
1621 if(info.zs + info.zm == info.mz) { // Rank on +Zeta boundary
1622 if(strcmp(data_layout, "Face-Centered")==0){
1623 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, K-DIR]: Idx %2d (Last Physical Face) = (%.3e)\n", rank, km_phys - 1,
1624 l_arr[km_phys - 1][j_mid][i_mid].z);
1625 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, K-DIR]: Idx %2d (Unused/Ghost Location) = (%.3e)\n", rank, km_phys,
1626 l_arr[km_phys][j_mid][i_mid].z);
1627 }else if(strcmp(data_layout, "Node-Centered")==0){// Node-Centered
1628 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, K-DIR]: Idx %2d (Last Physical Face/Node) = (%.3e, %.3e, %.3e)\n", rank, km_phys - 1,
1629 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);
1630 // Also show the value at the unused memory location for verification
1631 PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[Rank %d, K-DIR]: Idx %2d (Unused/Ghost Location) = (%.3e, %.3e, %.3e)\n", rank, km_phys,
1632 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);
1633 }
1634 }
1635
1636 ierr = DMDAVecRestoreArrayRead(dm, vec_local, (void*)&l_arr); CHKERRQ(ierr);
1637 }
1638 else {
1639 SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "LOG_FIELD_ANATOMY only supports fields with 1 or 3 components & certain data-layouts, but field '%s' has %D components and an unsupported data-layout %s \n", field_name, dof,data_layout);
1640 }
1641 ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT); CHKERRQ(ierr);
1642 ierr = PetscBarrier(NULL);
1643 PetscFunctionReturn(0);
1644}
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:943
const char * BCHandlerTypeToString(BCHandlerType handler_type)
Converts a BCHandlerType enum to its string representation.
Definition logging.c:732
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:948
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:770
#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
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:940
static PetscInt g_profiler_capacity
Definition logging.c:949
const char * ParticleInitializationToString(PetscInt ParticleInitialization)
Helper function to convert ParticleInitialization to a string representation.
Definition logging.c:675
PetscErrorCode ProfilingFinalize(SimCtx *simCtx)
the profiling excercise and build a profiling summary which is then printed to a log file.
Definition logging.c:1106
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:1273
static int gNumAllowed
Number of entries in the gAllowedFunctions array.
Definition logging.c:26
double total_time
Definition logging.c:938
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:1207
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:1435
static PetscErrorCode _FindOrCreateEntry(const char *func_name, PetscInt *idx)
Definition logging.c:952
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:1022
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:705
const char * ParticleLocationStatusToString(ParticleLocationStatus level)
A function that outputs the name of the current level in the ParticleLocation enum.
Definition logging.c:921
static int _CompareProfiledFunctions(const void *a, const void *b)
Definition logging.c:1089
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:802
PetscErrorCode LOG_CONTINUITY_METRICS(UserCtx *user)
Logs continuity metrics for a single block to a file.
Definition logging.c:862
double current_step_time
Definition logging.c:939
long long current_step_call_count
Definition logging.c:941
static ProfiledFunction * g_profiler_registry
Definition logging.c:947
PetscErrorCode ProfilingInitialize(SimCtx *simCtx)
Initializes the custom profiling system using configuration from SimCtx.
Definition logging.c:998
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
const char * LESFlagToString(PetscInt LESFlag)
Helper function to convert LES Flag to a string representation.
Definition logging.c:690
PetscErrorCode ProfilingResetTimestepCounters(void)
Definition logging.c:1037
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:942
PetscErrorCode ProfilingLogTimestepSummary(PetscInt step)
Logs the performance summary for the current timestep and resets timers.
Definition logging.c:1047
void _ProfilingStart(const char *func_name)
Definition logging.c:1015
const char * name
Definition logging.c:937
Logging utilities and macros for PETSc-based applications.
#define LOG_ALLOW_SYNC(scope, level, fmt,...)
----β€” DEBUG ---------------------------------------β€” #define LOG_ALLOW(scope, level,...
Definition logging.h:268
PetscBool log_to_console
Definition logging.h:59
#define LOCAL
Logging scope definitions for controlling message output.
Definition logging.h:46
#define GLOBAL
Scope for global logging across all processes.
Definition logging.h:47
#define LOG_ALLOW(scope, level, fmt,...)
Logging macro that checks both the log level and whether the calling function is in the allowed-funct...
Definition logging.h:201
PetscReal bnorm
Definition logging.h:60
PetscInt step
Definition logging.h:61
LogLevel
Enumeration of logging levels.
Definition logging.h:28
@ LOG_ERROR
Critical errors that may halt the program.
Definition logging.h:29
@ LOG_PROFILE
Exclusive log level for performance timing and profiling.
Definition logging.h:31
@ LOG_TRACE
Very fine-grained tracing information for in-depth debugging.
Definition logging.h:34
@ LOG_INFO
Informational messages about program execution.
Definition logging.h:32
@ LOG_WARNING
Non-critical issues that warrant attention.
Definition logging.h:30
@ LOG_DEBUG
Detailed debugging information.
Definition logging.h:33
@ LOG_VERBOSE
Extremely detailed logs, typically for development use only.
Definition logging.h:35
FILE * file_handle
Definition logging.h:58
PetscInt block_id
Definition logging.h:62
Context for a dual-purpose KSP monitor.
Definition logging.h:57
BCType
Defines the general mathematical/physical category of a boundary.
Definition variables.h:207
@ INLET
Definition variables.h:214
@ NOGRAD
Definition variables.h:224
@ INTERFACE
Definition variables.h:209
@ FARFIELD
Definition variables.h:215
@ OUTLET
Definition variables.h:213
@ PERIODIC
Definition variables.h:216
@ WALL
Definition variables.h:210
PetscMPIInt rank
Definition variables.h:541
SimCtx * simCtx
Back-pointer to the master simulation context.
Definition variables.h:664
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:602
PetscInt KM
Definition variables.h:670
BCHandlerType
Defines the specific computational "strategy" for a boundary handler.
Definition variables.h:228
@ BC_HANDLER_INLET_PULSATILE_FLUX
Definition variables.h:232
@ BC_HANDLER_INLET_PARABOLIC
Definition variables.h:231
@ BC_HANDLER_INLET_CONSTANT_VELOCITY
Definition variables.h:232
@ BC_HANDLER_INTERFACE_OVERSET
Definition variables.h:236
@ BC_HANDLER_WALL_MOVING
Definition variables.h:230
@ BC_HANDLER_NOGRAD_COPY_GHOST
Definition variables.h:229
@ BC_HANDLER_WALL_NOSLIP
Definition variables.h:230
@ BC_HANDLER_OUTLET_CONSERVATION
Definition variables.h:234
@ BC_HANDLER_PERIODIC
Definition variables.h:236
@ BC_HANDLER_FARFIELD_NONREFLECTING
Definition variables.h:235
@ BC_HANDLER_OUTLET_PRESSURE
Definition variables.h:234
@ BC_HANDLER_SYMMETRY_PLANE
Definition variables.h:231
@ BC_HANDLER_UNDEFINED
Definition variables.h:229
Vec lCellFieldAtCorner
Definition variables.h:693
PetscInt _this
Definition variables.h:674
Vec lPsi
Definition variables.h:730
char ** criticalFuncs
Definition variables.h:643
PetscInt StartStep
Definition variables.h:548
PetscScalar x
Definition variables.h:101
PetscReal MaxDiv
Definition variables.h:638
PetscInt MaxDivx
Definition variables.h:639
PetscInt MaxDivy
Definition variables.h:639
PetscInt MaxDivz
Definition variables.h:639
char log_dir[PETSC_MAX_PATH_LEN]
Definition variables.h:562
PetscInt MaxDivFlatArg
Definition variables.h:639
PetscReal FluxInSum
Definition variables.h:602
PetscScalar z
Definition variables.h:101
Vec Ucat
Definition variables.h:688
PetscInt JM
Definition variables.h:670
Vec lUcont
Definition variables.h:688
PetscInt step
Definition variables.h:546
Vec lUcat
Definition variables.h:688
PetscScalar y
Definition variables.h:101
PetscInt IM
Definition variables.h:670
@ 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
PetscInt nCriticalFuncs
Definition variables.h:644
PetscReal summationRHS
Definition variables.h:637
Cmpnts vertices[8]
Coordinates of the eight vertices of the cell.
Definition variables.h:161
Vec Psi
Definition variables.h:730
BCFace
Identifies the six logical faces of a structured computational block.
Definition variables.h:200
@ BC_FACE_NEG_X
Definition variables.h:201
@ BC_FACE_POS_Z
Definition variables.h:203
@ BC_FACE_POS_Y
Definition variables.h:202
@ BC_FACE_NEG_Z
Definition variables.h:203
@ BC_FACE_POS_X
Definition variables.h:201
@ BC_FACE_NEG_Y
Definition variables.h:202
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:538
User-defined context containing data specific to a single computational grid level.
Definition variables.h:661