PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
Loading...
Searching...
No Matches
test_logging.c
Go to the documentation of this file.
1/**
2 * @file test_logging.c
3 * @brief C unit tests for runtime log-level, allow-list, conversion, and profiling helpers.
4 */
5
6#include "test_support.h"
7
8#include "logging.h"
9#include "interpolation.h"
10#include "setup.h"
11
12#include <fcntl.h>
13#include <math.h>
14#include <stdio.h>
15#include <stdlib.h>
16#include <string.h>
17#include <unistd.h>
18/**
19 * @brief Asserts that one text file contains a required substring.
20 */
21
22static PetscErrorCode AssertFileContains(const char *path, const char *needle, const char *context)
23{
24 FILE *fp = NULL;
25 long file_size = 0;
26 char *buffer = NULL;
27
28 PetscFunctionBeginUser;
29 fp = fopen(path, "rb");
30 if (!fp) {
31 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Failed to open '%s' for assertion.", path);
32 }
33 if (fseek(fp, 0, SEEK_END) != 0) {
34 fclose(fp);
35 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Failed to seek '%s'.", path);
36 }
37 file_size = ftell(fp);
38 if (file_size < 0) {
39 fclose(fp);
40 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Failed to measure '%s'.", path);
41 }
42 if (fseek(fp, 0, SEEK_SET) != 0) {
43 fclose(fp);
44 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Failed to rewind '%s'.", path);
45 }
46
47 PetscCall(PetscMalloc1((size_t)file_size + 1, &buffer));
48 if (file_size > 0) {
49 size_t bytes_read = fread(buffer, 1, (size_t)file_size, fp);
50 if (bytes_read != (size_t)file_size) {
51 fclose(fp);
52 PetscCall(PetscFree(buffer));
53 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Failed to read '%s'.", path);
54 }
55 }
56 buffer[file_size] = '\0';
57 fclose(fp);
58
59 PetscCall(PicurvAssertBool((PetscBool)(strstr(buffer, needle) != NULL), context));
60 PetscCall(PetscFree(buffer));
61 PetscFunctionReturn(0);
62}
63/**
64 * @brief Asserts that one text file does not contain an excluded substring.
65 */
66
67static PetscErrorCode AssertFileNotContains(const char *path, const char *needle, const char *context)
68{
69 FILE *fp = NULL;
70 long file_size = 0;
71 char *buffer = NULL;
72
73 PetscFunctionBeginUser;
74 fp = fopen(path, "rb");
75 if (!fp) {
76 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Failed to open '%s' for assertion.", path);
77 }
78 if (fseek(fp, 0, SEEK_END) != 0) {
79 fclose(fp);
80 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Failed to seek '%s'.", path);
81 }
82 file_size = ftell(fp);
83 if (file_size < 0) {
84 fclose(fp);
85 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Failed to measure '%s'.", path);
86 }
87 if (fseek(fp, 0, SEEK_SET) != 0) {
88 fclose(fp);
89 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Failed to rewind '%s'.", path);
90 }
91
92 PetscCall(PetscMalloc1((size_t)file_size + 1, &buffer));
93 if (file_size > 0) {
94 size_t bytes_read = fread(buffer, 1, (size_t)file_size, fp);
95 if (bytes_read != (size_t)file_size) {
96 fclose(fp);
97 PetscCall(PetscFree(buffer));
98 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Failed to read '%s'.", path);
99 }
100 }
101 buffer[file_size] = '\0';
102 fclose(fp);
103
104 PetscCall(PicurvAssertBool((PetscBool)(strstr(buffer, needle) == NULL), context));
105 PetscCall(PetscFree(buffer));
106 PetscFunctionReturn(0);
107}
108
109/**
110 * @brief Reads the CSV header and the requested 1-based data row from a text file.
111 */
112static PetscErrorCode ReadCsvHeaderAndRow(const char *path,
113 PetscInt row_index,
114 char *header,
115 size_t header_len,
116 char *row,
117 size_t row_len)
118{
119 FILE *fp = NULL;
120 PetscInt current_row = 0;
121
122 PetscFunctionBeginUser;
123 PetscCheck(path != NULL, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "CSV path cannot be NULL.");
124 PetscCheck(header != NULL && header_len > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "CSV header buffer cannot be NULL or empty.");
125 PetscCheck(row != NULL && row_len > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "CSV row buffer cannot be NULL or empty.");
126 PetscCheck(row_index >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "CSV row index must be >= 1.");
127
128 fp = fopen(path, "r");
129 PetscCheck(fp != NULL, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Failed to open CSV '%s'.", path);
130 PetscCheck(fgets(header, (int)header_len, fp) != NULL, PETSC_COMM_SELF, PETSC_ERR_FILE_READ, "CSV header missing in '%s'.", path);
131
132 while (fgets(row, (int)row_len, fp) != NULL) {
133 current_row++;
134 if (current_row == row_index) {
135 fclose(fp);
136 PetscFunctionReturn(0);
137 }
138 }
139
140 fclose(fp);
141 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_READ, "CSV '%s' does not contain data row %d.", path, (int)row_index);
142}
143
144/**
145 * @brief Returns the zero-based column index for one CSV header field.
146 */
147static PetscErrorCode CsvFindColumnIndex(const char *header, const char *column_name, PetscInt *index_out)
148{
149 char local_header[4096];
150 char *saveptr = NULL;
151 char *token = NULL;
152 PetscInt index = 0;
153
154 PetscFunctionBeginUser;
155 PetscCheck(header != NULL, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "CSV header cannot be NULL.");
156 PetscCheck(column_name != NULL, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "CSV column name cannot be NULL.");
157 PetscCheck(index_out != NULL, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "CSV column index output cannot be NULL.");
158 PetscCheck(strlen(header) < sizeof(local_header), PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "CSV header is too long for the local parser buffer.");
159
160 PetscCall(PetscStrncpy(local_header, header, sizeof(local_header)));
161 token = strtok_r(local_header, ",\r\n", &saveptr);
162 while (token != NULL) {
163 if (strcmp(token, column_name) == 0) {
164 *index_out = index;
165 PetscFunctionReturn(0);
166 }
167 token = strtok_r(NULL, ",\r\n", &saveptr);
168 index++;
169 }
170
171 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "CSV column '%s' was not found.", column_name);
172}
173
174/**
175 * @brief Extracts one CSV cell as text by header name.
176 */
177static PetscErrorCode CsvGetColumnText(const char *header,
178 const char *row,
179 const char *column_name,
180 char *value,
181 size_t value_len)
182{
183 char local_row[4096];
184 char *saveptr = NULL;
185 char *token = NULL;
186 PetscInt target_index = -1;
187 PetscInt index = 0;
188
189 PetscFunctionBeginUser;
190 PetscCheck(row != NULL, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "CSV row cannot be NULL.");
191 PetscCheck(value != NULL && value_len > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "CSV value buffer cannot be NULL or empty.");
192 PetscCheck(strlen(row) < sizeof(local_row), PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "CSV row is too long for the local parser buffer.");
193
194 PetscCall(CsvFindColumnIndex(header, column_name, &target_index));
195 PetscCall(PetscStrncpy(local_row, row, sizeof(local_row)));
196
197 token = strtok_r(local_row, ",\r\n", &saveptr);
198 while (token != NULL) {
199 if (index == target_index) {
200 PetscCall(PetscStrncpy(value, token, value_len));
201 PetscFunctionReturn(0);
202 }
203 token = strtok_r(NULL, ",\r\n", &saveptr);
204 index++;
205 }
206
207 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "CSV row is missing column '%s'.", column_name);
208}
209
210/**
211 * @brief Extracts one CSV cell as an integer by header name.
212 */
213static PetscErrorCode CsvGetColumnInt(const char *header, const char *row, const char *column_name, PetscInt *value_out)
214{
215 char text[256];
216 char *endptr = NULL;
217 long parsed = 0;
218
219 PetscFunctionBeginUser;
220 PetscCheck(value_out != NULL, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "CSV integer output cannot be NULL.");
221 PetscCall(CsvGetColumnText(header, row, column_name, text, sizeof(text)));
222 parsed = strtol(text, &endptr, 10);
223 PetscCheck(endptr != text, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "CSV column '%s' did not contain an integer.", column_name);
224 *value_out = (PetscInt)parsed;
225 PetscFunctionReturn(0);
226}
227
228/**
229 * @brief Extracts one CSV cell as a real by header name.
230 */
231static PetscErrorCode CsvGetColumnReal(const char *header, const char *row, const char *column_name, PetscReal *value_out)
232{
233 char text[256];
234 char *endptr = NULL;
235 double parsed = 0.0;
236
237 PetscFunctionBeginUser;
238 PetscCheck(value_out != NULL, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "CSV real output cannot be NULL.");
239 PetscCall(CsvGetColumnText(header, row, column_name, text, sizeof(text)));
240 parsed = strtod(text, &endptr);
241 PetscCheck(endptr != text, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "CSV column '%s' did not contain a real value.", column_name);
242 *value_out = (PetscReal)parsed;
243 PetscFunctionReturn(0);
244}
245
246/**
247 * @brief Reads the column header and the requested 1-based data row from a
248 * pipe-delimited solution-convergence log file.
249 *
250 * Lines starting with '=' (banner) or '-' (separator) are skipped. The first
251 * non-skipped line is treated as the column header; subsequent non-skipped
252 * lines are data rows numbered from 1.
253 */
254static PetscErrorCode ReadLogHeaderAndRow(const char *path,
255 PetscInt row_index,
256 char *header,
257 size_t header_len,
258 char *row,
259 size_t row_len)
260{
261 FILE *fp = NULL;
262 PetscInt current_row = 0;
263 char line[8192];
264 PetscBool header_found = PETSC_FALSE;
265
266 PetscFunctionBeginUser;
267 PetscCheck(path != NULL, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Log path cannot be NULL.");
268 PetscCheck(header != NULL && header_len > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Log header buffer cannot be NULL.");
269 PetscCheck(row != NULL && row_len > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Log row buffer cannot be NULL.");
270 PetscCheck(row_index >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Log row index must be >= 1.");
271
272 fp = fopen(path, "r");
273 PetscCheck(fp != NULL, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Failed to open log '%s'.", path);
274
275 while (fgets(line, sizeof(line), fp) != NULL) {
276 if (line[0] == '=' || line[0] == '-' || line[0] == '\n' || line[0] == '\r') continue;
277 if (!header_found) {
278 PetscCall(PetscStrncpy(header, line, header_len));
279 header_found = PETSC_TRUE;
280 continue;
281 }
282 current_row++;
283 if (current_row == row_index) {
284 PetscCall(PetscStrncpy(row, line, row_len));
285 fclose(fp);
286 PetscFunctionReturn(0);
287 }
288 }
289
290 fclose(fp);
291 if (!header_found) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_READ, "Log '%s' has no column header.", path);
292 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_READ, "Log '%s' does not contain data row %d.", path, (int)row_index);
293}
294
295/**
296 * @brief Returns the zero-based column index for one pipe-delimited header field.
297 */
298static PetscErrorCode LogFindColumnIndex(const char *header, const char *column_name, PetscInt *index_out)
299{
300 char local_header[8192];
301 char *saveptr = NULL;
302 char *token = NULL;
303 PetscInt index = 0;
304
305 PetscFunctionBeginUser;
306 PetscCheck(header != NULL, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Log header cannot be NULL.");
307 PetscCheck(column_name != NULL, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Log column name cannot be NULL.");
308 PetscCheck(index_out != NULL, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Log column index output cannot be NULL.");
309
310 PetscCall(PetscStrncpy(local_header, header, sizeof(local_header)));
311 token = strtok_r(local_header, "|\r\n", &saveptr);
312 while (token != NULL) {
313 while (*token == ' ') token++;
314 char *end = token + strlen(token) - 1;
315 while (end > token && (*end == ' ' || *end == '\r' || *end == '\n')) end--;
316 *(end + 1) = '\0';
317 if (strcmp(token, column_name) == 0) {
318 *index_out = index;
319 PetscFunctionReturn(0);
320 }
321 token = strtok_r(NULL, "|\r\n", &saveptr);
322 index++;
323 }
324 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Log column '%s' was not found.", column_name);
325}
326
327/**
328 * @brief Extracts one pipe-delimited log cell as text by header name.
329 */
330static PetscErrorCode LogGetColumnText(const char *header,
331 const char *row,
332 const char *column_name,
333 char *value,
334 size_t value_len)
335{
336 char local_row[8192];
337 char *saveptr = NULL;
338 char *token = NULL;
339 PetscInt target_index = -1;
340 PetscInt index = 0;
341
342 PetscFunctionBeginUser;
343 PetscCheck(row != NULL, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Log row cannot be NULL.");
344 PetscCheck(value != NULL && value_len > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Log value buffer cannot be NULL.");
345
346 PetscCall(LogFindColumnIndex(header, column_name, &target_index));
347 PetscCall(PetscStrncpy(local_row, row, sizeof(local_row)));
348
349 token = strtok_r(local_row, "|\r\n", &saveptr);
350 while (token != NULL) {
351 if (index == target_index) {
352 while (*token == ' ') token++;
353 char *end = token + strlen(token) - 1;
354 while (end > token && (*end == ' ' || *end == '\r' || *end == '\n')) end--;
355 *(end + 1) = '\0';
356 PetscCall(PetscStrncpy(value, token, value_len));
357 PetscFunctionReturn(0);
358 }
359 token = strtok_r(NULL, "|\r\n", &saveptr);
360 index++;
361 }
362 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Log row is missing column '%s'.", column_name);
363}
364
365/**
366 * @brief Extracts one pipe-delimited log cell as an integer by header name.
367 */
368static PetscErrorCode LogGetColumnInt(const char *header, const char *row, const char *column_name, PetscInt *value_out)
369{
370 char text[256];
371 char *endptr = NULL;
372 long parsed = 0;
373
374 PetscFunctionBeginUser;
375 PetscCheck(value_out != NULL, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Log integer output cannot be NULL.");
376 PetscCall(LogGetColumnText(header, row, column_name, text, sizeof(text)));
377 parsed = strtol(text, &endptr, 10);
378 PetscCheck(endptr != text, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Log column '%s' did not contain an integer.", column_name);
379 *value_out = (PetscInt)parsed;
380 PetscFunctionReturn(0);
381}
382
383/**
384 * @brief Extracts one pipe-delimited log cell as a real by header name.
385 */
386static PetscErrorCode LogGetColumnReal(const char *header, const char *row, const char *column_name, PetscReal *value_out)
387{
388 char text[256];
389 char *endptr = NULL;
390 double parsed = 0.0;
391
392 PetscFunctionBeginUser;
393 PetscCheck(value_out != NULL, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Log real output cannot be NULL.");
394 PetscCall(LogGetColumnText(header, row, column_name, text, sizeof(text)));
395 parsed = strtod(text, &endptr);
396 PetscCheck(endptr != text, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Log column '%s' did not contain a real value.", column_name);
397 *value_out = (PetscReal)parsed;
398 PetscFunctionReturn(0);
399}
400
401typedef PetscErrorCode (*CapturedLoggingFn)(UserCtx *user, SimCtx *simCtx, void *ctx);
402
403typedef struct AnatomyCaptureCtx {
404 const char *field_name;
405 const char *stage_name;
407
408/**
409 * @brief Captures stdout emitted by one logging helper into a temporary file-backed buffer.
410 */
411static PetscErrorCode CaptureLoggingOutput(UserCtx *user,
412 SimCtx *simCtx,
414 void *ctx,
415 char *captured,
416 size_t captured_len)
417{
418 char tmpdir[PETSC_MAX_PATH_LEN];
419 char capture_path[PETSC_MAX_PATH_LEN];
420 FILE *capture_file = NULL;
421 int saved_stdout = -1;
422 int capture_fd = -1;
423 size_t bytes_read = 0;
424 PetscErrorCode ierr;
425
426 PetscFunctionBeginUser;
427 PetscCheck(fn != NULL, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Capture callback cannot be NULL.");
428 PetscCheck(captured != NULL, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Capture buffer cannot be NULL.");
429 PetscCheck(captured_len > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Capture buffer must be non-empty.");
430
431 PetscCall(PicurvMakeTempDir(tmpdir, sizeof(tmpdir)));
432 PetscCall(PetscSNPrintf(capture_path, sizeof(capture_path), "%s/logging.out", tmpdir));
433
434 fflush(stdout);
435 saved_stdout = dup(STDOUT_FILENO);
436 PetscCheck(saved_stdout >= 0, PETSC_COMM_SELF, PETSC_ERR_SYS, "dup(STDOUT_FILENO) failed.");
437 capture_fd = open(capture_path, O_CREAT | O_TRUNC | O_WRONLY, 0600);
438 PetscCheck(capture_fd >= 0, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Failed to open capture file '%s'.", capture_path);
439 PetscCheck(dup2(capture_fd, STDOUT_FILENO) >= 0, PETSC_COMM_SELF, PETSC_ERR_SYS, "dup2() failed while redirecting stdout.");
440 close(capture_fd);
441 capture_fd = -1;
442
443 ierr = fn(user, simCtx, ctx);
444 fflush(stdout);
445 PetscCheck(dup2(saved_stdout, STDOUT_FILENO) >= 0, PETSC_COMM_SELF, PETSC_ERR_SYS, "dup2() failed while restoring stdout.");
446 close(saved_stdout);
447 saved_stdout = -1;
448 PetscCall(ierr);
449
450 capture_file = fopen(capture_path, "r");
451 PetscCheck(capture_file != NULL, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Failed to read capture file '%s'.", capture_path);
452 bytes_read = fread(captured, 1, captured_len - 1, capture_file);
453 captured[bytes_read] = '\0';
454 fclose(capture_file);
455 PetscCall(PicurvRemoveTempDir(tmpdir));
456 PetscFunctionReturn(0);
457}
458
459/**
460 * @brief Adapts `LOG_PARTICLE_FIELDS()` to the generic stdout-capture callback shape.
461 */
462static PetscErrorCode InvokeParticleFieldLog(UserCtx *user, SimCtx *simCtx, void *ctx)
463{
464 PetscInt print_interval = *((PetscInt *)ctx);
465
466 PetscFunctionBeginUser;
467 (void)simCtx;
468 PetscCall(LOG_PARTICLE_FIELDS(user, print_interval));
469 PetscFunctionReturn(0);
470}
471
472/**
473 * @brief Adapts `EmitParticleConsoleSnapshot()` to the generic stdout-capture callback shape.
474 */
475static PetscErrorCode InvokeParticleConsoleSnapshot(UserCtx *user, SimCtx *simCtx, void *ctx)
476{
477 PetscInt step = *((PetscInt *)ctx);
478
479 PetscFunctionBeginUser;
480 PetscCall(EmitParticleConsoleSnapshot(user, simCtx, step));
481 PetscFunctionReturn(0);
482}
483
484/**
485 * @brief Adapts `LOG_FIELD_ANATOMY()` to the generic stdout-capture callback shape.
486 */
487static PetscErrorCode InvokeFieldAnatomyLog(UserCtx *user, SimCtx *simCtx, void *ctx)
488{
489 const AnatomyCaptureCtx *anatomy_ctx = (const AnatomyCaptureCtx *)ctx;
490
491 PetscFunctionBeginUser;
492 (void)simCtx;
493 PetscCall(LOG_FIELD_ANATOMY(user, anatomy_ctx->field_name, anatomy_ctx->stage_name));
494 PetscFunctionReturn(0);
495}
496
497/**
498 * @brief Creates a small particle-bearing runtime fixture used by logging tests.
499 */
500static PetscErrorCode SeedLoggingParticleFixture(SimCtx **simCtx_out, UserCtx **user_out)
501{
502 PetscReal *positions = NULL;
503 PetscReal *velocities = NULL;
504 PetscReal *weights = NULL;
505 PetscInt *cell_ids = NULL;
506 PetscInt *status = NULL;
507 PetscReal *psi = NULL;
508
509 PetscFunctionBeginUser;
510 PetscCall(PicurvCreateMinimalContexts(simCtx_out, user_out, 4, 4, 4));
511 PetscCall(PicurvCreateSwarmPair(*user_out, 2, "ske"));
512 (*simCtx_out)->np = 2;
513 (*simCtx_out)->particleConsoleOutputFreq = 2;
514 (*simCtx_out)->LoggingFrequency = 1;
515
516 PetscCall(DMSwarmGetField((*user_out)->swarm, "position", NULL, NULL, (void **)&positions));
517 PetscCall(DMSwarmGetField((*user_out)->swarm, "velocity", NULL, NULL, (void **)&velocities));
518 PetscCall(DMSwarmGetField((*user_out)->swarm, "weight", NULL, NULL, (void **)&weights));
519 PetscCall(DMSwarmGetField((*user_out)->swarm, "DMSwarm_CellID", NULL, NULL, (void **)&cell_ids));
520 PetscCall(DMSwarmGetField((*user_out)->swarm, "DMSwarm_location_status", NULL, NULL, (void **)&status));
521 PetscCall(DMSwarmGetField((*user_out)->swarm, "Psi", NULL, NULL, (void **)&psi));
522
523 positions[0] = 0.25; positions[1] = 0.50; positions[2] = 0.75;
524 positions[3] = 0.50; positions[4] = 0.50; positions[5] = 0.50;
525 velocities[0] = 1.0; velocities[1] = 2.0; velocities[2] = 3.0;
526 velocities[3] = 4.0; velocities[4] = 5.0; velocities[5] = 6.0;
527 weights[0] = 0.2; weights[1] = 0.3; weights[2] = 0.4;
528 weights[3] = 0.5; weights[4] = 0.5; weights[5] = 0.5;
529 cell_ids[0] = 0; cell_ids[1] = 0; cell_ids[2] = 0;
530 cell_ids[3] = 1; cell_ids[4] = 1; cell_ids[5] = 1;
531 status[0] = ACTIVE_AND_LOCATED;
532 status[1] = ACTIVE_AND_LOCATED;
533 psi[0] = 1.0;
534 psi[1] = 3.0;
535
536 PetscCall(DMSwarmRestoreField((*user_out)->swarm, "Psi", NULL, NULL, (void **)&psi));
537 PetscCall(DMSwarmRestoreField((*user_out)->swarm, "DMSwarm_location_status", NULL, NULL, (void **)&status));
538 PetscCall(DMSwarmRestoreField((*user_out)->swarm, "DMSwarm_CellID", NULL, NULL, (void **)&cell_ids));
539 PetscCall(DMSwarmRestoreField((*user_out)->swarm, "weight", NULL, NULL, (void **)&weights));
540 PetscCall(DMSwarmRestoreField((*user_out)->swarm, "velocity", NULL, NULL, (void **)&velocities));
541 PetscCall(DMSwarmRestoreField((*user_out)->swarm, "position", NULL, NULL, (void **)&positions));
542 PetscFunctionReturn(0);
543}
544
545/**
546 * @brief Fills one Cartesian velocity field with a uniform constant state.
547 */
548static PetscErrorCode SetUniformVelocityField(UserCtx *user, Vec field, PetscReal ux, PetscReal uy, PetscReal uz)
549{
550 Cmpnts ***arr = NULL;
551
552 PetscFunctionBeginUser;
553 PetscCheck(user != NULL, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "UserCtx cannot be NULL.");
554 PetscCheck(field != NULL, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Velocity field cannot be NULL.");
555
556 PetscCall(DMDAVecGetArray(user->fda, field, &arr));
557 for (PetscInt k = user->info.zs; k < user->info.zs + user->info.zm; ++k) {
558 for (PetscInt j = user->info.ys; j < user->info.ys + user->info.ym; ++j) {
559 for (PetscInt i = user->info.xs; i < user->info.xs + user->info.xm; ++i) {
560 arr[k][j][i].x = ux;
561 arr[k][j][i].y = uy;
562 arr[k][j][i].z = uz;
563 }
564 }
565 }
566 PetscCall(DMDAVecRestoreArray(user->fda, field, &arr));
567 PetscFunctionReturn(0);
568}
569
570/**
571 * @brief Fills one scalar field with a uniform constant state.
572 */
573static PetscErrorCode SetUniformScalarField(UserCtx *user, Vec field, PetscReal value)
574{
575 PetscReal ***arr = NULL;
576
577 PetscFunctionBeginUser;
578 PetscCheck(user != NULL, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "UserCtx cannot be NULL.");
579 PetscCheck(field != NULL, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Scalar field cannot be NULL.");
580
581 PetscCall(DMDAVecGetArray(user->da, field, &arr));
582 for (PetscInt k = user->info.zs; k < user->info.zs + user->info.zm; ++k) {
583 for (PetscInt j = user->info.ys; j < user->info.ys + user->info.ym; ++j) {
584 for (PetscInt i = user->info.xs; i < user->info.xs + user->info.xm; ++i) {
585 arr[k][j][i] = value;
586 }
587 }
588 }
589 PetscCall(DMDAVecRestoreArray(user->da, field, &arr));
590 PetscFunctionReturn(0);
591}
592/**
593 * @brief Tests string-conversion helpers for configured enums and unknown values.
594 */
595
596static PetscErrorCode TestStringConversionHelpers(void)
597{
598 PetscFunctionBeginUser;
599 PetscCall(PicurvAssertBool(strcmp(BCFaceToString(BC_FACE_NEG_X), "-Xi (I-Min)") == 0,
600 "BCFaceToString should report the negative-x face"));
601 PetscCall(PicurvAssertBool(strcmp(FieldInitializationToString(0), "Zero") == 0,
602 "FieldInitializationToString should report the zero mode"));
603 PetscCall(PicurvAssertBool(strcmp(FieldInitializationToString(99), "Unknown Field Initialization") == 0,
604 "FieldInitializationToString should reject unknown selectors"));
606 "ParticleInitializationToString should report the volume mode"));
607 PetscCall(PicurvAssertBool(strcmp(LESModelToString(CONSTANT_SMAGORINSKY), "Constant Smagorinsky") == 0,
608 "LESModelToString should report the constant model"));
609 PetscCall(PicurvAssertBool(strcmp(MomentumSolverTypeToString(MOMENTUM_SOLVER_EXPLICIT_RK), "Explicit 4 stage Runge-Kutta ") == 0,
610 "MomentumSolverTypeToString should report the explicit solver"));
611 PetscCall(PicurvAssertBool(strcmp(BCTypeToString(PERIODIC), "PERIODIC") == 0,
612 "BCTypeToString should report periodic boundaries"));
614 "BCHandlerTypeToString should report the driven periodic handler"));
615 PetscCall(PicurvAssertBool(strcmp(ParticleLocationStatusToString(LOST), "LOST") == 0,
616 "ParticleLocationStatusToString should report LOST state"));
617 PetscFunctionReturn(0);
618}
619/**
620 * @brief Tests that log level selection honors the environment variable.
621 */
622
623static PetscErrorCode TestGetLogLevelFromEnvironment(void)
624{
625 PetscFunctionBeginUser;
627 "get_log_level should honor LOG_LEVEL=INFO in this test binary"));
628 PetscCall(print_log_level());
629 PetscFunctionReturn(0);
630}
631/**
632 * @brief Tests the function allow-list filter used by the logging layer.
633 */
634
635static PetscErrorCode TestAllowedFunctionsFilter(void)
636{
637 const char *allow_list[] = {"ComputeSpecificKE", "WriteEulerianFile"};
638
639 PetscFunctionBeginUser;
640 set_allowed_functions(allow_list, 2);
641 PetscCall(PicurvAssertBool(is_function_allowed("ComputeSpecificKE"),
642 "Allowed list should include ComputeSpecificKE"));
643 PetscCall(PicurvAssertBool((PetscBool)!is_function_allowed("UnlistedFunction"),
644 "Allowed list should exclude unknown function names"));
645
646 set_allowed_functions(NULL, 0);
647 PetscCall(PicurvAssertBool(is_function_allowed("AnyFunction"),
648 "Empty allow-list should permit all functions"));
649 PetscFunctionReturn(0);
650}
651/**
652 * @brief Tests periodic particle console snapshot enablement and cadence.
653 */
654
655static PetscErrorCode TestParticleConsoleSnapshotCadence(void)
656{
657 SimCtx simCtx;
658
659 PetscFunctionBeginUser;
660 PetscCall(PetscMemzero(&simCtx, sizeof(simCtx)));
661 simCtx.np = 32;
662 simCtx.particleConsoleOutputFreq = 4;
663
665 "Particle snapshot contract should be enabled when particles and cadence are configured"));
667 "Snapshot should emit on cadence-aligned completed steps"));
668 PetscCall(PicurvAssertBool((PetscBool)!ShouldEmitPeriodicParticleConsoleSnapshot(&simCtx, 7),
669 "Snapshot should not emit off-cadence"));
670
671 simCtx.particleConsoleOutputFreq = 0;
672 PetscCall(PicurvAssertBool((PetscBool)!IsParticleConsoleSnapshotEnabled(&simCtx),
673 "Zero cadence should disable periodic particle snapshots"));
674 PetscCall(PicurvAssertBool((PetscBool)!ShouldEmitPeriodicParticleConsoleSnapshot(NULL, 4),
675 "NULL SimCtx should never emit periodic snapshots"));
676 PetscFunctionReturn(0);
677}
678/**
679 * @brief Tests logging-side file parsing, helper formatting, and progress utilities.
680 */
681
683{
684 char tmpdir[PETSC_MAX_PATH_LEN];
685 char allow_path[PETSC_MAX_PATH_LEN];
686 char dual_log_path[PETSC_MAX_PATH_LEN];
687 FILE *file = NULL;
688 char **funcs = NULL;
689 PetscInt nfuncs = 0;
690 Cell cell;
691 PetscReal distances[6] = {1.0, 2.0, 3.0, 4.0, 5.0, 6.0};
692 DualMonitorCtx *monctx = NULL;
693 void *ctx = NULL;
694
695 PetscFunctionBeginUser;
696 PetscCall(PetscMemzero(&cell, sizeof(cell)));
697 PetscCall(PicurvMakeTempDir(tmpdir, sizeof(tmpdir)));
698 PetscCall(PetscSNPrintf(allow_path, sizeof(allow_path), "%s/allowed_functions.txt", tmpdir));
699 PetscCall(PetscSNPrintf(dual_log_path, sizeof(dual_log_path), "%s/dual-monitor.log", tmpdir));
700
701 file = fopen(allow_path, "w");
702 PetscCheck(file != NULL, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Failed to create allow-list file '%s'.", allow_path);
703 fputs(" ComputeSpecificKE \n", file);
704 fputs("# comment-only line\n", file);
705 for (PetscInt i = 0; i < 17; ++i) {
706 fprintf(file, "Helper_%02d # trailing comment\n", (int)i);
707 }
708 fclose(file);
709 file = NULL;
710
711 PetscCall(LoadAllowedFunctionsFromFile(allow_path, &funcs, &nfuncs));
712 PetscCall(PicurvAssertIntEqual(18, nfuncs, "LoadAllowedFunctionsFromFile should trim comments and keep all identifiers"));
713 PetscCall(PicurvAssertBool((PetscBool)(strcmp(funcs[0], "ComputeSpecificKE") == 0),
714 "LoadAllowedFunctionsFromFile should trim leading and trailing whitespace"));
715 PetscCall(PicurvAssertBool((PetscBool)(strcmp(funcs[17], "Helper_16") == 0),
716 "LoadAllowedFunctionsFromFile should grow past the initial pointer capacity"));
717 PetscCall(FreeAllowedFunctions(funcs, nfuncs));
718
719 for (PetscInt i = 0; i < 8; ++i) {
720 cell.vertices[i].x = (PetscReal)i;
721 cell.vertices[i].y = (PetscReal)(i + 1);
722 cell.vertices[i].z = (PetscReal)(i + 2);
723 }
724 PetscCall(LOG_CELL_VERTICES(&cell, 0));
725 PetscCall(LOG_FACE_DISTANCES(distances));
726
727 PetscCall(PetscCalloc1(1, &monctx));
728 monctx->file_handle = fopen(dual_log_path, "w");
729 PetscCheck(monctx->file_handle != NULL, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Failed to create dual-monitor log '%s'.", dual_log_path);
730 ctx = monctx;
731 PetscCall(DualMonitorDestroy(&ctx));
732 PetscCall(PicurvAssertBool((PetscBool)(ctx == NULL),
733 "DualMonitorDestroy should clear the caller-owned context pointer"));
734
735 PrintProgressBar(0, 0, 4, 0.10);
736 PrintProgressBar(3, 0, 4, 0.40);
737 PrintProgressBar(0, 0, 0, 0.00);
738 PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
739
740 PetscCall(PicurvRemoveTempDir(tmpdir));
741 PetscFunctionReturn(0);
742}
743/**
744 * @brief Tests continuity, min/max, and anatomy logging helpers on minimal runtime fixtures.
745 */
746
748{
749 SimCtx *simCtx = NULL;
750 UserCtx *user = NULL;
751 char tmpdir[PETSC_MAX_PATH_LEN];
752 char continuity_path[PETSC_MAX_PATH_LEN];
753 PetscReal ***p = NULL;
754 Cmpnts ***ucat = NULL;
755 Cmpnts ***ucont = NULL;
756 PetscErrorCode ierr_minmax = 0;
757
758 PetscFunctionBeginUser;
759 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
760 PetscCall(PicurvMakeTempDir(tmpdir, sizeof(tmpdir)));
761 PetscCall(PetscStrncpy(simCtx->log_dir, tmpdir, sizeof(simCtx->log_dir)));
762
763 simCtx->StartStep = 0;
764 simCtx->step = 1;
765 simCtx->MaxDiv = 1.25;
766 simCtx->MaxDivx = 1;
767 simCtx->MaxDivy = 2;
768 simCtx->MaxDivz = 3;
769 simCtx->MaxDivFlatArg = 17;
770 simCtx->summationRHS = 8.5;
771 simCtx->FluxInSum = 5.0;
772 simCtx->FluxOutSum = 3.25;
773 PetscCall(LOG_CONTINUITY_METRICS(user));
774
775 simCtx->step = 2;
776 simCtx->MaxDiv = 0.75;
777 simCtx->summationRHS = 4.5;
778 simCtx->FluxInSum = 2.5;
779 simCtx->FluxOutSum = 1.0;
780 PetscCall(LOG_CONTINUITY_METRICS(user));
781
782 PetscCall(PetscSNPrintf(continuity_path, sizeof(continuity_path), "%s/Continuity_Metrics.log", simCtx->log_dir));
783 PetscCall(PicurvAssertFileExists(continuity_path, "continuity metrics log should be written"));
784 PetscCall(AssertFileContains(continuity_path, "Timestep", "continuity metrics log should include the header"));
785 PetscCall(AssertFileContains(continuity_path, "([3][2][1] = 17)", "continuity metrics log should include the divergence location"));
786 PetscCall(AssertFileContains(continuity_path, "2 | 0", "continuity metrics log should append later timesteps"));
787
788 PetscCall(DMDAVecGetArray(user->da, user->P, &p));
789 PetscCall(DMDAVecGetArray(user->fda, user->Ucat, &ucat));
790 PetscCall(DMDAVecGetArray(user->fda, user->Ucont, &ucont));
791 for (PetscInt k = user->info.zs; k < user->info.zs + user->info.zm; ++k) {
792 for (PetscInt j = user->info.ys; j < user->info.ys + user->info.ym; ++j) {
793 for (PetscInt i = user->info.xs; i < user->info.xs + user->info.xm; ++i) {
794 p[k][j][i] = (PetscReal)(i + j + k);
795 ucat[k][j][i].x = (PetscReal)i;
796 ucat[k][j][i].y = (PetscReal)(-j);
797 ucat[k][j][i].z = (PetscReal)(2 * k);
798 ucont[k][j][i].x = (PetscReal)(10 + i);
799 ucont[k][j][i].y = (PetscReal)(20 + j);
800 ucont[k][j][i].z = (PetscReal)(30 + k);
801 }
802 }
803 }
804 PetscCall(DMDAVecRestoreArray(user->fda, user->Ucont, &ucont));
805 PetscCall(DMDAVecRestoreArray(user->fda, user->Ucat, &ucat));
806 PetscCall(DMDAVecRestoreArray(user->da, user->P, &p));
807
808 PetscCall(DMGlobalToLocalBegin(user->da, user->P, INSERT_VALUES, user->lP));
809 PetscCall(DMGlobalToLocalEnd(user->da, user->P, INSERT_VALUES, user->lP));
810 PetscCall(DMGlobalToLocalBegin(user->fda, user->Ucat, INSERT_VALUES, user->lUcat));
811 PetscCall(DMGlobalToLocalEnd(user->fda, user->Ucat, INSERT_VALUES, user->lUcat));
812 PetscCall(DMGlobalToLocalBegin(user->fda, user->Ucont, INSERT_VALUES, user->lUcont));
813 PetscCall(DMGlobalToLocalEnd(user->fda, user->Ucont, INSERT_VALUES, user->lUcont));
814
815 PetscCall(LOG_FIELD_MIN_MAX(user, "P"));
816 PetscCall(LOG_FIELD_MIN_MAX(user, "Ucat"));
817 PetscCall(LOG_FIELD_MIN_MAX(user, "Coordinates"));
818 PetscCall(LOG_FIELD_MIN_MAX(user, "Ucont"));
819
820 PetscCall(PetscPushErrorHandler(PetscIgnoreErrorHandler, NULL));
821 ierr_minmax = LOG_FIELD_MIN_MAX(user, "NotARealField");
822 PetscCall(PetscPopErrorHandler());
823 PetscCall(PicurvAssertBool((PetscBool)(ierr_minmax != 0),
824 "LOG_FIELD_MIN_MAX should reject unknown field names"));
825
826 PetscCall(PicurvRemoveTempDir(tmpdir));
827 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
828 PetscFunctionReturn(0);
829}
830/**
831 * @brief Tests interpolation-error logging against an analytically matched particle field.
832 */
833
834static PetscErrorCode TestInterpolationErrorLogging(void)
835{
836 SimCtx *simCtx = NULL;
837 UserCtx *user = NULL;
838 PetscReal (*pos_arr)[3] = NULL;
839 PetscReal (*vel_arr)[3] = NULL;
840 Vec position_vec = NULL;
841 Vec analytical_vec = NULL;
842 const PetscScalar *analytical_arr = NULL;
843
844 PetscFunctionBeginUser;
845 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
846 PetscCall(PicurvCreateSwarmPair(user, 2, "ske"));
847 PetscCall(PetscStrncpy(simCtx->AnalyticalSolutionType, "TGV3D", sizeof(simCtx->AnalyticalSolutionType)));
848 simCtx->ren = 1.0;
849 simCtx->ti = 0.0;
850
851 PetscCall(DMSwarmGetField(user->swarm, "position", NULL, NULL, (void *)&pos_arr));
852 pos_arr[0][0] = 0.5 * PETSC_PI; pos_arr[0][1] = 0.0; pos_arr[0][2] = 0.0;
853 pos_arr[1][0] = 0.0; pos_arr[1][1] = 0.5 * PETSC_PI; pos_arr[1][2] = 0.0;
854 PetscCall(DMSwarmRestoreField(user->swarm, "position", NULL, NULL, (void *)&pos_arr));
855
856 PetscCall(DMSwarmCreateGlobalVectorFromField(user->swarm, "position", &position_vec));
857 PetscCall(VecDuplicate(position_vec, &analytical_vec));
858 PetscCall(VecCopy(position_vec, analytical_vec));
859 PetscCall(SetAnalyticalSolutionForParticles(analytical_vec, simCtx));
860
861 PetscCall(DMSwarmGetField(user->swarm, "velocity", NULL, NULL, (void *)&vel_arr));
862 PetscCall(VecGetArrayRead(analytical_vec, &analytical_arr));
863 for (PetscInt particle = 0; particle < 2; ++particle) {
864 vel_arr[particle][0] = PetscRealPart(analytical_arr[3 * particle + 0]);
865 vel_arr[particle][1] = PetscRealPart(analytical_arr[3 * particle + 1]);
866 vel_arr[particle][2] = PetscRealPart(analytical_arr[3 * particle + 2]);
867 }
868 PetscCall(VecRestoreArrayRead(analytical_vec, &analytical_arr));
869 PetscCall(DMSwarmRestoreField(user->swarm, "velocity", NULL, NULL, (void *)&vel_arr));
870
871 PetscCall(VecDestroy(&analytical_vec));
872 PetscCall(DMSwarmDestroyGlobalVectorFromField(user->swarm, "position", &position_vec));
873 PetscCall(LOG_INTERPOLATION_ERROR(user));
874 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
875 PetscFunctionReturn(0);
876}
877/**
878 * @brief Tests file-backed scatter metrics logging against a fully occupied constant field.
879 */
880
881static PetscErrorCode TestScatterMetricsLogging(void)
882{
883 SimCtx *simCtx = NULL;
884 UserCtx *user = NULL;
885 char tmpdir[PETSC_MAX_PATH_LEN];
886 char metrics_path[PETSC_MAX_PATH_LEN];
887 PetscReal *positions = NULL;
888 PetscReal *psi = NULL;
889 PetscInt *cell_ids = NULL;
890 PetscInt *status = NULL;
891 PetscInt particle = 0;
892
893 PetscFunctionBeginUser;
894 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
895 PetscCall(PicurvCreateSwarmPair(user, 27, "ske"));
896 PetscCall(PicurvMakeTempDir(tmpdir, sizeof(tmpdir)));
897 PetscCall(PetscStrncpy(simCtx->log_dir, tmpdir, sizeof(simCtx->log_dir)));
898 simCtx->np = 27;
899 simCtx->step = 3;
900 simCtx->ti = 0.3;
901 simCtx->verificationScalar.enabled = PETSC_TRUE;
902 PetscCall(PetscStrncpy(simCtx->verificationScalar.mode,
903 "analytical",
904 sizeof(simCtx->verificationScalar.mode)));
905 PetscCall(PetscStrncpy(simCtx->verificationScalar.profile,
906 "CONSTANT",
907 sizeof(simCtx->verificationScalar.profile)));
908 simCtx->verificationScalar.value = 2.0;
909
910 PetscCall(DMSwarmGetField(user->swarm, "position", NULL, NULL, (void **)&positions));
911 PetscCall(DMSwarmGetField(user->swarm, "DMSwarm_CellID", NULL, NULL, (void **)&cell_ids));
912 PetscCall(DMSwarmGetField(user->swarm, "DMSwarm_location_status", NULL, NULL, (void **)&status));
913 PetscCall(DMSwarmGetField(user->swarm, "Psi", NULL, NULL, (void **)&psi));
914
915 for (PetscInt k = 0; k < 3; ++k) {
916 for (PetscInt j = 0; j < 3; ++j) {
917 for (PetscInt i = 0; i < 3; ++i) {
918 positions[3 * particle + 0] = (i + 0.5) / 4.0;
919 positions[3 * particle + 1] = (j + 0.5) / 4.0;
920 positions[3 * particle + 2] = (k + 0.5) / 4.0;
921 cell_ids[3 * particle + 0] = i;
922 cell_ids[3 * particle + 1] = j;
923 cell_ids[3 * particle + 2] = k;
924 status[particle] = ACTIVE_AND_LOCATED;
925 psi[particle] = 2.0;
926 ++particle;
927 }
928 }
929 }
930
931 PetscCall(DMSwarmRestoreField(user->swarm, "Psi", NULL, NULL, (void **)&psi));
932 PetscCall(DMSwarmRestoreField(user->swarm, "DMSwarm_location_status", NULL, NULL, (void **)&status));
933 PetscCall(DMSwarmRestoreField(user->swarm, "DMSwarm_CellID", NULL, NULL, (void **)&cell_ids));
934 PetscCall(DMSwarmRestoreField(user->swarm, "position", NULL, NULL, (void **)&positions));
935
937 PetscCall(LOG_SCATTER_METRICS(user));
938
939 PetscCall(PetscSNPrintf(metrics_path, sizeof(metrics_path), "%s/scatter_metrics.csv", simCtx->log_dir));
940 PetscCall(PicurvAssertFileExists(metrics_path, "LOG_SCATTER_METRICS should write scatter_metrics.csv"));
941 PetscCall(AssertFileContains(metrics_path, "relative_L2_error",
942 "Scatter metrics CSV header should include relative_L2_error"));
943 PetscCall(AssertFileContains(metrics_path,
944 "3,3.000000e-01,27,27,27,1.000000e+00,1.000000e+00,5.400000e+01,5.400000e+01,0.000000e+00,0.000000e+00,0.000000e+00,0.000000e+00,0.000000e+00",
945 "Scatter metrics CSV should record the expected constant-field zero-error row"));
946
947 PetscCall(PicurvRemoveTempDir(tmpdir));
948 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
949 PetscFunctionReturn(0);
950}
951/**
952 * @brief Tests stdout particle-table logging on a production-like swarm fixture.
953 */
954
955static PetscErrorCode TestParticleFieldTableLogging(void)
956{
957 SimCtx *simCtx = NULL;
958 UserCtx *user = NULL;
959 PetscInt print_interval = 1;
960 char captured[8192];
961
962 PetscFunctionBeginUser;
963 PetscCall(SeedLoggingParticleFixture(&simCtx, &user));
964 PetscCall(CaptureLoggingOutput(user, simCtx, InvokeParticleFieldLog, &print_interval, captured, sizeof(captured)));
965 PetscCall(PicurvAssertBool((PetscBool)(strstr(captured, "Position (x,y,z)") != NULL),
966 "LOG_PARTICLE_FIELDS should print the particle table header"));
967 PetscCall(PicurvAssertBool((PetscBool)(strstr(captured, "Weights (a1,a2,a3)") != NULL),
968 "LOG_PARTICLE_FIELDS should print the weight-column header"));
969 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
970 PetscFunctionReturn(0);
971}
972/**
973 * @brief Tests console snapshot logging against the public periodic-snapshot helper.
974 */
975
976static PetscErrorCode TestParticleConsoleSnapshotLogging(void)
977{
978 SimCtx *simCtx = NULL;
979 UserCtx *user = NULL;
980 PetscInt step = 4;
981 char captured[8192];
982
983 PetscFunctionBeginUser;
984 PetscCall(SeedLoggingParticleFixture(&simCtx, &user));
985 PetscCall(CaptureLoggingOutput(user, simCtx, InvokeParticleConsoleSnapshot, &step, captured, sizeof(captured)));
986 PetscCall(PicurvAssertBool((PetscBool)(strstr(captured, "Particle states at step 4") != NULL),
987 "EmitParticleConsoleSnapshot should print the step banner"));
988 PetscCall(PicurvAssertBool((PetscBool)(strstr(captured, "Position (x,y,z)") != NULL),
989 "EmitParticleConsoleSnapshot should reuse the particle table output"));
990 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
991 PetscFunctionReturn(0);
992}
993/**
994 * @brief Tests file-backed particle metrics logging after derived metrics are computed.
995 */
996
997static PetscErrorCode TestParticleMetricsLogging(void)
998{
999 SimCtx *simCtx = NULL;
1000 UserCtx *user = NULL;
1001 char tmpdir[PETSC_MAX_PATH_LEN];
1002 char metrics_path[PETSC_MAX_PATH_LEN];
1003
1004 PetscFunctionBeginUser;
1005 PetscCall(SeedLoggingParticleFixture(&simCtx, &user));
1006 PetscCall(PicurvMakeTempDir(tmpdir, sizeof(tmpdir)));
1007 PetscCall(PetscStrncpy(simCtx->log_dir, tmpdir, sizeof(simCtx->log_dir)));
1008 simCtx->StartStep = 0;
1009 simCtx->step = 1;
1010 simCtx->particlesLostLastStep = 1;
1011 simCtx->particlesLostCumulative = 7;
1012 simCtx->particlesMigratedLastStep = 2;
1013 simCtx->migrationPassesLastStep = 3;
1014
1015 PetscCall(CalculateParticleCountPerCell(user));
1016 PetscCall(CalculateAdvancedParticleMetrics(user));
1017 PetscCall(LOG_PARTICLE_METRICS(user, "Timestep Metrics"));
1018
1019 PetscCall(PetscSNPrintf(metrics_path, sizeof(metrics_path), "%s/Particle_Metrics.log", simCtx->log_dir));
1020 PetscCall(PicurvAssertFileExists(metrics_path, "LOG_PARTICLE_METRICS should write Particle_Metrics.log"));
1021 PetscCall(AssertFileContains(metrics_path, "Timestep Metrics", "Particle metrics log should include the caller-provided stage label"));
1022 PetscCall(AssertFileContains(metrics_path, "Occupied Cells", "Particle metrics log should include the metrics table header"));
1023 PetscCall(AssertFileContains(metrics_path, "Lost Total", "Particle metrics log should include the cumulative-loss column"));
1024 PetscCall(AssertFileContains(metrics_path, "| 1 | 7 | 2", "Particle metrics log should record both per-step and cumulative loss values"));
1025 PetscCall(PicurvRemoveTempDir(tmpdir));
1026 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
1027 PetscFunctionReturn(0);
1028}
1029/**
1030 * @brief Tests file-backed search metrics logging with the compact CSV contract.
1031 */
1032
1033static PetscErrorCode TestSearchMetricsLogging(void)
1034{
1035 SimCtx *simCtx = NULL;
1036 UserCtx *user = NULL;
1037 char tmpdir[PETSC_MAX_PATH_LEN];
1038 char metrics_path[PETSC_MAX_PATH_LEN];
1039
1040 PetscFunctionBeginUser;
1041 PetscCall(SeedLoggingParticleFixture(&simCtx, &user));
1042 PetscCall(PicurvMakeTempDir(tmpdir, sizeof(tmpdir)));
1043 PetscCall(PetscStrncpy(simCtx->log_dir, tmpdir, sizeof(simCtx->log_dir)));
1044 simCtx->step = 2;
1045 simCtx->ti = 0.2;
1046 simCtx->particlesLostLastStep = 1;
1047 simCtx->particlesLostCumulative = 7;
1048 simCtx->particlesMigratedLastStep = 2;
1049 simCtx->migrationPassesLastStep = 3;
1050 simCtx->particleLoadImbalance = 1.5;
1051 simCtx->searchMetrics.searchAttempts = 4;
1052 simCtx->searchMetrics.searchPopulation = 2;
1054 simCtx->searchMetrics.searchLostCount = 1;
1055 simCtx->searchMetrics.traversalStepsSum = 10;
1056 simCtx->searchMetrics.reSearchCount = 2;
1057 simCtx->searchMetrics.maxTraversalSteps = 6;
1059 simCtx->searchMetrics.tieBreakCount = 1;
1064
1065 PetscCall(LOG_SEARCH_METRICS(user));
1066
1067 PetscCall(PetscSNPrintf(metrics_path, sizeof(metrics_path), "%s/search_metrics.csv", simCtx->log_dir));
1068 PetscCall(PicurvAssertFileExists(metrics_path, "LOG_SEARCH_METRICS should write search_metrics.csv"));
1069 PetscCall(AssertFileContains(metrics_path, "search_attempts", "Search metrics CSV header should include search_attempts"));
1070 PetscCall(AssertFileContains(metrics_path, "search_population", "Search metrics CSV header should include search_population"));
1071 PetscCall(AssertFileContains(metrics_path, "max_particle_pass_depth", "Search metrics CSV header should include max_particle_pass_depth"));
1072 PetscCall(AssertFileContains(metrics_path, "search_work_index", "Search metrics CSV header should include search_work_index"));
1073 PetscCall(AssertFileContains(metrics_path, "re_search_fraction", "Search metrics CSV header should include re_search_fraction"));
1074 PetscCall(AssertFileContains(metrics_path, "lost_cumulative", "Search metrics CSV header should include the cumulative-loss column"));
1075 PetscCall(AssertFileContains(metrics_path, "2,2.000000e-01,2,1,7,2,3,4,2.500000e+00,6,1,2,3,1,3,1.500000e+00,2,1,1,10,2,1,5.000000e-01,5.000000e+00,1.000000e+00", "Search metrics CSV should record the V2 raw and derived search metrics"));
1076 PetscCall(PicurvRemoveTempDir(tmpdir));
1077 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
1078 PetscFunctionReturn(0);
1079}
1080/**
1081 * @brief Tests stdout field-anatomy logging on the corrected production-like DM fixture.
1082 */
1083
1084static PetscErrorCode TestFieldAnatomyLogging(void)
1085{
1086 SimCtx *simCtx = NULL;
1087 UserCtx *user = NULL;
1088 char captured[8192];
1089 AnatomyCaptureCtx anatomy_ctx = {"P", "unit-test"};
1090
1091 PetscFunctionBeginUser;
1092 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
1093 PetscCall(VecSet(user->P, 7.0));
1094 PetscCall(DMGlobalToLocalBegin(user->da, user->P, INSERT_VALUES, user->lP));
1095 PetscCall(DMGlobalToLocalEnd(user->da, user->P, INSERT_VALUES, user->lP));
1096 PetscCall(CaptureLoggingOutput(user, simCtx, InvokeFieldAnatomyLog, &anatomy_ctx, captured, sizeof(captured)));
1097 PetscCall(PicurvAssertBool((PetscBool)(strstr(captured, "Field Anatomy Log: [P]") != NULL),
1098 "LOG_FIELD_ANATOMY should print the requested field name"));
1099 PetscCall(PicurvAssertBool((PetscBool)(strstr(captured, "Layout: [Cell-Centered]") != NULL),
1100 "LOG_FIELD_ANATOMY should report the inferred data layout"));
1101 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
1102 PetscFunctionReturn(0);
1103}
1104/**
1105 * @brief Tests profiling helper lifecycle logging for timestep and final-summary outputs.
1106 */
1107
1108static PetscErrorCode TestProfilingLifecycleHelpers(void)
1109{
1110 SimCtx simCtx;
1111 char tmpdir[PETSC_MAX_PATH_LEN];
1112 char timestep_path[PETSC_MAX_PATH_LEN];
1113 char summary_path[PETSC_MAX_PATH_LEN];
1114 static char selected_name[] = "FlowSolver";
1115 char *selected_funcs[] = {selected_name};
1116
1117 PetscFunctionBeginUser;
1118 PetscCall(PetscMemzero(&simCtx, sizeof(simCtx)));
1119 PetscCall(PicurvMakeTempDir(tmpdir, sizeof(tmpdir)));
1120 PetscCall(PetscStrncpy(simCtx.log_dir, tmpdir, sizeof(simCtx.log_dir)));
1121 PetscCall(PetscStrncpy(simCtx.profilingTimestepMode, "selected", sizeof(simCtx.profilingTimestepMode)));
1122 PetscCall(PetscStrncpy(simCtx.profilingTimestepFile, "Profiling_Timestep_Summary.csv", sizeof(simCtx.profilingTimestepFile)));
1123 simCtx.rank = 0;
1124 simCtx.exec_mode = EXEC_MODE_SOLVER;
1125 simCtx.StartStep = 0;
1126 simCtx.nProfilingSelectedFuncs = 1;
1127 simCtx.profilingSelectedFuncs = selected_funcs;
1128 simCtx.profilingFinalSummary = PETSC_TRUE;
1129
1130 PetscCall(ProfilingInitialize(&simCtx));
1131
1132 _ProfilingStart("FlowSolver");
1133 _ProfilingEnd("FlowSolver");
1134 _ProfilingStart("UnselectedHelper");
1135 _ProfilingEnd("UnselectedHelper");
1136 PetscCall(ProfilingLogTimestepSummary(&simCtx, 1));
1137
1138 _ProfilingStart("FlowSolver");
1139 _ProfilingEnd("FlowSolver");
1140 PetscCall(ProfilingResetTimestepCounters());
1141 PetscCall(ProfilingLogTimestepSummary(&simCtx, 2));
1142 PetscCall(ProfilingFinalize(&simCtx));
1143
1144 PetscCall(PetscSNPrintf(timestep_path, sizeof(timestep_path), "%s/%s", simCtx.log_dir, simCtx.profilingTimestepFile));
1145 PetscCall(PetscSNPrintf(summary_path, sizeof(summary_path), "%s/ProfilingSummary_Solver.log", simCtx.log_dir));
1146 PetscCall(PicurvAssertFileExists(timestep_path, "profiling timestep summary should be written"));
1147 PetscCall(PicurvAssertFileExists(summary_path, "profiling final summary should be written"));
1148 PetscCall(AssertFileContains(timestep_path, "step,function,calls,step_time_s",
1149 "profiling timestep summary should contain the CSV header"));
1150 PetscCall(AssertFileContains(timestep_path, "1,FlowSolver,1,",
1151 "profiling timestep summary should log selected functions"));
1152 PetscCall(AssertFileNotContains(timestep_path, "UnselectedHelper",
1153 "profiling timestep summary should omit unselected functions in selected mode"));
1154 PetscCall(AssertFileContains(summary_path, "FINAL PROFILING SUMMARY",
1155 "profiling final summary should include its table banner"));
1156 PetscCall(AssertFileContains(summary_path, "FlowSolver",
1157 "profiling final summary should include selected functions"));
1158 PetscCall(AssertFileContains(summary_path, "UnselectedHelper",
1159 "profiling final summary should include total-time entries for unselected functions"));
1160 PetscCall(PicurvRemoveTempDir(tmpdir));
1161 PetscFunctionReturn(0);
1162}
1163
1164/**
1165 * @brief Tests runtime memory log header, step rows, final rows, and disabled mode.
1166 */
1167static PetscErrorCode TestRuntimeMemoryLogHelpers(void)
1168{
1169 SimCtx simCtx;
1170 char tmpdir[PETSC_MAX_PATH_LEN];
1171 char memory_path[PETSC_MAX_PATH_LEN];
1172
1173 PetscFunctionBeginUser;
1174 PetscCall(PetscMemzero(&simCtx, sizeof(simCtx)));
1175 PetscCall(PicurvMakeTempDir(tmpdir, sizeof(tmpdir)));
1176 PetscCall(PetscStrncpy(simCtx.log_dir, tmpdir, sizeof(simCtx.log_dir)));
1177 PetscCall(PetscStrncpy(simCtx.runtimeMemoryLogFile, "Runtime_Memory.log", sizeof(simCtx.runtimeMemoryLogFile)));
1178 simCtx.rank = 0;
1179 simCtx.runtimeMemoryLogEnabled = PETSC_TRUE;
1180 simCtx.runtimeMemoryLogStarted = PETSC_FALSE;
1181 simCtx.runtimeMemoryLogHasPrevious = PETSC_FALSE;
1182 simCtx.continueMode = PETSC_FALSE;
1183 simCtx.StartStep = 0;
1184 PetscCall(PetscMemorySetGetMaximumUsage());
1185
1186 PetscCall(RuntimeMemoryLogSample(&simCtx, 1, "Step", "-"));
1187 PetscCall(RuntimeMemoryLogSample(&simCtx, 1, "Final", "Complete"));
1188 PetscCall(PetscSNPrintf(memory_path, sizeof(memory_path), "%s/%s", simCtx.log_dir, simCtx.runtimeMemoryLogFile));
1189 PetscCall(PicurvAssertFileExists(memory_path, "runtime memory log should be written"));
1190 PetscCall(AssertFileContains(memory_path, "Process Current MB Max",
1191 "runtime memory log should contain readable column labels"));
1192 PetscCall(AssertFileContains(memory_path, "Step",
1193 "runtime memory log should contain a step row"));
1194 PetscCall(AssertFileContains(memory_path, "Final",
1195 "runtime memory log should contain a final row"));
1196 PetscCall(AssertFileContains(memory_path, "Complete",
1197 "runtime memory log should record final reason"));
1198
1199 simCtx.runtimeMemoryLogEnabled = PETSC_FALSE;
1200 PetscCall(PetscStrncpy(simCtx.runtimeMemoryLogFile, "Runtime_Memory_Disabled.log", sizeof(simCtx.runtimeMemoryLogFile)));
1201 PetscCall(RuntimeMemoryLogSample(&simCtx, 2, "Step", "-"));
1202 PetscCall(PetscSNPrintf(memory_path, sizeof(memory_path), "%s/%s", simCtx.log_dir, simCtx.runtimeMemoryLogFile));
1203 PetscCall(PicurvAssertBool((PetscBool)(access(memory_path, F_OK) != 0),
1204 "disabled runtime memory log should not create a file"));
1205 PetscCall(PicurvRemoveTempDir(tmpdir));
1206 PetscFunctionReturn(0);
1207}
1208
1209/**
1210 * @brief Tests steady solution-convergence log output, IBM masking, and gauge-invariant pressure drift.
1211 */
1212static PetscErrorCode TestSolutionConvergenceSteadyLogging(void)
1213{
1214 SimCtx *simCtx = NULL;
1215 UserCtx *user = NULL;
1216 char tmpdir[PETSC_MAX_PATH_LEN];
1217 char log_path[PETSC_MAX_PATH_LEN];
1218 char header[4096];
1219 char row1[4096];
1220 char row2[4096];
1221 char mode[128];
1222 Cmpnts ***ucat = NULL;
1223 PetscReal ***nvert = NULL;
1224 PetscReal mean_speed = NAN;
1225 PetscReal mean_speed_ref = NAN;
1226 PetscReal mean_speed_abs = NAN;
1227 PetscReal p_abs = NAN;
1228 PetscInt has_reference_1 = 0;
1229 PetscInt has_reference_2 = 0;
1230
1231 PetscFunctionBeginUser;
1232 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
1233 PetscCall(PicurvMakeTempDir(tmpdir, sizeof(tmpdir)));
1234 PetscCall(PetscStrncpy(simCtx->log_dir, tmpdir, sizeof(simCtx->log_dir)));
1236
1237 PetscCall(InitializeSolutionConvergenceState(simCtx));
1238
1239 simCtx->step = 1;
1240 simCtx->ti = 0.1;
1241 PetscCall(SetUniformVelocityField(user, user->Ucat, 1.0, 0.0, 0.0));
1242 PetscCall(SetUniformScalarField(user, user->P, 11.0));
1243 PetscCall(DMDAVecGetArray(user->da, user->Nvert, &nvert));
1244 nvert[1][1][1] = 1.0;
1245 PetscCall(DMDAVecRestoreArray(user->da, user->Nvert, &nvert));
1246 PetscCall(DMDAVecGetArray(user->fda, user->Ucat, &ucat));
1247 ucat[1][1][1].x = 999.0;
1248 PetscCall(DMDAVecRestoreArray(user->fda, user->Ucat, &ucat));
1249 PetscCall(LOG_SOLUTION_CONVERGENCE(simCtx));
1250
1251 simCtx->step = 2;
1252 simCtx->ti = 0.2;
1253 PetscCall(SetUniformVelocityField(user, user->Ucat, 1.0, 0.0, 0.0));
1254 PetscCall(SetUniformVelocityField(user, user->Ucat_o, 0.5, 0.0, 0.0));
1255 PetscCall(SetUniformScalarField(user, user->P, 11.0));
1256 PetscCall(SetUniformScalarField(user, user->P_o, 7.0));
1257 PetscCall(DMDAVecGetArray(user->fda, user->Ucat, &ucat));
1258 ucat[1][1][1].x = 999.0;
1259 PetscCall(DMDAVecRestoreArray(user->fda, user->Ucat, &ucat));
1260 PetscCall(LOG_SOLUTION_CONVERGENCE(simCtx));
1261
1262 PetscCall(PetscSNPrintf(log_path, sizeof(log_path), "%s/solution_convergence.log", simCtx->log_dir));
1263 PetscCall(PicurvAssertFileExists(log_path, "steady solution-convergence logging should write the log"));
1264 PetscCall(ReadLogHeaderAndRow(log_path, 1, header, sizeof(header), row1, sizeof(row1)));
1265 PetscCall(ReadLogHeaderAndRow(log_path, 2, header, sizeof(header), row2, sizeof(row2)));
1266 PetscCall(LogGetColumnText(header, row1, "mode", mode, sizeof(mode)));
1267 PetscCall(LogGetColumnInt(header, row1, "ref", &has_reference_1));
1268 PetscCall(LogGetColumnInt(header, row2, "ref", &has_reference_2));
1269 PetscCall(LogGetColumnReal(header, row2, "mean_speed", &mean_speed));
1270 PetscCall(LogGetColumnReal(header, row2, "spd_ref", &mean_speed_ref));
1271 PetscCall(LogGetColumnReal(header, row2, "spd_abs", &mean_speed_abs));
1272 PetscCall(LogGetColumnReal(header, row2, "p_abs_l2", &p_abs));
1273
1274 PetscCall(PicurvAssertBool((PetscBool)(strcmp(mode, "steady_deterministic") == 0),
1275 "steady solution convergence row should record the mode name"));
1276 PetscCall(PicurvAssertIntEqual(0, has_reference_1,
1277 "the first steady solution-convergence row should be a warmup row"));
1278 PetscCall(PicurvAssertIntEqual(1, has_reference_2,
1279 "steady solution convergence should compare against the previous solved step"));
1280 PetscCall(PicurvAssertRealNear(1.0, mean_speed, 1.0e-12,
1281 "steady solution convergence should mask IBM-marked solid cells"));
1282 PetscCall(PicurvAssertRealNear(0.5, mean_speed_ref, 1.0e-12,
1283 "steady solution convergence should report the previous-step mean speed"));
1284 PetscCall(PicurvAssertRealNear(0.5, mean_speed_abs, 1.0e-12,
1285 "steady solution convergence should report mean-speed drift"));
1286 PetscCall(PicurvAssertRealNear(0.0, p_abs, 1.0e-12,
1287 "steady solution convergence pressure drift should be gauge invariant"));
1288
1289 PetscCall(PicurvRemoveTempDir(tmpdir));
1290 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
1291 PetscFunctionReturn(0);
1292}
1293
1294/**
1295 * @brief Tests periodic solution-convergence warmup and phase-aligned reference reuse.
1296 */
1298{
1299 SimCtx *simCtx = NULL;
1300 UserCtx *user = NULL;
1301 char tmpdir[PETSC_MAX_PATH_LEN];
1302 char log_path[PETSC_MAX_PATH_LEN];
1303 char header[4096];
1304 char row1[4096];
1305 char row2[4096];
1306 char row3[4096];
1307 PetscInt has_reference_1 = 0;
1308 PetscInt has_reference_2 = 0;
1309 PetscInt has_reference_3 = 0;
1310 PetscInt phase_step_1 = -1;
1311 PetscInt phase_step_2 = -1;
1312 PetscInt phase_step_3 = -1;
1313 PetscReal mean_speed_ref_2 = NAN;
1314 PetscReal mean_speed_abs_2 = NAN;
1315
1316 PetscFunctionBeginUser;
1317 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
1318 PetscCall(PicurvMakeTempDir(tmpdir, sizeof(tmpdir)));
1319 PetscCall(PetscStrncpy(simCtx->log_dir, tmpdir, sizeof(simCtx->log_dir)));
1322
1323 PetscCall(InitializeSolutionConvergenceState(simCtx));
1324
1325 simCtx->step = 1;
1326 simCtx->ti = 0.1;
1327 PetscCall(SetUniformVelocityField(user, user->Ucat, 1.0, 0.0, 0.0));
1328 PetscCall(SetUniformScalarField(user, user->P, 2.0));
1329 PetscCall(LOG_SOLUTION_CONVERGENCE(simCtx));
1330
1331 simCtx->step = 2;
1332 simCtx->ti = 0.2;
1333 PetscCall(SetUniformVelocityField(user, user->Ucat, 2.0, 0.0, 0.0));
1334 PetscCall(SetUniformScalarField(user, user->P, 5.0));
1335 PetscCall(LOG_SOLUTION_CONVERGENCE(simCtx));
1336
1337 simCtx->step = 3;
1338 simCtx->ti = 0.3;
1339 PetscCall(SetUniformVelocityField(user, user->Ucat, 1.25, 0.0, 0.0));
1340 PetscCall(SetUniformScalarField(user, user->P, 9.0));
1341 PetscCall(LOG_SOLUTION_CONVERGENCE(simCtx));
1342
1343 PetscCall(PetscSNPrintf(log_path, sizeof(log_path), "%s/solution_convergence.log", simCtx->log_dir));
1344 PetscCall(PicurvAssertFileExists(log_path, "periodic solution-convergence logging should write the log"));
1345 PetscCall(ReadLogHeaderAndRow(log_path, 1, header, sizeof(header), row1, sizeof(row1)));
1346 PetscCall(ReadLogHeaderAndRow(log_path, 2, header, sizeof(header), row2, sizeof(row2)));
1347 PetscCall(ReadLogHeaderAndRow(log_path, 3, header, sizeof(header), row3, sizeof(row3)));
1348 PetscCall(LogGetColumnInt(header, row1, "ref", &has_reference_1));
1349 PetscCall(LogGetColumnInt(header, row2, "ref", &has_reference_2));
1350 PetscCall(LogGetColumnInt(header, row3, "ref", &has_reference_3));
1351 PetscCall(LogGetColumnInt(header, row1, "ph", &phase_step_1));
1352 PetscCall(LogGetColumnInt(header, row2, "ph", &phase_step_2));
1353 PetscCall(LogGetColumnInt(header, row3, "ph", &phase_step_3));
1354 PetscCall(LogGetColumnReal(header, row3, "spd_ref", &mean_speed_ref_2));
1355 PetscCall(LogGetColumnReal(header, row3, "spd_abs", &mean_speed_abs_2));
1356
1357 PetscCall(PicurvAssertIntEqual(0, has_reference_1,
1358 "the first periodic phase visit should log warmup without a reference"));
1359 PetscCall(PicurvAssertIntEqual(0, has_reference_2,
1360 "the first periodic cycle should fully warm up before comparisons begin"));
1361 PetscCall(PicurvAssertIntEqual(1, has_reference_3,
1362 "the repeated periodic phase visit should compare against the stored reference"));
1363 PetscCall(PicurvAssertIntEqual(1, phase_step_1,
1364 "periodic solution convergence should log the current phase slot"));
1365 PetscCall(PicurvAssertIntEqual(0, phase_step_2,
1366 "periodic solution convergence should log distinct phase slots during warmup"));
1367 PetscCall(PicurvAssertIntEqual(1, phase_step_3,
1368 "periodic solution convergence should reuse the same phase slot on later cycles"));
1369 PetscCall(PicurvAssertRealNear(1.0, mean_speed_ref_2, 1.0e-12,
1370 "periodic solution convergence should report the stored phase-aligned reference"));
1371 PetscCall(PicurvAssertRealNear(0.25, mean_speed_abs_2, 1.0e-12,
1372 "periodic solution convergence should report the phase-aligned drift"));
1373
1374 PetscCall(PicurvRemoveTempDir(tmpdir));
1375 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
1376 PetscFunctionReturn(0);
1377}
1378
1379/**
1380 * @brief Tests statistical solution-convergence sliding-window metrics.
1381 */
1383{
1384 SimCtx *simCtx = NULL;
1385 UserCtx *user = NULL;
1386 char tmpdir[PETSC_MAX_PATH_LEN];
1387 char log_path[PETSC_MAX_PATH_LEN];
1388 char header[4096];
1389 char row[4096];
1390 PetscReal mean_speed_window = NAN;
1391 PetscReal mean_speed_window_prev = NAN;
1392 PetscReal mean_speed_window_abs = NAN;
1393 PetscReal mean_speed_rms_window = NAN;
1394 PetscReal mean_ke_window = NAN;
1395 PetscReal mean_ke_window_prev = NAN;
1396 PetscReal mean_ke_window_abs = NAN;
1397 PetscReal mean_ke_rms_window_abs = NAN;
1398 PetscInt has_reference = 0;
1399 const PetscReal speed_samples[4] = {1.0, 3.0, 2.0, 4.0};
1400
1401 PetscFunctionBeginUser;
1402 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
1403 PetscCall(PicurvMakeTempDir(tmpdir, sizeof(tmpdir)));
1404 PetscCall(PetscStrncpy(simCtx->log_dir, tmpdir, sizeof(simCtx->log_dir)));
1407
1408 PetscCall(InitializeSolutionConvergenceState(simCtx));
1409 for (PetscInt step = 0; step < 4; ++step) {
1410 simCtx->step = step + 1;
1411 simCtx->ti = 0.1 * (PetscReal)(step + 1);
1412 PetscCall(SetUniformVelocityField(user, user->Ucat, speed_samples[step], 0.0, 0.0));
1413 PetscCall(LOG_SOLUTION_CONVERGENCE(simCtx));
1414 }
1415
1416 PetscCall(PetscSNPrintf(log_path, sizeof(log_path), "%s/solution_convergence.log", simCtx->log_dir));
1417 PetscCall(PicurvAssertFileExists(log_path, "statistical solution-convergence logging should write the log"));
1418 PetscCall(ReadLogHeaderAndRow(log_path, 4, header, sizeof(header), row, sizeof(row)));
1419 PetscCall(LogGetColumnInt(header, row, "ref", &has_reference));
1420 PetscCall(LogGetColumnReal(header, row, "spd_win", &mean_speed_window));
1421 PetscCall(LogGetColumnReal(header, row, "spd_win_prev", &mean_speed_window_prev));
1422 PetscCall(LogGetColumnReal(header, row, "spd_win_abs", &mean_speed_window_abs));
1423 PetscCall(LogGetColumnReal(header, row, "spd_rms_win", &mean_speed_rms_window));
1424 PetscCall(LogGetColumnReal(header, row, "ke_win", &mean_ke_window));
1425 PetscCall(LogGetColumnReal(header, row, "ke_win_prev", &mean_ke_window_prev));
1426 PetscCall(LogGetColumnReal(header, row, "ke_win_abs", &mean_ke_window_abs));
1427 PetscCall(LogGetColumnReal(header, row, "ke_rms_abs", &mean_ke_rms_window_abs));
1428
1429 PetscCall(PicurvAssertIntEqual(1, has_reference,
1430 "statistical solution convergence should emit adjacent-window drift once two windows exist"));
1431 PetscCall(PicurvAssertRealNear(3.0, mean_speed_window, 1.0e-12,
1432 "statistical solution convergence should report the current mean-speed window"));
1433 PetscCall(PicurvAssertRealNear(2.0, mean_speed_window_prev, 1.0e-12,
1434 "statistical solution convergence should report the previous mean-speed window"));
1435 PetscCall(PicurvAssertRealNear(1.0, mean_speed_window_abs, 1.0e-12,
1436 "statistical solution convergence should report mean-speed window drift"));
1437 PetscCall(PicurvAssertRealNear(1.0, mean_speed_rms_window, 1.0e-12,
1438 "statistical solution convergence should report current mean-speed RMS"));
1439 PetscCall(PicurvAssertRealNear(5.0, mean_ke_window, 1.0e-12,
1440 "statistical solution convergence should report the current kinetic-energy window"));
1441 PetscCall(PicurvAssertRealNear(2.5, mean_ke_window_prev, 1.0e-12,
1442 "statistical solution convergence should report the previous kinetic-energy window"));
1443 PetscCall(PicurvAssertRealNear(2.5, mean_ke_window_abs, 1.0e-12,
1444 "statistical solution convergence should report kinetic-energy window drift"));
1445 PetscCall(PicurvAssertRealNear(1.0, mean_ke_rms_window_abs, 1.0e-12,
1446 "statistical solution convergence should report RMS kinetic-energy drift"));
1447
1448 PetscCall(PicurvRemoveTempDir(tmpdir));
1449 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
1450 PetscFunctionReturn(0);
1451}
1452/**
1453 * @brief Runs the unit-logging PETSc test binary.
1454 */
1455
1456int main(int argc, char **argv)
1457{
1458 PetscErrorCode ierr;
1459 const PicurvTestCase cases[] = {
1460 {"string-conversion-helpers", TestStringConversionHelpers},
1461 {"get-log-level-from-environment", TestGetLogLevelFromEnvironment},
1462 {"allowed-functions-filter", TestAllowedFunctionsFilter},
1463 {"particle-console-snapshot-cadence", TestParticleConsoleSnapshotCadence},
1464 {"logging-file-parsing-and-formatting-helpers", TestLoggingFileParsingAndFormattingHelpers},
1465 {"logging-continuity-and-field-diagnostics", TestLoggingContinuityAndFieldDiagnostics},
1466 {"interpolation-error-logging", TestInterpolationErrorLogging},
1467 {"scatter-metrics-logging", TestScatterMetricsLogging},
1468 {"particle-field-table-logging", TestParticleFieldTableLogging},
1469 {"particle-console-snapshot-logging", TestParticleConsoleSnapshotLogging},
1470 {"particle-metrics-logging", TestParticleMetricsLogging},
1471 {"search-metrics-logging", TestSearchMetricsLogging},
1472 {"field-anatomy-logging", TestFieldAnatomyLogging},
1473 {"profiling-lifecycle-helpers", TestProfilingLifecycleHelpers},
1474 {"runtime-memory-log-helpers", TestRuntimeMemoryLogHelpers},
1475 {"solution-convergence-steady-logging", TestSolutionConvergenceSteadyLogging},
1476 {"solution-convergence-periodic-logging", TestSolutionConvergencePeriodicLogging},
1477 {"solution-convergence-statistical-logging", TestSolutionConvergenceStatisticalLogging},
1478 };
1479
1480 (void)setenv("LOG_LEVEL", "INFO", 1);
1481
1482 ierr = PetscInitialize(&argc, &argv, NULL, "PICurv logging tests");
1483 if (ierr) {
1484 return (int)ierr;
1485 }
1486
1487 ierr = PicurvRunTests("unit-logging", cases, sizeof(cases) / sizeof(cases[0]));
1488 if (ierr) {
1489 PetscFinalize();
1490 return (int)ierr;
1491 }
1492
1493 set_allowed_functions(NULL, 0);
1494
1495 ierr = PetscFinalize();
1496 return (int)ierr;
1497}
PetscErrorCode SetAnalyticalSolutionForParticles(Vec tempVec, SimCtx *simCtx)
Applies the analytical solution to particle velocity vector.
PetscErrorCode CalculateParticleCountPerCell(UserCtx *user)
Counts particles in each cell of the DMDA 'da' and stores the result in user->ParticleCount.
PetscErrorCode ScatterAllParticleFieldsToEulerFields(UserCtx *user)
Scatters a predefined set of particle fields to their corresponding Eulerian fields.
Logging utilities and macros for PETSc-based applications.
void set_allowed_functions(const char **functionList, int count)
Sets the global list of function names that are allowed to log.
Definition logging.c:152
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:3182
const char * BCHandlerTypeToString(BCHandlerType handler_type)
Converts a BCHandlerType enum to its string representation.
Definition logging.c:771
PetscBool is_function_allowed(const char *functionName)
Checks if a given function is in the allow-list.
Definition logging.c:183
const char * FieldInitializationToString(PetscInt FieldInitialization)
Helper function to convert FieldInitialization to a string representation.
Definition logging.c:687
PetscErrorCode DualMonitorDestroy(void **ctx)
Destroys the DualMonitorCtx.
Definition logging.c:809
PetscErrorCode LOG_INTERPOLATION_ERROR(UserCtx *user)
Logs the interpolation error between the analytical and computed solutions.
Definition logging.c:2719
PetscBool ShouldEmitPeriodicParticleConsoleSnapshot(const SimCtx *simCtx, PetscInt completed_step)
Returns whether a particle console snapshot should be emitted for the.
Definition logging.c:542
const char * BCFaceToString(BCFace face)
Helper function to convert BCFace enum to a string representation.
Definition logging.c:669
PetscErrorCode FreeAllowedFunctions(char **funcs, PetscInt n)
Free an array previously returned by LoadAllowedFunctionsFromFile().
Definition logging.c:650
PetscBool IsParticleConsoleSnapshotEnabled(const SimCtx *simCtx)
Returns whether periodic particle console snapshots are enabled.
Definition logging.c:525
PetscErrorCode print_log_level(void)
Prints the current logging level to the console.
Definition logging.c:116
PetscErrorCode EmitParticleConsoleSnapshot(UserCtx *user, SimCtx *simCtx, PetscInt step)
Emits one particle console snapshot into the main solver log.
Definition logging.c:556
PetscErrorCode ProfilingFinalize(SimCtx *simCtx)
the profiling excercise and build a profiling summary which is then printed to a log file.
Definition logging.c:2148
PetscErrorCode LoadAllowedFunctionsFromFile(const char filename[], char ***funcsOut, PetscInt *nOut)
Load function names from a text file.
Definition logging.c:596
PetscErrorCode LOG_FIELD_MIN_MAX(UserCtx *user, const char *fieldName)
Computes and logs the local and global min/max values of a 3-component vector field.
Definition logging.c:2301
void PrintProgressBar(PetscInt step, PetscInt startStep, PetscInt totalSteps, PetscReal currentTime)
Prints a progress bar to the console.
Definition logging.c:2254
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:2449
PetscErrorCode RuntimeMemoryLogSample(SimCtx *simCtx, PetscInt step, const char *event, const char *reason)
Append a reduced runtime memory sample to the configured memory log.
Definition logging.c:2037
LogLevel get_log_level()
Retrieves the current logging level from the environment variable LOG_LEVEL.
Definition logging.c:84
PetscErrorCode ProfilingLogTimestepSummary(SimCtx *simCtx, PetscInt step)
Logs the performance summary for the current timestep and resets timers.
Definition logging.c:1959
PetscErrorCode LOG_FACE_DISTANCES(PetscReal *d)
Prints the signed distances to each face of the cell.
Definition logging.c:230
PetscErrorCode LOG_PARTICLE_FIELDS(UserCtx *user, PetscInt printInterval)
Prints particle fields in a table that automatically adjusts its column widths.
Definition logging.c:397
void _ProfilingEnd(const char *func_name)
Internal profiling hook invoked by PROFILE_FUNCTION_END.
Definition logging.c:1920
const char * BCTypeToString(BCType type)
Helper function to convert BCType enum to a string representation.
Definition logging.c:751
PetscErrorCode CalculateAdvancedParticleMetrics(UserCtx *user)
Computes advanced particle statistics and stores them in SimCtx.
Definition logging.c:3128
const char * ParticleLocationStatusToString(ParticleLocationStatus level)
A function that outputs the name of the current level in the ParticleLocation enum.
Definition logging.c:1811
PetscErrorCode LOG_SCATTER_METRICS(UserCtx *user)
Logs particle-to-grid scatter verification metrics for the prescribed scalar truth path.
Definition logging.c:2795
PetscErrorCode LOG_SOLUTION_CONVERGENCE(SimCtx *simCtx)
Logs physical solution-convergence metrics once per completed timestep.
Definition logging.c:1560
PetscErrorCode LOG_CONTINUITY_METRICS(UserCtx *user)
Logs continuity metrics for a single block to a file.
Definition logging.c:1752
PetscErrorCode LOG_SEARCH_METRICS(UserCtx *user)
Writes compact runtime search metrics to CSV and optionally to console.
Definition logging.c:2977
PetscErrorCode ProfilingInitialize(SimCtx *simCtx)
Initializes the custom profiling system using configuration from SimCtx.
Definition logging.c:1882
@ LOG_INFO
Informational messages about program execution.
Definition logging.h:30
const char * LESModelToString(LESModelType LESFlag)
Helper function to convert LES Flag to a string representation.
Definition logging.c:720
PetscErrorCode LOG_CELL_VERTICES(const Cell *cell, PetscMPIInt rank)
Prints the coordinates of a cell's vertices.
Definition logging.c:205
PetscErrorCode ProfilingResetTimestepCounters(void)
Resets per-timestep profiling counters for the next solver step.
Definition logging.c:1942
const char * MomentumSolverTypeToString(MomentumSolverType SolverFlag)
Helper function to convert Momentum Solver flag to a string representation.
Definition logging.c:736
FILE * file_handle
Definition logging.h:56
const char * ParticleInitializationToString(ParticleInitializationType ParticleInitialization)
Helper function to convert ParticleInitialization to a string representation.
Definition logging.c:703
void _ProfilingStart(const char *func_name)
Internal profiling hook invoked by PROFILE_FUNCTION_BEGIN.
Definition logging.c:1906
Context for a dual-purpose KSP monitor.
Definition logging.h:55
PetscErrorCode InitializeSolutionConvergenceState(SimCtx *simCtx)
Allocates any runtime storage required by solution-convergence logging.
Definition setup.c:47
static PetscErrorCode TestParticleConsoleSnapshotCadence(void)
Tests periodic particle console snapshot enablement and cadence.
static PetscErrorCode SetUniformScalarField(UserCtx *user, Vec field, PetscReal value)
Fills one scalar field with a uniform constant state.
static PetscErrorCode TestFieldAnatomyLogging(void)
Tests stdout field-anatomy logging on the corrected production-like DM fixture.
static PetscErrorCode LogGetColumnText(const char *header, const char *row, const char *column_name, char *value, size_t value_len)
Extracts one pipe-delimited log cell as text by header name.
static PetscErrorCode InvokeParticleConsoleSnapshot(UserCtx *user, SimCtx *simCtx, void *ctx)
Adapts EmitParticleConsoleSnapshot() to the generic stdout-capture callback shape.
static PetscErrorCode TestParticleConsoleSnapshotLogging(void)
Tests console snapshot logging against the public periodic-snapshot helper.
static PetscErrorCode CsvGetColumnText(const char *header, const char *row, const char *column_name, char *value, size_t value_len)
Extracts one CSV cell as text by header name.
static PetscErrorCode CaptureLoggingOutput(UserCtx *user, SimCtx *simCtx, CapturedLoggingFn fn, void *ctx, char *captured, size_t captured_len)
Captures stdout emitted by one logging helper into a temporary file-backed buffer.
static PetscErrorCode ReadCsvHeaderAndRow(const char *path, PetscInt row_index, char *header, size_t header_len, char *row, size_t row_len)
Reads the CSV header and the requested 1-based data row from a text file.
static PetscErrorCode TestGetLogLevelFromEnvironment(void)
Tests that log level selection honors the environment variable.
int main(int argc, char **argv)
Runs the unit-logging PETSc test binary.
static PetscErrorCode TestLoggingContinuityAndFieldDiagnostics(void)
Tests continuity, min/max, and anatomy logging helpers on minimal runtime fixtures.
PetscErrorCode(* CapturedLoggingFn)(UserCtx *user, SimCtx *simCtx, void *ctx)
static PetscErrorCode TestSolutionConvergenceStatisticalLogging(void)
Tests statistical solution-convergence sliding-window metrics.
static PetscErrorCode AssertFileNotContains(const char *path, const char *needle, const char *context)
Asserts that one text file does not contain an excluded substring.
static PetscErrorCode InvokeFieldAnatomyLog(UserCtx *user, SimCtx *simCtx, void *ctx)
Adapts LOG_FIELD_ANATOMY() to the generic stdout-capture callback shape.
static PetscErrorCode TestSolutionConvergencePeriodicLogging(void)
Tests periodic solution-convergence warmup and phase-aligned reference reuse.
static PetscErrorCode TestStringConversionHelpers(void)
Tests string-conversion helpers for configured enums and unknown values.
static PetscErrorCode TestLoggingFileParsingAndFormattingHelpers(void)
Tests logging-side file parsing, helper formatting, and progress utilities.
static PetscErrorCode TestAllowedFunctionsFilter(void)
Tests the function allow-list filter used by the logging layer.
static PetscErrorCode TestInterpolationErrorLogging(void)
Tests interpolation-error logging against an analytically matched particle field.
static PetscErrorCode SeedLoggingParticleFixture(SimCtx **simCtx_out, UserCtx **user_out)
Creates a small particle-bearing runtime fixture used by logging tests.
static PetscErrorCode TestScatterMetricsLogging(void)
Tests file-backed scatter metrics logging against a fully occupied constant field.
static PetscErrorCode TestParticleMetricsLogging(void)
Tests file-backed particle metrics logging after derived metrics are computed.
static PetscErrorCode CsvGetColumnReal(const char *header, const char *row, const char *column_name, PetscReal *value_out)
Extracts one CSV cell as a real by header name.
static PetscErrorCode TestSearchMetricsLogging(void)
Tests file-backed search metrics logging with the compact CSV contract.
static PetscErrorCode TestParticleFieldTableLogging(void)
Tests stdout particle-table logging on a production-like swarm fixture.
static PetscErrorCode CsvFindColumnIndex(const char *header, const char *column_name, PetscInt *index_out)
Returns the zero-based column index for one CSV header field.
static PetscErrorCode LogFindColumnIndex(const char *header, const char *column_name, PetscInt *index_out)
Returns the zero-based column index for one pipe-delimited header field.
static PetscErrorCode SetUniformVelocityField(UserCtx *user, Vec field, PetscReal ux, PetscReal uy, PetscReal uz)
Fills one Cartesian velocity field with a uniform constant state.
static PetscErrorCode ReadLogHeaderAndRow(const char *path, PetscInt row_index, char *header, size_t header_len, char *row, size_t row_len)
Reads the column header and the requested 1-based data row from a pipe-delimited solution-convergence...
static PetscErrorCode AssertFileContains(const char *path, const char *needle, const char *context)
Asserts that one text file contains a required substring.
static PetscErrorCode CsvGetColumnInt(const char *header, const char *row, const char *column_name, PetscInt *value_out)
Extracts one CSV cell as an integer by header name.
static PetscErrorCode LogGetColumnReal(const char *header, const char *row, const char *column_name, PetscReal *value_out)
Extracts one pipe-delimited log cell as a real by header name.
static PetscErrorCode TestSolutionConvergenceSteadyLogging(void)
Tests steady solution-convergence log output, IBM masking, and gauge-invariant pressure drift.
static PetscErrorCode TestProfilingLifecycleHelpers(void)
Tests profiling helper lifecycle logging for timestep and final-summary outputs.
static PetscErrorCode LogGetColumnInt(const char *header, const char *row, const char *column_name, PetscInt *value_out)
Extracts one pipe-delimited log cell as an integer by header name.
const char * stage_name
const char * field_name
static PetscErrorCode InvokeParticleFieldLog(UserCtx *user, SimCtx *simCtx, void *ctx)
Adapts LOG_PARTICLE_FIELDS() to the generic stdout-capture callback shape.
static PetscErrorCode TestRuntimeMemoryLogHelpers(void)
Tests runtime memory log header, step rows, final rows, and disabled mode.
PetscErrorCode PicurvMakeTempDir(char *path, size_t path_len)
Creates a unique temporary directory for one test case.
PetscErrorCode PicurvCreateMinimalContexts(SimCtx **simCtx_out, UserCtx **user_out, PetscInt mx, PetscInt my, PetscInt mz)
Builds minimal SimCtx and UserCtx fixtures for C unit tests.
PetscErrorCode PicurvAssertRealNear(PetscReal expected, PetscReal actual, PetscReal tol, const char *context)
Asserts that two real values agree within tolerance.
PetscErrorCode PicurvDestroyMinimalContexts(SimCtx **simCtx_ptr, UserCtx **user_ptr)
Destroys minimal SimCtx/UserCtx fixtures and all owned PETSc objects.
PetscErrorCode PicurvCreateSwarmPair(UserCtx *user, PetscInt nlocal, const char *post_field_name)
Creates matched solver and post-processing swarms for tests.
PetscErrorCode PicurvRunTests(const char *suite_name, const PicurvTestCase *cases, size_t case_count)
Runs a named C test suite and prints pass/fail progress markers.
PetscErrorCode PicurvAssertFileExists(const char *path, const char *context)
Asserts that a filesystem path exists as a readable file.
PetscErrorCode PicurvAssertIntEqual(PetscInt expected, PetscInt actual, const char *context)
Asserts that two integer values are equal.
PetscErrorCode PicurvAssertBool(PetscBool value, const char *context)
Asserts that one boolean condition is true.
PetscErrorCode PicurvRemoveTempDir(const char *path)
Recursively removes a temporary directory created by PicurvMakeTempDir.
Shared declarations for the PICurv C test fixture and assertion layer.
Named test case descriptor consumed by PicurvRunTests.
@ CONSTANT_SMAGORINSKY
Definition variables.h:490
@ PERIODIC
Definition variables.h:260
PetscBool continueMode
Definition variables.h:668
PetscBool profilingFinalSummary
Definition variables.h:794
PetscMPIInt rank
Definition variables.h:654
char profilingTimestepFile[PETSC_MAX_PATH_LEN]
Definition variables.h:793
PetscInt64 searchLocatedCount
Definition variables.h:224
PetscInt64 searchLostCount
Definition variables.h:225
@ PARTICLE_INIT_VOLUME
Random volumetric distribution across the domain.
Definition variables.h:518
@ LOST
Definition variables.h:139
@ ACTIVE_AND_LOCATED
Definition variables.h:137
PetscReal FluxOutSum
Definition variables.h:735
PetscBool runtimeMemoryLogEnabled
Enable the rank-reduced runtime memory log.
Definition variables.h:809
PetscInt64 boundaryClampCount
Definition variables.h:231
PetscInt particlesLostLastStep
Definition variables.h:760
PetscInt64 traversalStepsSum
Definition variables.h:226
@ BC_HANDLER_PERIODIC_DRIVEN_CONSTANT_FLUX
Definition variables.h:286
PetscReal ren
Definition variables.h:699
PetscInt64 searchPopulation
Definition variables.h:223
char runtimeMemoryLogFile[PETSC_MAX_PATH_LEN]
File name written under log_dir.
Definition variables.h:810
PetscBool runtimeMemoryLogStarted
True after rank 0 writes the log header.
Definition variables.h:811
char profilingTimestepMode[32]
Definition variables.h:792
PetscInt np
Definition variables.h:753
Vec Ucont
Definition variables.h:856
PetscInt StartStep
Definition variables.h:661
@ MOMENTUM_SOLVER_EXPLICIT_RK
Definition variables.h:503
PetscInt solutionConvergencePeriodSteps
Definition variables.h:708
PetscScalar x
Definition variables.h:101
PetscInt64 reSearchCount
Definition variables.h:227
PetscReal MaxDiv
Definition variables.h:785
PetscInt64 bboxGuessFallbackCount
Definition variables.h:233
VerificationScalarConfig verificationScalar
Definition variables.h:714
Vec Ucat_o
Definition variables.h:863
PetscInt MaxDivx
Definition variables.h:786
PetscInt MaxDivy
Definition variables.h:786
PetscInt64 bboxGuessSuccessCount
Definition variables.h:232
PetscInt MaxDivz
Definition variables.h:786
char log_dir[PETSC_MAX_PATH_LEN]
Definition variables.h:676
PetscInt MaxDivFlatArg
Definition variables.h:786
PetscReal FluxInSum
Definition variables.h:735
PetscInt64 maxParticlePassDepth
Definition variables.h:234
PetscInt64 maxTraversalSteps
Definition variables.h:228
PetscScalar z
Definition variables.h:101
Vec Ucat
Definition variables.h:856
PetscBool runtimeMemoryLogHasPrevious
True after the first process-memory sample.
Definition variables.h:812
char ** profilingSelectedFuncs
Definition variables.h:790
PetscInt solutionConvergenceWindowSteps
Definition variables.h:709
PetscInt particlesLostCumulative
Definition variables.h:761
PetscInt nProfilingSelectedFuncs
Definition variables.h:791
PetscInt particlesMigratedLastStep
Definition variables.h:763
char AnalyticalSolutionType[PETSC_MAX_PATH_LEN]
Definition variables.h:684
PetscInt particleConsoleOutputFreq
Definition variables.h:664
SearchMetricsState searchMetrics
Definition variables.h:766
Vec lUcont
Definition variables.h:856
PetscInt step
Definition variables.h:659
DMDALocalInfo info
Definition variables.h:837
Vec lUcat
Definition variables.h:856
PetscInt migrationPassesLastStep
Definition variables.h:762
PetscScalar y
Definition variables.h:101
@ EXEC_MODE_SOLVER
Definition variables.h:624
Vec Nvert
Definition variables.h:856
@ SOLUTION_CONVERGENCE_PERIODIC_DETERMINISTIC
Definition variables.h:510
@ SOLUTION_CONVERGENCE_STATISTICAL_STEADY
Definition variables.h:511
@ SOLUTION_CONVERGENCE_STEADY_DETERMINISTIC
Definition variables.h:509
SolutionConvergenceMode solutionConvergenceMode
Definition variables.h:707
PetscInt64 searchAttempts
Definition variables.h:222
ExecutionMode exec_mode
Definition variables.h:670
PetscInt64 tieBreakCount
Definition variables.h:230
PetscReal ti
Definition variables.h:660
PetscReal summationRHS
Definition variables.h:784
PetscInt64 maxTraversalFailCount
Definition variables.h:229
Cmpnts vertices[8]
Coordinates of the eight vertices of the cell.
Definition variables.h:161
PetscReal particleLoadImbalance
Definition variables.h:765
Vec P_o
Definition variables.h:863
@ BC_FACE_NEG_X
Definition variables.h:245
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:651
User-defined context containing data specific to a single computational grid level.
Definition variables.h:830