20static PetscErrorCode
AssertFileContains(
const char *path,
const char *needle,
const char *context)
26 PetscFunctionBeginUser;
27 fp = fopen(path,
"rb");
29 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN,
"Failed to open '%s' for assertion.", path);
31 if (fseek(fp, 0, SEEK_END) != 0) {
33 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN,
"Failed to seek '%s'.", path);
35 file_size = ftell(fp);
38 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN,
"Failed to measure '%s'.", path);
40 if (fseek(fp, 0, SEEK_SET) != 0) {
42 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN,
"Failed to rewind '%s'.", path);
45 PetscCall(PetscMalloc1((
size_t)file_size + 1, &buffer));
47 size_t bytes_read = fread(buffer, 1, (
size_t)file_size, fp);
48 if (bytes_read != (
size_t)file_size) {
50 PetscCall(PetscFree(buffer));
51 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN,
"Failed to read '%s'.", path);
54 buffer[file_size] =
'\0';
57 PetscCall(
PicurvAssertBool((PetscBool)(strstr(buffer, needle) != NULL), context));
58 PetscCall(PetscFree(buffer));
59 PetscFunctionReturn(0);
71 PetscFunctionBeginUser;
72 fp = fopen(path,
"rb");
74 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN,
"Failed to open '%s' for assertion.", path);
76 if (fseek(fp, 0, SEEK_END) != 0) {
78 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN,
"Failed to seek '%s'.", path);
80 file_size = ftell(fp);
83 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN,
"Failed to measure '%s'.", path);
85 if (fseek(fp, 0, SEEK_SET) != 0) {
87 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN,
"Failed to rewind '%s'.", path);
90 PetscCall(PetscMalloc1((
size_t)file_size + 1, &buffer));
92 size_t bytes_read = fread(buffer, 1, (
size_t)file_size, fp);
93 if (bytes_read != (
size_t)file_size) {
95 PetscCall(PetscFree(buffer));
96 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN,
"Failed to read '%s'.", path);
99 buffer[file_size] =
'\0';
102 PetscCall(
PicurvAssertBool((PetscBool)(strstr(buffer, needle) == NULL), context));
103 PetscCall(PetscFree(buffer));
104 PetscFunctionReturn(0);
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;
129 size_t bytes_read = 0;
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.");
138 PetscCall(PetscSNPrintf(capture_path,
sizeof(capture_path),
"%s/logging.out", tmpdir));
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.");
149 ierr = fn(user, simCtx, ctx);
151 PetscCheck(dup2(saved_stdout, STDOUT_FILENO) >= 0, PETSC_COMM_SELF, PETSC_ERR_SYS,
"dup2() failed while restoring stdout.");
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);
162 PetscFunctionReturn(0);
170 PetscInt print_interval = *((PetscInt *)ctx);
172 PetscFunctionBeginUser;
175 PetscFunctionReturn(0);
183 PetscInt step = *((PetscInt *)ctx);
185 PetscFunctionBeginUser;
187 PetscFunctionReturn(0);
197 PetscFunctionBeginUser;
200 PetscFunctionReturn(0);
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;
215 PetscFunctionBeginUser;
218 (*simCtx_out)->np = 2;
219 (*simCtx_out)->particleConsoleOutputFreq = 2;
220 (*simCtx_out)->LoggingFrequency = 1;
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));
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;
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);
256 PetscFunctionBeginUser;
258 "BCFaceToString should report the negative-x face"));
260 "FieldInitializationToString should report the zero mode"));
262 "FieldInitializationToString should reject unknown selectors"));
264 "ParticleInitializationToString should report the volume mode"));
266 "LESModelToString should report the constant model"));
268 "MomentumSolverTypeToString should report the explicit solver"));
270 "BCTypeToString should report periodic boundaries"));
272 "BCHandlerTypeToString should report the driven periodic handler"));
274 "ParticleLocationStatusToString should report LOST state"));
275 PetscFunctionReturn(0);
283 PetscFunctionBeginUser;
285 "get_log_level should honor LOG_LEVEL=INFO in this test binary"));
287 PetscFunctionReturn(0);
295 const char *allow_list[] = {
"ComputeSpecificKE",
"WriteEulerianFile"};
297 PetscFunctionBeginUser;
300 "Allowed list should include ComputeSpecificKE"));
302 "Allowed list should exclude unknown function names"));
306 "Empty allow-list should permit all functions"));
307 PetscFunctionReturn(0);
317 PetscFunctionBeginUser;
318 PetscCall(PetscMemzero(&simCtx,
sizeof(simCtx)));
323 "Particle snapshot contract should be enabled when particles and cadence are configured"));
325 "Snapshot should emit on cadence-aligned completed steps"));
327 "Snapshot should not emit off-cadence"));
331 "Zero cadence should disable periodic particle snapshots"));
333 "NULL SimCtx should never emit periodic snapshots"));
334 PetscFunctionReturn(0);
342 char tmpdir[PETSC_MAX_PATH_LEN];
343 char allow_path[PETSC_MAX_PATH_LEN];
344 char dual_log_path[PETSC_MAX_PATH_LEN];
349 PetscReal distances[6] = {1.0, 2.0, 3.0, 4.0, 5.0, 6.0};
353 PetscFunctionBeginUser;
354 PetscCall(PetscMemzero(&cell,
sizeof(cell)));
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));
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);
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"));
377 for (PetscInt i = 0; i < 8; ++i) {
385 PetscCall(PetscCalloc1(1, &monctx));
387 PetscCheck(monctx->
file_handle != NULL, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN,
"Failed to create dual-monitor log '%s'.", dual_log_path);
391 "DualMonitorDestroy should clear the caller-owned context pointer"));
396 PetscCall(PetscPrintf(PETSC_COMM_SELF,
"\n"));
399 PetscFunctionReturn(0);
409 char tmpdir[PETSC_MAX_PATH_LEN];
410 char continuity_path[PETSC_MAX_PATH_LEN];
411 PetscReal ***p = NULL;
414 PetscErrorCode ierr_minmax = 0;
416 PetscFunctionBeginUser;
419 PetscCall(PetscStrncpy(simCtx->
log_dir, tmpdir,
sizeof(simCtx->
log_dir)));
440 PetscCall(PetscSNPrintf(continuity_path,
sizeof(continuity_path),
"%s/Continuity_Metrics.log", simCtx->
log_dir));
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"));
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);
462 PetscCall(DMDAVecRestoreArray(user->
fda, user->
Ucont, &ucont));
463 PetscCall(DMDAVecRestoreArray(user->
fda, user->
Ucat, &ucat));
464 PetscCall(DMDAVecRestoreArray(user->
da, user->
P, &p));
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));
478 PetscCall(PetscPushErrorHandler(PetscIgnoreErrorHandler, NULL));
480 PetscCall(PetscPopErrorHandler());
482 "LOG_FIELD_MIN_MAX should reject unknown field names"));
486 PetscFunctionReturn(0);
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;
502 PetscFunctionBeginUser;
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));
514 PetscCall(DMSwarmCreateGlobalVectorFromField(user->
swarm,
"position", &position_vec));
515 PetscCall(VecDuplicate(position_vec, &analytical_vec));
516 PetscCall(VecCopy(position_vec, analytical_vec));
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]);
526 PetscCall(VecRestoreArrayRead(analytical_vec, &analytical_arr));
527 PetscCall(DMSwarmRestoreField(user->
swarm,
"velocity", NULL, NULL, (
void *)&vel_arr));
529 PetscCall(VecDestroy(&analytical_vec));
530 PetscCall(DMSwarmDestroyGlobalVectorFromField(user->
swarm,
"position", &position_vec));
533 PetscFunctionReturn(0);
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;
551 PetscFunctionBeginUser;
555 PetscCall(PetscStrncpy(simCtx->
log_dir, tmpdir,
sizeof(simCtx->
log_dir)));
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));
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;
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));
597 PetscCall(PetscSNPrintf(metrics_path,
sizeof(metrics_path),
"%s/scatter_metrics.csv", simCtx->
log_dir));
600 "Scatter metrics CSV header should include relative_L2_error"));
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"));
607 PetscFunctionReturn(0);
617 PetscInt print_interval = 1;
620 PetscFunctionBeginUser;
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"));
628 PetscFunctionReturn(0);
641 PetscFunctionBeginUser;
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"));
649 PetscFunctionReturn(0);
659 char tmpdir[PETSC_MAX_PATH_LEN];
660 char metrics_path[PETSC_MAX_PATH_LEN];
662 PetscFunctionBeginUser;
665 PetscCall(PetscStrncpy(simCtx->
log_dir, tmpdir,
sizeof(simCtx->
log_dir)));
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"));
685 PetscFunctionReturn(0);
695 char tmpdir[PETSC_MAX_PATH_LEN];
696 char metrics_path[PETSC_MAX_PATH_LEN];
698 PetscFunctionBeginUser;
701 PetscCall(PetscStrncpy(simCtx->
log_dir, tmpdir,
sizeof(simCtx->
log_dir)));
725 PetscCall(PetscSNPrintf(metrics_path,
sizeof(metrics_path),
"%s/search_metrics.csv", simCtx->
log_dir));
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"));
736 PetscFunctionReturn(0);
749 PetscFunctionBeginUser;
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));
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"));
760 PetscFunctionReturn(0);
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};
775 PetscFunctionBeginUser;
776 PetscCall(PetscMemzero(&simCtx,
sizeof(simCtx)));
778 PetscCall(PetscStrncpy(simCtx.
log_dir, tmpdir,
sizeof(simCtx.
log_dir)));
803 PetscCall(PetscSNPrintf(summary_path,
sizeof(summary_path),
"%s/ProfilingSummary_Solver.log", simCtx.
log_dir));
807 "profiling timestep summary should contain the CSV header"));
809 "profiling timestep summary should log selected functions"));
811 "profiling timestep summary should omit unselected functions in selected mode"));
813 "profiling final summary should include its table banner"));
815 "profiling final summary should include selected functions"));
817 "profiling final summary should include total-time entries for unselected functions"));
819 PetscFunctionReturn(0);
845 (void)setenv(
"LOG_LEVEL",
"INFO", 1);
847 ierr = PetscInitialize(&argc, &argv, NULL,
"PICurv logging tests");
852 ierr =
PicurvRunTests(
"unit-logging", cases,
sizeof(cases) /
sizeof(cases[0]));
860 ierr = PetscFinalize();
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.
PetscErrorCode LOG_PARTICLE_METRICS(UserCtx *user, const char *stageName)
Logs particle swarm metrics, adapting its behavior based on a boolean flag in SimCtx.
const char * BCHandlerTypeToString(BCHandlerType handler_type)
Converts a BCHandlerType enum to its string representation.
PetscBool is_function_allowed(const char *functionName)
Checks if a given function is in the allow-list.
const char * FieldInitializationToString(PetscInt FieldInitialization)
Helper function to convert FieldInitialization to a string representation.
PetscErrorCode DualMonitorDestroy(void **ctx)
Destroys the DualMonitorCtx.
PetscErrorCode LOG_INTERPOLATION_ERROR(UserCtx *user)
Logs the interpolation error between the analytical and computed solutions.
PetscBool ShouldEmitPeriodicParticleConsoleSnapshot(const SimCtx *simCtx, PetscInt completed_step)
Returns whether a particle console snapshot should be emitted for the.
const char * BCFaceToString(BCFace face)
Helper function to convert BCFace enum to a string representation.
PetscErrorCode FreeAllowedFunctions(char **funcs, PetscInt n)
Free an array previously returned by LoadAllowedFunctionsFromFile().
PetscBool IsParticleConsoleSnapshotEnabled(const SimCtx *simCtx)
Returns whether periodic particle console snapshots are enabled.
PetscErrorCode print_log_level(void)
Prints the current logging level to the console.
PetscErrorCode EmitParticleConsoleSnapshot(UserCtx *user, SimCtx *simCtx, PetscInt step)
Emits one particle console snapshot into the main solver log.
PetscErrorCode ProfilingFinalize(SimCtx *simCtx)
the profiling excercise and build a profiling summary which is then printed to a log file.
PetscErrorCode LoadAllowedFunctionsFromFile(const char filename[], char ***funcsOut, PetscInt *nOut)
Load function names from a text file.
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.
void PrintProgressBar(PetscInt step, PetscInt startStep, PetscInt totalSteps, PetscReal currentTime)
Prints a progress bar to the console.
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...
LogLevel get_log_level()
Retrieves the current logging level from the environment variable LOG_LEVEL.
PetscErrorCode ProfilingLogTimestepSummary(SimCtx *simCtx, PetscInt step)
Logs the performance summary for the current timestep and resets timers.
PetscErrorCode LOG_FACE_DISTANCES(PetscReal *d)
Prints the signed distances to each face of the cell.
PetscErrorCode LOG_PARTICLE_FIELDS(UserCtx *user, PetscInt printInterval)
Prints particle fields in a table that automatically adjusts its column widths.
void _ProfilingEnd(const char *func_name)
Internal profiling hook invoked by PROFILE_FUNCTION_END.
const char * BCTypeToString(BCType type)
Helper function to convert BCType enum to a string representation.
PetscErrorCode CalculateAdvancedParticleMetrics(UserCtx *user)
Computes advanced particle statistics and stores them in SimCtx.
const char * ParticleLocationStatusToString(ParticleLocationStatus level)
A function that outputs the name of the current level in the ParticleLocation enum.
PetscErrorCode LOG_SCATTER_METRICS(UserCtx *user)
Logs particle-to-grid scatter verification metrics for the prescribed scalar truth path.
PetscErrorCode LOG_CONTINUITY_METRICS(UserCtx *user)
Logs continuity metrics for a single block to a file.
PetscErrorCode LOG_SEARCH_METRICS(UserCtx *user)
Writes compact runtime search metrics to CSV and optionally to console.
PetscErrorCode ProfilingInitialize(SimCtx *simCtx)
Initializes the custom profiling system using configuration from SimCtx.
@ LOG_INFO
Informational messages about program execution.
const char * LESModelToString(LESModelType LESFlag)
Helper function to convert LES Flag to a string representation.
PetscErrorCode LOG_CELL_VERTICES(const Cell *cell, PetscMPIInt rank)
Prints the coordinates of a cell's vertices.
PetscErrorCode ProfilingResetTimestepCounters(void)
Resets per-timestep profiling counters for the next solver step.
const char * MomentumSolverTypeToString(MomentumSolverType SolverFlag)
Helper function to convert Momentum Solver flag to a string representation.
const char * ParticleInitializationToString(ParticleInitializationType ParticleInitialization)
Helper function to convert ParticleInitialization to a string representation.
void _ProfilingStart(const char *func_name)
Internal profiling hook invoked by PROFILE_FUNCTION_BEGIN.
Context for a dual-purpose KSP monitor.
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.
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.
PetscBool profilingFinalSummary
char profilingTimestepFile[PETSC_MAX_PATH_LEN]
PetscInt64 searchLocatedCount
PetscInt64 searchLostCount
@ PARTICLE_INIT_VOLUME
Random volumetric distribution across the domain.
PetscInt64 boundaryClampCount
PetscInt particlesLostLastStep
PetscInt64 traversalStepsSum
@ BC_HANDLER_PERIODIC_DRIVEN_CONSTANT_FLUX
PetscInt64 searchPopulation
char profilingTimestepMode[32]
@ MOMENTUM_SOLVER_EXPLICIT_RK
PetscInt64 bboxGuessFallbackCount
VerificationScalarConfig verificationScalar
PetscInt64 bboxGuessSuccessCount
char log_dir[PETSC_MAX_PATH_LEN]
PetscInt64 maxParticlePassDepth
PetscInt64 maxTraversalSteps
char ** profilingSelectedFuncs
PetscInt particlesLostCumulative
PetscInt nProfilingSelectedFuncs
PetscInt particlesMigratedLastStep
char AnalyticalSolutionType[PETSC_MAX_PATH_LEN]
PetscInt particleConsoleOutputFreq
SearchMetricsState searchMetrics
PetscInt migrationPassesLastStep
PetscInt64 searchAttempts
PetscInt64 maxTraversalFailCount
Cmpnts vertices[8]
Coordinates of the eight vertices of the cell.
PetscReal particleLoadImbalance
Defines the vertices of a single hexahedral grid cell.
A 3D point or vector with PetscScalar components.
The master context for the entire simulation.
User-defined context containing data specific to a single computational grid level.