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
11#include <fcntl.h>
12#include <stdio.h>
13#include <stdlib.h>
14#include <string.h>
15#include <unistd.h>
16/**
17 * @brief Asserts that one text file contains a required substring.
18 */
19
20static PetscErrorCode AssertFileContains(const char *path, const char *needle, const char *context)
21{
22 FILE *fp = NULL;
23 long file_size = 0;
24 char *buffer = NULL;
25
26 PetscFunctionBeginUser;
27 fp = fopen(path, "rb");
28 if (!fp) {
29 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Failed to open '%s' for assertion.", path);
30 }
31 if (fseek(fp, 0, SEEK_END) != 0) {
32 fclose(fp);
33 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Failed to seek '%s'.", path);
34 }
35 file_size = ftell(fp);
36 if (file_size < 0) {
37 fclose(fp);
38 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Failed to measure '%s'.", path);
39 }
40 if (fseek(fp, 0, SEEK_SET) != 0) {
41 fclose(fp);
42 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Failed to rewind '%s'.", path);
43 }
44
45 PetscCall(PetscMalloc1((size_t)file_size + 1, &buffer));
46 if (file_size > 0) {
47 size_t bytes_read = fread(buffer, 1, (size_t)file_size, fp);
48 if (bytes_read != (size_t)file_size) {
49 fclose(fp);
50 PetscCall(PetscFree(buffer));
51 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Failed to read '%s'.", path);
52 }
53 }
54 buffer[file_size] = '\0';
55 fclose(fp);
56
57 PetscCall(PicurvAssertBool((PetscBool)(strstr(buffer, needle) != NULL), context));
58 PetscCall(PetscFree(buffer));
59 PetscFunctionReturn(0);
60}
61/**
62 * @brief Asserts that one text file does not contain an excluded substring.
63 */
64
65static PetscErrorCode AssertFileNotContains(const char *path, const char *needle, const char *context)
66{
67 FILE *fp = NULL;
68 long file_size = 0;
69 char *buffer = NULL;
70
71 PetscFunctionBeginUser;
72 fp = fopen(path, "rb");
73 if (!fp) {
74 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Failed to open '%s' for assertion.", path);
75 }
76 if (fseek(fp, 0, SEEK_END) != 0) {
77 fclose(fp);
78 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Failed to seek '%s'.", path);
79 }
80 file_size = ftell(fp);
81 if (file_size < 0) {
82 fclose(fp);
83 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Failed to measure '%s'.", path);
84 }
85 if (fseek(fp, 0, SEEK_SET) != 0) {
86 fclose(fp);
87 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Failed to rewind '%s'.", path);
88 }
89
90 PetscCall(PetscMalloc1((size_t)file_size + 1, &buffer));
91 if (file_size > 0) {
92 size_t bytes_read = fread(buffer, 1, (size_t)file_size, fp);
93 if (bytes_read != (size_t)file_size) {
94 fclose(fp);
95 PetscCall(PetscFree(buffer));
96 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Failed to read '%s'.", path);
97 }
98 }
99 buffer[file_size] = '\0';
100 fclose(fp);
101
102 PetscCall(PicurvAssertBool((PetscBool)(strstr(buffer, needle) == NULL), context));
103 PetscCall(PetscFree(buffer));
104 PetscFunctionReturn(0);
105}
106
107typedef PetscErrorCode (*CapturedLoggingFn)(UserCtx *user, SimCtx *simCtx, void *ctx);
108
109typedef struct AnatomyCaptureCtx {
110 const char *field_name;
111 const char *stage_name;
113
114/**
115 * @brief Captures stdout emitted by one logging helper into a temporary file-backed buffer.
116 */
117static PetscErrorCode CaptureLoggingOutput(UserCtx *user,
118 SimCtx *simCtx,
120 void *ctx,
121 char *captured,
122 size_t captured_len)
123{
124 char tmpdir[PETSC_MAX_PATH_LEN];
125 char capture_path[PETSC_MAX_PATH_LEN];
126 FILE *capture_file = NULL;
127 int saved_stdout = -1;
128 int capture_fd = -1;
129 size_t bytes_read = 0;
130 PetscErrorCode ierr;
131
132 PetscFunctionBeginUser;
133 PetscCheck(fn != NULL, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Capture callback cannot be NULL.");
134 PetscCheck(captured != NULL, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Capture buffer cannot be NULL.");
135 PetscCheck(captured_len > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Capture buffer must be non-empty.");
136
137 PetscCall(PicurvMakeTempDir(tmpdir, sizeof(tmpdir)));
138 PetscCall(PetscSNPrintf(capture_path, sizeof(capture_path), "%s/logging.out", tmpdir));
139
140 fflush(stdout);
141 saved_stdout = dup(STDOUT_FILENO);
142 PetscCheck(saved_stdout >= 0, PETSC_COMM_SELF, PETSC_ERR_SYS, "dup(STDOUT_FILENO) failed.");
143 capture_fd = open(capture_path, O_CREAT | O_TRUNC | O_WRONLY, 0600);
144 PetscCheck(capture_fd >= 0, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Failed to open capture file '%s'.", capture_path);
145 PetscCheck(dup2(capture_fd, STDOUT_FILENO) >= 0, PETSC_COMM_SELF, PETSC_ERR_SYS, "dup2() failed while redirecting stdout.");
146 close(capture_fd);
147 capture_fd = -1;
148
149 ierr = fn(user, simCtx, ctx);
150 fflush(stdout);
151 PetscCheck(dup2(saved_stdout, STDOUT_FILENO) >= 0, PETSC_COMM_SELF, PETSC_ERR_SYS, "dup2() failed while restoring stdout.");
152 close(saved_stdout);
153 saved_stdout = -1;
154 PetscCall(ierr);
155
156 capture_file = fopen(capture_path, "r");
157 PetscCheck(capture_file != NULL, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Failed to read capture file '%s'.", capture_path);
158 bytes_read = fread(captured, 1, captured_len - 1, capture_file);
159 captured[bytes_read] = '\0';
160 fclose(capture_file);
161 PetscCall(PicurvRemoveTempDir(tmpdir));
162 PetscFunctionReturn(0);
163}
164
165/**
166 * @brief Adapts `LOG_PARTICLE_FIELDS()` to the generic stdout-capture callback shape.
167 */
168static PetscErrorCode InvokeParticleFieldLog(UserCtx *user, SimCtx *simCtx, void *ctx)
169{
170 PetscInt print_interval = *((PetscInt *)ctx);
171
172 PetscFunctionBeginUser;
173 (void)simCtx;
174 PetscCall(LOG_PARTICLE_FIELDS(user, print_interval));
175 PetscFunctionReturn(0);
176}
177
178/**
179 * @brief Adapts `EmitParticleConsoleSnapshot()` to the generic stdout-capture callback shape.
180 */
181static PetscErrorCode InvokeParticleConsoleSnapshot(UserCtx *user, SimCtx *simCtx, void *ctx)
182{
183 PetscInt step = *((PetscInt *)ctx);
184
185 PetscFunctionBeginUser;
186 PetscCall(EmitParticleConsoleSnapshot(user, simCtx, step));
187 PetscFunctionReturn(0);
188}
189
190/**
191 * @brief Adapts `LOG_FIELD_ANATOMY()` to the generic stdout-capture callback shape.
192 */
193static PetscErrorCode InvokeFieldAnatomyLog(UserCtx *user, SimCtx *simCtx, void *ctx)
194{
195 const AnatomyCaptureCtx *anatomy_ctx = (const AnatomyCaptureCtx *)ctx;
196
197 PetscFunctionBeginUser;
198 (void)simCtx;
199 PetscCall(LOG_FIELD_ANATOMY(user, anatomy_ctx->field_name, anatomy_ctx->stage_name));
200 PetscFunctionReturn(0);
201}
202
203/**
204 * @brief Creates a small particle-bearing runtime fixture used by logging tests.
205 */
206static PetscErrorCode SeedLoggingParticleFixture(SimCtx **simCtx_out, UserCtx **user_out)
207{
208 PetscReal *positions = NULL;
209 PetscReal *velocities = NULL;
210 PetscReal *weights = NULL;
211 PetscInt *cell_ids = NULL;
212 PetscInt *status = NULL;
213 PetscReal *psi = NULL;
214
215 PetscFunctionBeginUser;
216 PetscCall(PicurvCreateMinimalContexts(simCtx_out, user_out, 4, 4, 4));
217 PetscCall(PicurvCreateSwarmPair(*user_out, 2, "ske"));
218 (*simCtx_out)->np = 2;
219 (*simCtx_out)->particleConsoleOutputFreq = 2;
220 (*simCtx_out)->LoggingFrequency = 1;
221
222 PetscCall(DMSwarmGetField((*user_out)->swarm, "position", NULL, NULL, (void **)&positions));
223 PetscCall(DMSwarmGetField((*user_out)->swarm, "velocity", NULL, NULL, (void **)&velocities));
224 PetscCall(DMSwarmGetField((*user_out)->swarm, "weight", NULL, NULL, (void **)&weights));
225 PetscCall(DMSwarmGetField((*user_out)->swarm, "DMSwarm_CellID", NULL, NULL, (void **)&cell_ids));
226 PetscCall(DMSwarmGetField((*user_out)->swarm, "DMSwarm_location_status", NULL, NULL, (void **)&status));
227 PetscCall(DMSwarmGetField((*user_out)->swarm, "Psi", NULL, NULL, (void **)&psi));
228
229 positions[0] = 0.25; positions[1] = 0.50; positions[2] = 0.75;
230 positions[3] = 0.50; positions[4] = 0.50; positions[5] = 0.50;
231 velocities[0] = 1.0; velocities[1] = 2.0; velocities[2] = 3.0;
232 velocities[3] = 4.0; velocities[4] = 5.0; velocities[5] = 6.0;
233 weights[0] = 0.2; weights[1] = 0.3; weights[2] = 0.4;
234 weights[3] = 0.5; weights[4] = 0.5; weights[5] = 0.5;
235 cell_ids[0] = 0; cell_ids[1] = 0; cell_ids[2] = 0;
236 cell_ids[3] = 1; cell_ids[4] = 1; cell_ids[5] = 1;
237 status[0] = ACTIVE_AND_LOCATED;
238 status[1] = ACTIVE_AND_LOCATED;
239 psi[0] = 1.0;
240 psi[1] = 3.0;
241
242 PetscCall(DMSwarmRestoreField((*user_out)->swarm, "Psi", NULL, NULL, (void **)&psi));
243 PetscCall(DMSwarmRestoreField((*user_out)->swarm, "DMSwarm_location_status", NULL, NULL, (void **)&status));
244 PetscCall(DMSwarmRestoreField((*user_out)->swarm, "DMSwarm_CellID", NULL, NULL, (void **)&cell_ids));
245 PetscCall(DMSwarmRestoreField((*user_out)->swarm, "weight", NULL, NULL, (void **)&weights));
246 PetscCall(DMSwarmRestoreField((*user_out)->swarm, "velocity", NULL, NULL, (void **)&velocities));
247 PetscCall(DMSwarmRestoreField((*user_out)->swarm, "position", NULL, NULL, (void **)&positions));
248 PetscFunctionReturn(0);
249}
250/**
251 * @brief Tests string-conversion helpers for configured enums and unknown values.
252 */
253
254static PetscErrorCode TestStringConversionHelpers(void)
255{
256 PetscFunctionBeginUser;
257 PetscCall(PicurvAssertBool(strcmp(BCFaceToString(BC_FACE_NEG_X), "-Xi (I-Min)") == 0,
258 "BCFaceToString should report the negative-x face"));
259 PetscCall(PicurvAssertBool(strcmp(FieldInitializationToString(0), "Zero") == 0,
260 "FieldInitializationToString should report the zero mode"));
261 PetscCall(PicurvAssertBool(strcmp(FieldInitializationToString(99), "Unknown Field Initialization") == 0,
262 "FieldInitializationToString should reject unknown selectors"));
264 "ParticleInitializationToString should report the volume mode"));
265 PetscCall(PicurvAssertBool(strcmp(LESModelToString(CONSTANT_SMAGORINSKY), "Constant Smagorinsky") == 0,
266 "LESModelToString should report the constant model"));
267 PetscCall(PicurvAssertBool(strcmp(MomentumSolverTypeToString(MOMENTUM_SOLVER_EXPLICIT_RK), "Explicit 4 stage Runge-Kutta ") == 0,
268 "MomentumSolverTypeToString should report the explicit solver"));
269 PetscCall(PicurvAssertBool(strcmp(BCTypeToString(PERIODIC), "PERIODIC") == 0,
270 "BCTypeToString should report periodic boundaries"));
272 "BCHandlerTypeToString should report the driven periodic handler"));
273 PetscCall(PicurvAssertBool(strcmp(ParticleLocationStatusToString(LOST), "LOST") == 0,
274 "ParticleLocationStatusToString should report LOST state"));
275 PetscFunctionReturn(0);
276}
277/**
278 * @brief Tests that log level selection honors the environment variable.
279 */
280
281static PetscErrorCode TestGetLogLevelFromEnvironment(void)
282{
283 PetscFunctionBeginUser;
285 "get_log_level should honor LOG_LEVEL=INFO in this test binary"));
286 PetscCall(print_log_level());
287 PetscFunctionReturn(0);
288}
289/**
290 * @brief Tests the function allow-list filter used by the logging layer.
291 */
292
293static PetscErrorCode TestAllowedFunctionsFilter(void)
294{
295 const char *allow_list[] = {"ComputeSpecificKE", "WriteEulerianFile"};
296
297 PetscFunctionBeginUser;
298 set_allowed_functions(allow_list, 2);
299 PetscCall(PicurvAssertBool(is_function_allowed("ComputeSpecificKE"),
300 "Allowed list should include ComputeSpecificKE"));
301 PetscCall(PicurvAssertBool((PetscBool)!is_function_allowed("UnlistedFunction"),
302 "Allowed list should exclude unknown function names"));
303
304 set_allowed_functions(NULL, 0);
305 PetscCall(PicurvAssertBool(is_function_allowed("AnyFunction"),
306 "Empty allow-list should permit all functions"));
307 PetscFunctionReturn(0);
308}
309/**
310 * @brief Tests periodic particle console snapshot enablement and cadence.
311 */
312
313static PetscErrorCode TestParticleConsoleSnapshotCadence(void)
314{
315 SimCtx simCtx;
316
317 PetscFunctionBeginUser;
318 PetscCall(PetscMemzero(&simCtx, sizeof(simCtx)));
319 simCtx.np = 32;
320 simCtx.particleConsoleOutputFreq = 4;
321
323 "Particle snapshot contract should be enabled when particles and cadence are configured"));
325 "Snapshot should emit on cadence-aligned completed steps"));
326 PetscCall(PicurvAssertBool((PetscBool)!ShouldEmitPeriodicParticleConsoleSnapshot(&simCtx, 7),
327 "Snapshot should not emit off-cadence"));
328
329 simCtx.particleConsoleOutputFreq = 0;
330 PetscCall(PicurvAssertBool((PetscBool)!IsParticleConsoleSnapshotEnabled(&simCtx),
331 "Zero cadence should disable periodic particle snapshots"));
332 PetscCall(PicurvAssertBool((PetscBool)!ShouldEmitPeriodicParticleConsoleSnapshot(NULL, 4),
333 "NULL SimCtx should never emit periodic snapshots"));
334 PetscFunctionReturn(0);
335}
336/**
337 * @brief Tests logging-side file parsing, helper formatting, and progress utilities.
338 */
339
341{
342 char tmpdir[PETSC_MAX_PATH_LEN];
343 char allow_path[PETSC_MAX_PATH_LEN];
344 char dual_log_path[PETSC_MAX_PATH_LEN];
345 FILE *file = NULL;
346 char **funcs = NULL;
347 PetscInt nfuncs = 0;
348 Cell cell;
349 PetscReal distances[6] = {1.0, 2.0, 3.0, 4.0, 5.0, 6.0};
350 DualMonitorCtx *monctx = NULL;
351 void *ctx = NULL;
352
353 PetscFunctionBeginUser;
354 PetscCall(PetscMemzero(&cell, sizeof(cell)));
355 PetscCall(PicurvMakeTempDir(tmpdir, sizeof(tmpdir)));
356 PetscCall(PetscSNPrintf(allow_path, sizeof(allow_path), "%s/allowed_functions.txt", tmpdir));
357 PetscCall(PetscSNPrintf(dual_log_path, sizeof(dual_log_path), "%s/dual-monitor.log", tmpdir));
358
359 file = fopen(allow_path, "w");
360 PetscCheck(file != NULL, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Failed to create allow-list file '%s'.", allow_path);
361 fputs(" ComputeSpecificKE \n", file);
362 fputs("# comment-only line\n", file);
363 for (PetscInt i = 0; i < 17; ++i) {
364 fprintf(file, "Helper_%02d # trailing comment\n", (int)i);
365 }
366 fclose(file);
367 file = NULL;
368
369 PetscCall(LoadAllowedFunctionsFromFile(allow_path, &funcs, &nfuncs));
370 PetscCall(PicurvAssertIntEqual(18, nfuncs, "LoadAllowedFunctionsFromFile should trim comments and keep all identifiers"));
371 PetscCall(PicurvAssertBool((PetscBool)(strcmp(funcs[0], "ComputeSpecificKE") == 0),
372 "LoadAllowedFunctionsFromFile should trim leading and trailing whitespace"));
373 PetscCall(PicurvAssertBool((PetscBool)(strcmp(funcs[17], "Helper_16") == 0),
374 "LoadAllowedFunctionsFromFile should grow past the initial pointer capacity"));
375 PetscCall(FreeAllowedFunctions(funcs, nfuncs));
376
377 for (PetscInt i = 0; i < 8; ++i) {
378 cell.vertices[i].x = (PetscReal)i;
379 cell.vertices[i].y = (PetscReal)(i + 1);
380 cell.vertices[i].z = (PetscReal)(i + 2);
381 }
382 PetscCall(LOG_CELL_VERTICES(&cell, 0));
383 PetscCall(LOG_FACE_DISTANCES(distances));
384
385 PetscCall(PetscCalloc1(1, &monctx));
386 monctx->file_handle = fopen(dual_log_path, "w");
387 PetscCheck(monctx->file_handle != NULL, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Failed to create dual-monitor log '%s'.", dual_log_path);
388 ctx = monctx;
389 PetscCall(DualMonitorDestroy(&ctx));
390 PetscCall(PicurvAssertBool((PetscBool)(ctx == NULL),
391 "DualMonitorDestroy should clear the caller-owned context pointer"));
392
393 PrintProgressBar(0, 0, 4, 0.10);
394 PrintProgressBar(3, 0, 4, 0.40);
395 PrintProgressBar(0, 0, 0, 0.00);
396 PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
397
398 PetscCall(PicurvRemoveTempDir(tmpdir));
399 PetscFunctionReturn(0);
400}
401/**
402 * @brief Tests continuity, min/max, and anatomy logging helpers on minimal runtime fixtures.
403 */
404
406{
407 SimCtx *simCtx = NULL;
408 UserCtx *user = NULL;
409 char tmpdir[PETSC_MAX_PATH_LEN];
410 char continuity_path[PETSC_MAX_PATH_LEN];
411 PetscReal ***p = NULL;
412 Cmpnts ***ucat = NULL;
413 Cmpnts ***ucont = NULL;
414 PetscErrorCode ierr_minmax = 0;
415
416 PetscFunctionBeginUser;
417 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
418 PetscCall(PicurvMakeTempDir(tmpdir, sizeof(tmpdir)));
419 PetscCall(PetscStrncpy(simCtx->log_dir, tmpdir, sizeof(simCtx->log_dir)));
420
421 simCtx->StartStep = 0;
422 simCtx->step = 1;
423 simCtx->MaxDiv = 1.25;
424 simCtx->MaxDivx = 1;
425 simCtx->MaxDivy = 2;
426 simCtx->MaxDivz = 3;
427 simCtx->MaxDivFlatArg = 17;
428 simCtx->summationRHS = 8.5;
429 simCtx->FluxInSum = 5.0;
430 simCtx->FluxOutSum = 3.25;
431 PetscCall(LOG_CONTINUITY_METRICS(user));
432
433 simCtx->step = 2;
434 simCtx->MaxDiv = 0.75;
435 simCtx->summationRHS = 4.5;
436 simCtx->FluxInSum = 2.5;
437 simCtx->FluxOutSum = 1.0;
438 PetscCall(LOG_CONTINUITY_METRICS(user));
439
440 PetscCall(PetscSNPrintf(continuity_path, sizeof(continuity_path), "%s/Continuity_Metrics.log", simCtx->log_dir));
441 PetscCall(PicurvAssertFileExists(continuity_path, "continuity metrics log should be written"));
442 PetscCall(AssertFileContains(continuity_path, "Timestep", "continuity metrics log should include the header"));
443 PetscCall(AssertFileContains(continuity_path, "([3][2][1] = 17)", "continuity metrics log should include the divergence location"));
444 PetscCall(AssertFileContains(continuity_path, "2 | 0", "continuity metrics log should append later timesteps"));
445
446 PetscCall(DMDAVecGetArray(user->da, user->P, &p));
447 PetscCall(DMDAVecGetArray(user->fda, user->Ucat, &ucat));
448 PetscCall(DMDAVecGetArray(user->fda, user->Ucont, &ucont));
449 for (PetscInt k = user->info.zs; k < user->info.zs + user->info.zm; ++k) {
450 for (PetscInt j = user->info.ys; j < user->info.ys + user->info.ym; ++j) {
451 for (PetscInt i = user->info.xs; i < user->info.xs + user->info.xm; ++i) {
452 p[k][j][i] = (PetscReal)(i + j + k);
453 ucat[k][j][i].x = (PetscReal)i;
454 ucat[k][j][i].y = (PetscReal)(-j);
455 ucat[k][j][i].z = (PetscReal)(2 * k);
456 ucont[k][j][i].x = (PetscReal)(10 + i);
457 ucont[k][j][i].y = (PetscReal)(20 + j);
458 ucont[k][j][i].z = (PetscReal)(30 + k);
459 }
460 }
461 }
462 PetscCall(DMDAVecRestoreArray(user->fda, user->Ucont, &ucont));
463 PetscCall(DMDAVecRestoreArray(user->fda, user->Ucat, &ucat));
464 PetscCall(DMDAVecRestoreArray(user->da, user->P, &p));
465
466 PetscCall(DMGlobalToLocalBegin(user->da, user->P, INSERT_VALUES, user->lP));
467 PetscCall(DMGlobalToLocalEnd(user->da, user->P, INSERT_VALUES, user->lP));
468 PetscCall(DMGlobalToLocalBegin(user->fda, user->Ucat, INSERT_VALUES, user->lUcat));
469 PetscCall(DMGlobalToLocalEnd(user->fda, user->Ucat, INSERT_VALUES, user->lUcat));
470 PetscCall(DMGlobalToLocalBegin(user->fda, user->Ucont, INSERT_VALUES, user->lUcont));
471 PetscCall(DMGlobalToLocalEnd(user->fda, user->Ucont, INSERT_VALUES, user->lUcont));
472
473 PetscCall(LOG_FIELD_MIN_MAX(user, "P"));
474 PetscCall(LOG_FIELD_MIN_MAX(user, "Ucat"));
475 PetscCall(LOG_FIELD_MIN_MAX(user, "Coordinates"));
476 PetscCall(LOG_FIELD_MIN_MAX(user, "Ucont"));
477
478 PetscCall(PetscPushErrorHandler(PetscIgnoreErrorHandler, NULL));
479 ierr_minmax = LOG_FIELD_MIN_MAX(user, "NotARealField");
480 PetscCall(PetscPopErrorHandler());
481 PetscCall(PicurvAssertBool((PetscBool)(ierr_minmax != 0),
482 "LOG_FIELD_MIN_MAX should reject unknown field names"));
483
484 PetscCall(PicurvRemoveTempDir(tmpdir));
485 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
486 PetscFunctionReturn(0);
487}
488/**
489 * @brief Tests interpolation-error logging against an analytically matched particle field.
490 */
491
492static PetscErrorCode TestInterpolationErrorLogging(void)
493{
494 SimCtx *simCtx = NULL;
495 UserCtx *user = NULL;
496 PetscReal (*pos_arr)[3] = NULL;
497 PetscReal (*vel_arr)[3] = NULL;
498 Vec position_vec = NULL;
499 Vec analytical_vec = NULL;
500 const PetscScalar *analytical_arr = NULL;
501
502 PetscFunctionBeginUser;
503 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
504 PetscCall(PicurvCreateSwarmPair(user, 2, "ske"));
505 PetscCall(PetscStrncpy(simCtx->AnalyticalSolutionType, "TGV3D", sizeof(simCtx->AnalyticalSolutionType)));
506 simCtx->ren = 1.0;
507 simCtx->ti = 0.0;
508
509 PetscCall(DMSwarmGetField(user->swarm, "position", NULL, NULL, (void *)&pos_arr));
510 pos_arr[0][0] = 0.5 * PETSC_PI; pos_arr[0][1] = 0.0; pos_arr[0][2] = 0.0;
511 pos_arr[1][0] = 0.0; pos_arr[1][1] = 0.5 * PETSC_PI; pos_arr[1][2] = 0.0;
512 PetscCall(DMSwarmRestoreField(user->swarm, "position", NULL, NULL, (void *)&pos_arr));
513
514 PetscCall(DMSwarmCreateGlobalVectorFromField(user->swarm, "position", &position_vec));
515 PetscCall(VecDuplicate(position_vec, &analytical_vec));
516 PetscCall(VecCopy(position_vec, analytical_vec));
517 PetscCall(SetAnalyticalSolutionForParticles(analytical_vec, simCtx));
518
519 PetscCall(DMSwarmGetField(user->swarm, "velocity", NULL, NULL, (void *)&vel_arr));
520 PetscCall(VecGetArrayRead(analytical_vec, &analytical_arr));
521 for (PetscInt particle = 0; particle < 2; ++particle) {
522 vel_arr[particle][0] = PetscRealPart(analytical_arr[3 * particle + 0]);
523 vel_arr[particle][1] = PetscRealPart(analytical_arr[3 * particle + 1]);
524 vel_arr[particle][2] = PetscRealPart(analytical_arr[3 * particle + 2]);
525 }
526 PetscCall(VecRestoreArrayRead(analytical_vec, &analytical_arr));
527 PetscCall(DMSwarmRestoreField(user->swarm, "velocity", NULL, NULL, (void *)&vel_arr));
528
529 PetscCall(VecDestroy(&analytical_vec));
530 PetscCall(DMSwarmDestroyGlobalVectorFromField(user->swarm, "position", &position_vec));
531 PetscCall(LOG_INTERPOLATION_ERROR(user));
532 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
533 PetscFunctionReturn(0);
534}
535/**
536 * @brief Tests file-backed scatter metrics logging against a fully occupied constant field.
537 */
538
539static PetscErrorCode TestScatterMetricsLogging(void)
540{
541 SimCtx *simCtx = NULL;
542 UserCtx *user = NULL;
543 char tmpdir[PETSC_MAX_PATH_LEN];
544 char metrics_path[PETSC_MAX_PATH_LEN];
545 PetscReal *positions = NULL;
546 PetscReal *psi = NULL;
547 PetscInt *cell_ids = NULL;
548 PetscInt *status = NULL;
549 PetscInt particle = 0;
550
551 PetscFunctionBeginUser;
552 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
553 PetscCall(PicurvCreateSwarmPair(user, 27, "ske"));
554 PetscCall(PicurvMakeTempDir(tmpdir, sizeof(tmpdir)));
555 PetscCall(PetscStrncpy(simCtx->log_dir, tmpdir, sizeof(simCtx->log_dir)));
556 simCtx->np = 27;
557 simCtx->step = 3;
558 simCtx->ti = 0.3;
559 simCtx->verificationScalar.enabled = PETSC_TRUE;
560 PetscCall(PetscStrncpy(simCtx->verificationScalar.mode,
561 "analytical",
562 sizeof(simCtx->verificationScalar.mode)));
563 PetscCall(PetscStrncpy(simCtx->verificationScalar.profile,
564 "CONSTANT",
565 sizeof(simCtx->verificationScalar.profile)));
566 simCtx->verificationScalar.value = 2.0;
567
568 PetscCall(DMSwarmGetField(user->swarm, "position", NULL, NULL, (void **)&positions));
569 PetscCall(DMSwarmGetField(user->swarm, "DMSwarm_CellID", NULL, NULL, (void **)&cell_ids));
570 PetscCall(DMSwarmGetField(user->swarm, "DMSwarm_location_status", NULL, NULL, (void **)&status));
571 PetscCall(DMSwarmGetField(user->swarm, "Psi", NULL, NULL, (void **)&psi));
572
573 for (PetscInt k = 0; k < 3; ++k) {
574 for (PetscInt j = 0; j < 3; ++j) {
575 for (PetscInt i = 0; i < 3; ++i) {
576 positions[3 * particle + 0] = (i + 0.5) / 4.0;
577 positions[3 * particle + 1] = (j + 0.5) / 4.0;
578 positions[3 * particle + 2] = (k + 0.5) / 4.0;
579 cell_ids[3 * particle + 0] = i;
580 cell_ids[3 * particle + 1] = j;
581 cell_ids[3 * particle + 2] = k;
582 status[particle] = ACTIVE_AND_LOCATED;
583 psi[particle] = 2.0;
584 ++particle;
585 }
586 }
587 }
588
589 PetscCall(DMSwarmRestoreField(user->swarm, "Psi", NULL, NULL, (void **)&psi));
590 PetscCall(DMSwarmRestoreField(user->swarm, "DMSwarm_location_status", NULL, NULL, (void **)&status));
591 PetscCall(DMSwarmRestoreField(user->swarm, "DMSwarm_CellID", NULL, NULL, (void **)&cell_ids));
592 PetscCall(DMSwarmRestoreField(user->swarm, "position", NULL, NULL, (void **)&positions));
593
595 PetscCall(LOG_SCATTER_METRICS(user));
596
597 PetscCall(PetscSNPrintf(metrics_path, sizeof(metrics_path), "%s/scatter_metrics.csv", simCtx->log_dir));
598 PetscCall(PicurvAssertFileExists(metrics_path, "LOG_SCATTER_METRICS should write scatter_metrics.csv"));
599 PetscCall(AssertFileContains(metrics_path, "relative_L2_error",
600 "Scatter metrics CSV header should include relative_L2_error"));
601 PetscCall(AssertFileContains(metrics_path,
602 "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",
603 "Scatter metrics CSV should record the expected constant-field zero-error row"));
604
605 PetscCall(PicurvRemoveTempDir(tmpdir));
606 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
607 PetscFunctionReturn(0);
608}
609/**
610 * @brief Tests stdout particle-table logging on a production-like swarm fixture.
611 */
612
613static PetscErrorCode TestParticleFieldTableLogging(void)
614{
615 SimCtx *simCtx = NULL;
616 UserCtx *user = NULL;
617 PetscInt print_interval = 1;
618 char captured[8192];
619
620 PetscFunctionBeginUser;
621 PetscCall(SeedLoggingParticleFixture(&simCtx, &user));
622 PetscCall(CaptureLoggingOutput(user, simCtx, InvokeParticleFieldLog, &print_interval, captured, sizeof(captured)));
623 PetscCall(PicurvAssertBool((PetscBool)(strstr(captured, "Position (x,y,z)") != NULL),
624 "LOG_PARTICLE_FIELDS should print the particle table header"));
625 PetscCall(PicurvAssertBool((PetscBool)(strstr(captured, "Weights (a1,a2,a3)") != NULL),
626 "LOG_PARTICLE_FIELDS should print the weight-column header"));
627 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
628 PetscFunctionReturn(0);
629}
630/**
631 * @brief Tests console snapshot logging against the public periodic-snapshot helper.
632 */
633
634static PetscErrorCode TestParticleConsoleSnapshotLogging(void)
635{
636 SimCtx *simCtx = NULL;
637 UserCtx *user = NULL;
638 PetscInt step = 4;
639 char captured[8192];
640
641 PetscFunctionBeginUser;
642 PetscCall(SeedLoggingParticleFixture(&simCtx, &user));
643 PetscCall(CaptureLoggingOutput(user, simCtx, InvokeParticleConsoleSnapshot, &step, captured, sizeof(captured)));
644 PetscCall(PicurvAssertBool((PetscBool)(strstr(captured, "Particle states at step 4") != NULL),
645 "EmitParticleConsoleSnapshot should print the step banner"));
646 PetscCall(PicurvAssertBool((PetscBool)(strstr(captured, "Position (x,y,z)") != NULL),
647 "EmitParticleConsoleSnapshot should reuse the particle table output"));
648 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
649 PetscFunctionReturn(0);
650}
651/**
652 * @brief Tests file-backed particle metrics logging after derived metrics are computed.
653 */
654
655static PetscErrorCode TestParticleMetricsLogging(void)
656{
657 SimCtx *simCtx = NULL;
658 UserCtx *user = NULL;
659 char tmpdir[PETSC_MAX_PATH_LEN];
660 char metrics_path[PETSC_MAX_PATH_LEN];
661
662 PetscFunctionBeginUser;
663 PetscCall(SeedLoggingParticleFixture(&simCtx, &user));
664 PetscCall(PicurvMakeTempDir(tmpdir, sizeof(tmpdir)));
665 PetscCall(PetscStrncpy(simCtx->log_dir, tmpdir, sizeof(simCtx->log_dir)));
666 simCtx->StartStep = 0;
667 simCtx->step = 1;
668 simCtx->particlesLostLastStep = 1;
669 simCtx->particlesLostCumulative = 7;
670 simCtx->particlesMigratedLastStep = 2;
671 simCtx->migrationPassesLastStep = 3;
672
673 PetscCall(CalculateParticleCountPerCell(user));
674 PetscCall(CalculateAdvancedParticleMetrics(user));
675 PetscCall(LOG_PARTICLE_METRICS(user, "Timestep Metrics"));
676
677 PetscCall(PetscSNPrintf(metrics_path, sizeof(metrics_path), "%s/Particle_Metrics.log", simCtx->log_dir));
678 PetscCall(PicurvAssertFileExists(metrics_path, "LOG_PARTICLE_METRICS should write Particle_Metrics.log"));
679 PetscCall(AssertFileContains(metrics_path, "Timestep Metrics", "Particle metrics log should include the caller-provided stage label"));
680 PetscCall(AssertFileContains(metrics_path, "Occupied Cells", "Particle metrics log should include the metrics table header"));
681 PetscCall(AssertFileContains(metrics_path, "Lost Total", "Particle metrics log should include the cumulative-loss column"));
682 PetscCall(AssertFileContains(metrics_path, "| 1 | 7 | 2", "Particle metrics log should record both per-step and cumulative loss values"));
683 PetscCall(PicurvRemoveTempDir(tmpdir));
684 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
685 PetscFunctionReturn(0);
686}
687/**
688 * @brief Tests file-backed search metrics logging with the compact CSV contract.
689 */
690
691static PetscErrorCode TestSearchMetricsLogging(void)
692{
693 SimCtx *simCtx = NULL;
694 UserCtx *user = NULL;
695 char tmpdir[PETSC_MAX_PATH_LEN];
696 char metrics_path[PETSC_MAX_PATH_LEN];
697
698 PetscFunctionBeginUser;
699 PetscCall(SeedLoggingParticleFixture(&simCtx, &user));
700 PetscCall(PicurvMakeTempDir(tmpdir, sizeof(tmpdir)));
701 PetscCall(PetscStrncpy(simCtx->log_dir, tmpdir, sizeof(simCtx->log_dir)));
702 simCtx->step = 2;
703 simCtx->ti = 0.2;
704 simCtx->particlesLostLastStep = 1;
705 simCtx->particlesLostCumulative = 7;
706 simCtx->particlesMigratedLastStep = 2;
707 simCtx->migrationPassesLastStep = 3;
708 simCtx->particleLoadImbalance = 1.5;
709 simCtx->searchMetrics.searchAttempts = 4;
712 simCtx->searchMetrics.searchLostCount = 1;
713 simCtx->searchMetrics.traversalStepsSum = 10;
714 simCtx->searchMetrics.reSearchCount = 2;
717 simCtx->searchMetrics.tieBreakCount = 1;
722
723 PetscCall(LOG_SEARCH_METRICS(user));
724
725 PetscCall(PetscSNPrintf(metrics_path, sizeof(metrics_path), "%s/search_metrics.csv", simCtx->log_dir));
726 PetscCall(PicurvAssertFileExists(metrics_path, "LOG_SEARCH_METRICS should write search_metrics.csv"));
727 PetscCall(AssertFileContains(metrics_path, "search_attempts", "Search metrics CSV header should include search_attempts"));
728 PetscCall(AssertFileContains(metrics_path, "search_population", "Search metrics CSV header should include search_population"));
729 PetscCall(AssertFileContains(metrics_path, "max_particle_pass_depth", "Search metrics CSV header should include max_particle_pass_depth"));
730 PetscCall(AssertFileContains(metrics_path, "search_work_index", "Search metrics CSV header should include search_work_index"));
731 PetscCall(AssertFileContains(metrics_path, "re_search_fraction", "Search metrics CSV header should include re_search_fraction"));
732 PetscCall(AssertFileContains(metrics_path, "lost_cumulative", "Search metrics CSV header should include the cumulative-loss column"));
733 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"));
734 PetscCall(PicurvRemoveTempDir(tmpdir));
735 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
736 PetscFunctionReturn(0);
737}
738/**
739 * @brief Tests stdout field-anatomy logging on the corrected production-like DM fixture.
740 */
741
742static PetscErrorCode TestFieldAnatomyLogging(void)
743{
744 SimCtx *simCtx = NULL;
745 UserCtx *user = NULL;
746 char captured[8192];
747 AnatomyCaptureCtx anatomy_ctx = {"P", "unit-test"};
748
749 PetscFunctionBeginUser;
750 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
751 PetscCall(VecSet(user->P, 7.0));
752 PetscCall(DMGlobalToLocalBegin(user->da, user->P, INSERT_VALUES, user->lP));
753 PetscCall(DMGlobalToLocalEnd(user->da, user->P, INSERT_VALUES, user->lP));
754 PetscCall(CaptureLoggingOutput(user, simCtx, InvokeFieldAnatomyLog, &anatomy_ctx, captured, sizeof(captured)));
755 PetscCall(PicurvAssertBool((PetscBool)(strstr(captured, "Field Anatomy Log: [P]") != NULL),
756 "LOG_FIELD_ANATOMY should print the requested field name"));
757 PetscCall(PicurvAssertBool((PetscBool)(strstr(captured, "Layout: [Cell-Centered]") != NULL),
758 "LOG_FIELD_ANATOMY should report the inferred data layout"));
759 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
760 PetscFunctionReturn(0);
761}
762/**
763 * @brief Tests profiling helper lifecycle logging for timestep and final-summary outputs.
764 */
765
766static PetscErrorCode TestProfilingLifecycleHelpers(void)
767{
768 SimCtx simCtx;
769 char tmpdir[PETSC_MAX_PATH_LEN];
770 char timestep_path[PETSC_MAX_PATH_LEN];
771 char summary_path[PETSC_MAX_PATH_LEN];
772 static char selected_name[] = "FlowSolver";
773 char *selected_funcs[] = {selected_name};
774
775 PetscFunctionBeginUser;
776 PetscCall(PetscMemzero(&simCtx, sizeof(simCtx)));
777 PetscCall(PicurvMakeTempDir(tmpdir, sizeof(tmpdir)));
778 PetscCall(PetscStrncpy(simCtx.log_dir, tmpdir, sizeof(simCtx.log_dir)));
779 PetscCall(PetscStrncpy(simCtx.profilingTimestepMode, "selected", sizeof(simCtx.profilingTimestepMode)));
780 PetscCall(PetscStrncpy(simCtx.profilingTimestepFile, "Profiling_Timestep_Summary.csv", sizeof(simCtx.profilingTimestepFile)));
781 simCtx.rank = 0;
783 simCtx.StartStep = 0;
784 simCtx.nProfilingSelectedFuncs = 1;
785 simCtx.profilingSelectedFuncs = selected_funcs;
786 simCtx.profilingFinalSummary = PETSC_TRUE;
787
788 PetscCall(ProfilingInitialize(&simCtx));
789
790 _ProfilingStart("FlowSolver");
791 _ProfilingEnd("FlowSolver");
792 _ProfilingStart("UnselectedHelper");
793 _ProfilingEnd("UnselectedHelper");
794 PetscCall(ProfilingLogTimestepSummary(&simCtx, 1));
795
796 _ProfilingStart("FlowSolver");
797 _ProfilingEnd("FlowSolver");
799 PetscCall(ProfilingLogTimestepSummary(&simCtx, 2));
800 PetscCall(ProfilingFinalize(&simCtx));
801
802 PetscCall(PetscSNPrintf(timestep_path, sizeof(timestep_path), "%s/%s", simCtx.log_dir, simCtx.profilingTimestepFile));
803 PetscCall(PetscSNPrintf(summary_path, sizeof(summary_path), "%s/ProfilingSummary_Solver.log", simCtx.log_dir));
804 PetscCall(PicurvAssertFileExists(timestep_path, "profiling timestep summary should be written"));
805 PetscCall(PicurvAssertFileExists(summary_path, "profiling final summary should be written"));
806 PetscCall(AssertFileContains(timestep_path, "step,function,calls,step_time_s",
807 "profiling timestep summary should contain the CSV header"));
808 PetscCall(AssertFileContains(timestep_path, "1,FlowSolver,1,",
809 "profiling timestep summary should log selected functions"));
810 PetscCall(AssertFileNotContains(timestep_path, "UnselectedHelper",
811 "profiling timestep summary should omit unselected functions in selected mode"));
812 PetscCall(AssertFileContains(summary_path, "FINAL PROFILING SUMMARY",
813 "profiling final summary should include its table banner"));
814 PetscCall(AssertFileContains(summary_path, "FlowSolver",
815 "profiling final summary should include selected functions"));
816 PetscCall(AssertFileContains(summary_path, "UnselectedHelper",
817 "profiling final summary should include total-time entries for unselected functions"));
818 PetscCall(PicurvRemoveTempDir(tmpdir));
819 PetscFunctionReturn(0);
820}
821/**
822 * @brief Runs the unit-logging PETSc test binary.
823 */
824
825int main(int argc, char **argv)
826{
827 PetscErrorCode ierr;
828 const PicurvTestCase cases[] = {
829 {"string-conversion-helpers", TestStringConversionHelpers},
830 {"get-log-level-from-environment", TestGetLogLevelFromEnvironment},
831 {"allowed-functions-filter", TestAllowedFunctionsFilter},
832 {"particle-console-snapshot-cadence", TestParticleConsoleSnapshotCadence},
833 {"logging-file-parsing-and-formatting-helpers", TestLoggingFileParsingAndFormattingHelpers},
834 {"logging-continuity-and-field-diagnostics", TestLoggingContinuityAndFieldDiagnostics},
835 {"interpolation-error-logging", TestInterpolationErrorLogging},
836 {"scatter-metrics-logging", TestScatterMetricsLogging},
837 {"particle-field-table-logging", TestParticleFieldTableLogging},
838 {"particle-console-snapshot-logging", TestParticleConsoleSnapshotLogging},
839 {"particle-metrics-logging", TestParticleMetricsLogging},
840 {"search-metrics-logging", TestSearchMetricsLogging},
841 {"field-anatomy-logging", TestFieldAnatomyLogging},
842 {"profiling-lifecycle-helpers", TestProfilingLifecycleHelpers},
843 };
844
845 (void)setenv("LOG_LEVEL", "INFO", 1);
846
847 ierr = PetscInitialize(&argc, &argv, NULL, "PICurv logging tests");
848 if (ierr) {
849 return (int)ierr;
850 }
851
852 ierr = PicurvRunTests("unit-logging", cases, sizeof(cases) / sizeof(cases[0]));
853 if (ierr) {
854 PetscFinalize();
855 return (int)ierr;
856 }
857
858 set_allowed_functions(NULL, 0);
859
860 ierr = PetscFinalize();
861 return (int)ierr;
862}
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:2248
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:808
PetscErrorCode LOG_INTERPOLATION_ERROR(UserCtx *user)
Logs the interpolation error between the analytical and computed solutions.
Definition logging.c:1785
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:1214
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:1367
void PrintProgressBar(PetscInt step, PetscInt startStep, PetscInt totalSteps, PetscReal currentTime)
Prints a progress bar to the console.
Definition logging.c:1320
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:1515
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:1121
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:1082
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:2194
const char * ParticleLocationStatusToString(ParticleLocationStatus level)
A function that outputs the name of the current level in the ParticleLocation enum.
Definition logging.c:973
PetscErrorCode LOG_SCATTER_METRICS(UserCtx *user)
Logs particle-to-grid scatter verification metrics for the prescribed scalar truth path.
Definition logging.c:1861
PetscErrorCode LOG_CONTINUITY_METRICS(UserCtx *user)
Logs continuity metrics for a single block to a file.
Definition logging.c:914
PetscErrorCode LOG_SEARCH_METRICS(UserCtx *user)
Writes compact runtime search metrics to CSV and optionally to console.
Definition logging.c:2043
PetscErrorCode ProfilingInitialize(SimCtx *simCtx)
Initializes the custom profiling system using configuration from SimCtx.
Definition logging.c:1044
@ 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:1104
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:1068
Context for a dual-purpose KSP monitor.
Definition logging.h:55
static PetscErrorCode TestParticleConsoleSnapshotCadence(void)
Tests periodic particle console snapshot enablement and cadence.
static PetscErrorCode TestFieldAnatomyLogging(void)
Tests stdout field-anatomy logging on the corrected production-like DM fixture.
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 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 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 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 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 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 AssertFileContains(const char *path, const char *needle, const char *context)
Asserts that one text file contains a required substring.
static PetscErrorCode TestProfilingLifecycleHelpers(void)
Tests profiling helper lifecycle logging for timestep and final-summary outputs.
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.
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 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 profilingFinalSummary
Definition variables.h:780
PetscMPIInt rank
Definition variables.h:646
char profilingTimestepFile[PETSC_MAX_PATH_LEN]
Definition variables.h:779
PetscInt64 searchLocatedCount
Definition variables.h:224
PetscInt64 searchLostCount
Definition variables.h:225
@ PARTICLE_INIT_VOLUME
Random volumetric distribution across the domain.
Definition variables.h:510
@ LOST
Definition variables.h:139
@ ACTIVE_AND_LOCATED
Definition variables.h:137
PetscReal FluxOutSum
Definition variables.h:721
PetscInt64 boundaryClampCount
Definition variables.h:231
PetscInt particlesLostLastStep
Definition variables.h:746
PetscInt64 traversalStepsSum
Definition variables.h:226
@ BC_HANDLER_PERIODIC_DRIVEN_CONSTANT_FLUX
Definition variables.h:286
PetscReal ren
Definition variables.h:691
PetscInt64 searchPopulation
Definition variables.h:223
char profilingTimestepMode[32]
Definition variables.h:778
PetscInt np
Definition variables.h:739
Vec Ucont
Definition variables.h:837
PetscInt StartStep
Definition variables.h:653
@ MOMENTUM_SOLVER_EXPLICIT_RK
Definition variables.h:503
PetscScalar x
Definition variables.h:101
PetscInt64 reSearchCount
Definition variables.h:227
PetscReal MaxDiv
Definition variables.h:771
PetscInt64 bboxGuessFallbackCount
Definition variables.h:233
VerificationScalarConfig verificationScalar
Definition variables.h:700
PetscInt MaxDivx
Definition variables.h:772
PetscInt MaxDivy
Definition variables.h:772
PetscInt64 bboxGuessSuccessCount
Definition variables.h:232
PetscInt MaxDivz
Definition variables.h:772
char log_dir[PETSC_MAX_PATH_LEN]
Definition variables.h:668
PetscInt MaxDivFlatArg
Definition variables.h:772
PetscReal FluxInSum
Definition variables.h:721
PetscInt64 maxParticlePassDepth
Definition variables.h:234
PetscInt64 maxTraversalSteps
Definition variables.h:228
PetscScalar z
Definition variables.h:101
Vec Ucat
Definition variables.h:837
char ** profilingSelectedFuncs
Definition variables.h:776
PetscInt particlesLostCumulative
Definition variables.h:747
PetscInt nProfilingSelectedFuncs
Definition variables.h:777
PetscInt particlesMigratedLastStep
Definition variables.h:749
char AnalyticalSolutionType[PETSC_MAX_PATH_LEN]
Definition variables.h:676
PetscInt particleConsoleOutputFreq
Definition variables.h:656
SearchMetricsState searchMetrics
Definition variables.h:752
Vec lUcont
Definition variables.h:837
PetscInt step
Definition variables.h:651
DMDALocalInfo info
Definition variables.h:818
Vec lUcat
Definition variables.h:837
PetscInt migrationPassesLastStep
Definition variables.h:748
PetscScalar y
Definition variables.h:101
@ EXEC_MODE_SOLVER
Definition variables.h:616
PetscInt64 searchAttempts
Definition variables.h:222
ExecutionMode exec_mode
Definition variables.h:662
PetscInt64 tieBreakCount
Definition variables.h:230
PetscReal ti
Definition variables.h:652
PetscReal summationRHS
Definition variables.h:770
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:751
@ 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:643
User-defined context containing data specific to a single computational grid level.
Definition variables.h:811