25 PetscFunctionBeginUser;
26 PetscCall(DMGlobalToLocalBegin(user->
da, user->
P, INSERT_VALUES, user->
lP));
27 PetscCall(DMGlobalToLocalEnd(user->
da, user->
P, INSERT_VALUES, user->
lP));
28 PetscCall(DMGlobalToLocalBegin(user->
da, user->
Psi, INSERT_VALUES, user->
lPsi));
29 PetscCall(DMGlobalToLocalEnd(user->
da, user->
Psi, INSERT_VALUES, user->
lPsi));
32 PetscCall(DMGlobalToLocalBegin(user->
fda, user->
Ucat, INSERT_VALUES, user->
lUcat));
33 PetscCall(DMGlobalToLocalEnd(user->
fda, user->
Ucat, INSERT_VALUES, user->
lUcat));
34 PetscCall(DMGlobalToLocalBegin(user->
fda, user->
Ucont, INSERT_VALUES, user->
lUcont));
35 PetscCall(DMGlobalToLocalEnd(user->
fda, user->
Ucont, INSERT_VALUES, user->
lUcont));
38 PetscFunctionReturn(0);
54 PetscInt status_value)
56 PetscReal *positions = NULL;
57 PetscReal *weights = NULL;
58 PetscInt *cell_ids = NULL;
59 PetscInt *status = NULL;
61 PetscFunctionBeginUser;
62 PetscCall(DMSwarmGetField(user->
swarm,
"position", NULL, NULL, (
void **)&positions));
63 PetscCall(DMSwarmGetField(user->
swarm,
"weight", NULL, NULL, (
void **)&weights));
64 PetscCall(DMSwarmGetField(user->
swarm,
"DMSwarm_CellID", NULL, NULL, (
void **)&cell_ids));
65 PetscCall(DMSwarmGetField(user->
swarm,
"DMSwarm_location_status", NULL, NULL, (
void **)&status));
76 status[0] = status_value;
78 PetscCall(DMSwarmRestoreField(user->
swarm,
"DMSwarm_location_status", NULL, NULL, (
void **)&status));
79 PetscCall(DMSwarmRestoreField(user->
swarm,
"DMSwarm_CellID", NULL, NULL, (
void **)&cell_ids));
80 PetscCall(DMSwarmRestoreField(user->
swarm,
"weight", NULL, NULL, (
void **)&weights));
81 PetscCall(DMSwarmRestoreField(user->
swarm,
"position", NULL, NULL, (
void **)&positions));
82 PetscFunctionReturn(0);
90 PetscInt particles_per_rank = -1;
91 PetscInt remainder = -1;
93 PetscFunctionBeginUser;
95 PetscCall(
PicurvAssertIntEqual(4, particles_per_rank,
"rank 0 should receive one remainder particle"));
99 PetscCall(
PicurvAssertIntEqual(3, particles_per_rank,
"last rank should receive base particle count"));
100 PetscCall(
PicurvAssertIntEqual(1, remainder,
"remainder should remain unchanged across ranks"));
101 PetscFunctionReturn(0);
112 PetscFunctionBeginUser;
113 PetscCall(PetscMemzero(&bbox,
sizeof(bbox)));
114 PetscCall(PetscMemzero(&particle,
sizeof(particle)));
123 particle.
loc.
x = 0.25;
124 particle.
loc.
y = 1.0;
125 particle.
loc.
z = 2.5;
128 particle.
loc.
x = 1.5;
130 PetscFunctionReturn(0);
139 PetscReal distances[
NUM_FACES] = {1.0, 3.0, 2.0, 2.0, 4.0, 1.0};
140 PetscReal clamped[
NUM_FACES] = {0.0, 0.0, 0.0, 0.0, 0.0, 0.0};
142 PetscFunctionBeginUser;
143 PetscCall(PetscMemzero(&particle,
sizeof(particle)));
153 PetscFunctionReturn(0);
165 PetscFunctionBeginUser;
169 PetscCall(PetscMemzero(&particle,
sizeof(particle)));
170 particle.
loc.
x = 1.0;
171 particle.
loc.
y = -2.0;
172 particle.
loc.
z = 3.0;
173 particle.
vel.
x = 0.5;
174 particle.
vel.
y = -1.0;
175 particle.
vel.
z = 2.0;
187 PetscFunctionReturn(0);
199 PetscFunctionBeginUser;
203 PetscCall(PetscMemzero(&particle,
sizeof(particle)));
204 particle.
loc.
x = 0.25;
205 particle.
loc.
y = 0.5;
206 particle.
loc.
z = 0.75;
207 particle.
vel.
x = 0.0;
208 particle.
vel.
y = 0.0;
209 particle.
vel.
z = 0.0;
221 PetscFunctionReturn(0);
230 const PetscReal dt = 0.5;
231 const PetscReal diffusivity = 0.2;
232 const PetscReal mean_val = 3.0;
233 const PetscReal cell_vol = 8.0;
234 const PetscReal c_model = 2.0;
235 PetscReal unchanged = 7.0;
237 PetscFunctionBeginUser;
240 mean_val + (1.0 - mean_val) * PetscExpReal(-(c_model * diffusivity / PetscPowReal(cell_vol, 0.6666667)) * dt),
243 "IEM update should match analytical relaxation"));
245 PetscCall(
UpdateParticleField(
"UnrelatedField", dt, &unchanged, diffusivity, mean_val, cell_vol, c_model));
246 PetscCall(
PicurvAssertRealNear(7.0, unchanged, 1.0e-12,
"unknown field should remain unchanged"));
247 PetscFunctionReturn(0);
258 PetscFunctionBeginUser;
260 PetscCall(VecSet(user->
Ucont, 7.0));
266 PetscFunctionReturn(0);
278 PetscFunctionBeginUser;
286 PetscCall(VecSet(user->
Ucont, 0.0));
289 PetscCall(DMDAVecGetArrayRead(user->
fda, user->
Ucont, &ucont));
290 PetscCall(
PicurvAssertRealNear(0.0, ucont[1][1][1].x, 1.0e-12,
"x contravariant component remains zero"));
291 PetscCall(
PicurvAssertRealNear(0.0, ucont[1][1][1].y, 1.0e-12,
"y contravariant component remains zero"));
292 PetscCall(
PicurvAssertRealNear(2.0, ucont[1][1][1].z, 1.0e-12,
"z contravariant component should match constant inlet magnitude"));
293 PetscCall(
PicurvAssertRealNear(0.0, ucont[0][1][1].z, 1.0e-12,
"boundary/ghost cell should remain unchanged"));
294 PetscCall(DMDAVecRestoreArrayRead(user->
fda, user->
Ucont, &ucont));
297 PetscFunctionReturn(0);
308 PetscReal *velocity = NULL;
309 PetscReal *diffusivity = NULL;
310 PetscReal *diffusivity_gradient = NULL;
312 PetscFunctionBeginUser;
315 PetscCall(VecSet(user->
Ucat, 2.0));
319 for (PetscInt k = user->
info.zs; k < user->
info.zs + user->
info.zm; ++k) {
320 for (PetscInt j = user->
info.ys; j < user->
info.ys + user->
info.ym; ++j) {
321 for (PetscInt i = user->
info.xs; i < user->
info.xs + user->
info.xm; ++i) {
322 grad[k][j][i].
x = 0.1;
323 grad[k][j][i].
y = 0.2;
324 grad[k][j][i].
z = 0.3;
330 PetscCall(
SeedSingleParticle(user, 0, 0, 0, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5,
ACTIVE_AND_LOCATED));
334 PetscCall(DMSwarmGetField(user->
swarm,
"velocity", NULL, NULL, (
void **)&velocity));
335 PetscCall(DMSwarmGetField(user->
swarm,
"Diffusivity", NULL, NULL, (
void **)&diffusivity));
336 PetscCall(DMSwarmGetField(user->
swarm,
"DiffusivityGradient", NULL, NULL, (
void **)&diffusivity_gradient));
337 PetscCall(
PicurvAssertRealNear(2.0, velocity[0], 1.0e-12,
"Interpolated velocity x should match constant Eulerian field"));
338 PetscCall(
PicurvAssertRealNear(2.0, velocity[1], 1.0e-12,
"Interpolated velocity y should match constant Eulerian field"));
339 PetscCall(
PicurvAssertRealNear(2.0, velocity[2], 1.0e-12,
"Interpolated velocity z should match constant Eulerian field"));
340 PetscCall(
PicurvAssertRealNear(0.25, diffusivity[0], 1.0e-12,
"Interpolated scalar diffusivity should match constant Eulerian field"));
341 PetscCall(
PicurvAssertRealNear(0.1, diffusivity_gradient[0], 1.0e-12,
"Interpolated diffusivity-gradient x component"));
342 PetscCall(
PicurvAssertRealNear(0.2, diffusivity_gradient[1], 1.0e-12,
"Interpolated diffusivity-gradient y component"));
343 PetscCall(
PicurvAssertRealNear(0.3, diffusivity_gradient[2], 1.0e-12,
"Interpolated diffusivity-gradient z component"));
344 PetscCall(DMSwarmRestoreField(user->
swarm,
"DiffusivityGradient", NULL, NULL, (
void **)&diffusivity_gradient));
345 PetscCall(DMSwarmRestoreField(user->
swarm,
"Diffusivity", NULL, NULL, (
void **)&diffusivity));
346 PetscCall(DMSwarmRestoreField(user->
swarm,
"velocity", NULL, NULL, (
void **)&velocity));
349 PetscFunctionReturn(0);
360 PetscReal *velocity = NULL;
361 PetscReal *diffusivity = NULL;
362 PetscReal *diffusivity_gradient = NULL;
364 PetscFunctionBeginUser;
368 PetscCall(VecSet(user->
Ucat, 2.0));
372 for (PetscInt k = user->
info.zs; k < user->
info.zs + user->
info.zm; ++k) {
373 for (PetscInt j = user->
info.ys; j < user->
info.ys + user->
info.ym; ++j) {
374 for (PetscInt i = user->
info.xs; i < user->
info.xs + user->
info.xm; ++i) {
375 grad[k][j][i].
x = 0.1;
376 grad[k][j][i].
y = 0.2;
377 grad[k][j][i].
z = 0.3;
383 PetscCall(
SeedSingleParticle(user, 0, 0, 0, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5,
ACTIVE_AND_LOCATED));
387 PetscCall(DMSwarmGetField(user->
swarm,
"velocity", NULL, NULL, (
void **)&velocity));
388 PetscCall(DMSwarmGetField(user->
swarm,
"Diffusivity", NULL, NULL, (
void **)&diffusivity));
389 PetscCall(DMSwarmGetField(user->
swarm,
"DiffusivityGradient", NULL, NULL, (
void **)&diffusivity_gradient));
390 PetscCall(
PicurvAssertRealNear(2.0, velocity[0], 1.0e-12,
"CornerAveraged: interpolated velocity x should match constant Eulerian field"));
391 PetscCall(
PicurvAssertRealNear(2.0, velocity[1], 1.0e-12,
"CornerAveraged: interpolated velocity y should match constant Eulerian field"));
392 PetscCall(
PicurvAssertRealNear(2.0, velocity[2], 1.0e-12,
"CornerAveraged: interpolated velocity z should match constant Eulerian field"));
393 PetscCall(
PicurvAssertRealNear(0.25, diffusivity[0], 1.0e-12,
"CornerAveraged: interpolated scalar diffusivity should match constant Eulerian field"));
394 PetscCall(
PicurvAssertRealNear(0.1, diffusivity_gradient[0], 1.0e-12,
"CornerAveraged: interpolated diffusivity-gradient x component"));
395 PetscCall(
PicurvAssertRealNear(0.2, diffusivity_gradient[1], 1.0e-12,
"CornerAveraged: interpolated diffusivity-gradient y component"));
396 PetscCall(
PicurvAssertRealNear(0.3, diffusivity_gradient[2], 1.0e-12,
"CornerAveraged: interpolated diffusivity-gradient z component"));
397 PetscCall(DMSwarmRestoreField(user->
swarm,
"DiffusivityGradient", NULL, NULL, (
void **)&diffusivity_gradient));
398 PetscCall(DMSwarmRestoreField(user->
swarm,
"Diffusivity", NULL, NULL, (
void **)&diffusivity));
399 PetscCall(DMSwarmRestoreField(user->
swarm,
"velocity", NULL, NULL, (
void **)&velocity));
402 PetscFunctionReturn(0);
412 PetscInt *cell_ids = NULL;
413 PetscReal *psi = NULL;
414 PetscReal ***psi_grid = NULL;
416 PetscFunctionBeginUser;
420 PetscCall(DMSwarmGetField(user->
swarm,
"DMSwarm_CellID", NULL, NULL, (
void **)&cell_ids));
421 cell_ids[0] = 0; cell_ids[1] = 0; cell_ids[2] = 0;
422 cell_ids[3] = 0; cell_ids[4] = 0; cell_ids[5] = 0;
423 PetscCall(DMSwarmRestoreField(user->
swarm,
"DMSwarm_CellID", NULL, NULL, (
void **)&cell_ids));
425 PetscCall(DMSwarmGetField(user->
swarm,
"Psi", NULL, NULL, (
void **)&psi));
428 PetscCall(DMSwarmRestoreField(user->
swarm,
"Psi", NULL, NULL, (
void **)&psi));
432 PetscCall(DMDAVecGetArrayRead(user->
da, user->
Psi, &psi_grid));
433 PetscCall(
PicurvAssertRealNear(2.0, psi_grid[1][1][1], 1.0e-12,
"Scatter should average particle Psi values into the owning cell"));
434 PetscCall(DMDAVecRestoreArrayRead(user->
da, user->
Psi, &psi_grid));
437 PetscFunctionReturn(0);
447 PetscInt *cell_ids = NULL;
448 PetscReal ***counts = NULL;
450 PetscFunctionBeginUser;
454 PetscCall(DMSwarmGetField(user->
swarm,
"DMSwarm_CellID", NULL, NULL, (
void **)&cell_ids));
455 cell_ids[0] = 0; cell_ids[1] = 0; cell_ids[2] = 0;
456 cell_ids[3] = 0; cell_ids[4] = 0; cell_ids[5] = 0;
457 cell_ids[6] = 1; cell_ids[7] = 0; cell_ids[8] = 0;
458 PetscCall(DMSwarmRestoreField(user->
swarm,
"DMSwarm_CellID", NULL, NULL, (
void **)&cell_ids));
462 PetscCall(DMDAVecGetArrayRead(user->
da, user->
ParticleCount, &counts));
463 PetscCall(
PicurvAssertRealNear(2.0, counts[1][1][1], 1.0e-12,
"Two particles in cell (0,0,0) should accumulate at shifted index (1,1,1)"));
464 PetscCall(
PicurvAssertRealNear(1.0, counts[1][1][2], 1.0e-12,
"One particle in cell (1,0,0) should accumulate at shifted index (2,1,1)"));
465 PetscCall(DMDAVecRestoreArrayRead(user->
da, user->
ParticleCount, &counts));
468 PetscFunctionReturn(0);
478 PetscInt *status = NULL;
480 PetscFunctionBeginUser;
484 PetscCall(DMSwarmGetField(user->
swarm,
"DMSwarm_location_status", NULL, NULL, (
void **)&status));
488 PetscCall(DMSwarmRestoreField(user->
swarm,
"DMSwarm_location_status", NULL, NULL, (
void **)&status));
492 PetscCall(DMSwarmGetField(user->
swarm,
"DMSwarm_location_status", NULL, NULL, (
void **)&status));
496 PetscCall(DMSwarmRestoreField(user->
swarm,
"DMSwarm_location_status", NULL, NULL, (
void **)&status));
499 PetscFunctionReturn(0);
509 PetscReal *positions = NULL;
510 PetscInt removed_local = 0;
511 PetscInt removed_global = 0;
514 PetscFunctionBeginUser;
518 PetscCall(DMSwarmGetField(user->
swarm,
"position", NULL, NULL, (
void **)&positions));
519 positions[0] = 0.5; positions[1] = 0.5; positions[2] = 0.5;
520 positions[3] = 9.0; positions[4] = 9.0; positions[5] = 9.0;
521 PetscCall(DMSwarmRestoreField(user->
swarm,
"position", NULL, NULL, (
void **)&positions));
524 PetscCall(DMSwarmGetLocalSize(user->
swarm, &nlocal));
525 PetscCall(
PicurvAssertIntEqual(1, removed_local,
"Exactly one particle should be removed as out-of-bounds on a single rank"));
526 PetscCall(
PicurvAssertIntEqual(1, removed_global,
"Global out-of-bounds removal count should match the local single-rank result"));
527 PetscCall(
PicurvAssertIntEqual(1, nlocal,
"One in-bounds particle should remain after out-of-bounds removal"));
530 PetscFunctionReturn(0);
540 PetscInt *status = NULL;
541 PetscInt removed_local = 0;
542 PetscInt removed_global = 0;
545 PetscFunctionBeginUser;
549 PetscCall(DMSwarmGetField(user->
swarm,
"DMSwarm_location_status", NULL, NULL, (
void **)&status));
553 PetscCall(DMSwarmRestoreField(user->
swarm,
"DMSwarm_location_status", NULL, NULL, (
void **)&status));
556 PetscCall(DMSwarmGetLocalSize(user->
swarm, &nlocal));
557 PetscCall(
PicurvAssertIntEqual(2, removed_local,
"Two LOST particles should be removed locally"));
558 PetscCall(
PicurvAssertIntEqual(2, removed_global,
"Global LOST-particle removal count should match the local single-rank result"));
562 PetscFunctionReturn(0);
572 char tmpdir[PETSC_MAX_PATH_LEN];
576 PetscFunctionBeginUser;
580 "runtime setup path should initialize the Brownian RNG"));
590 PetscCall(
PicurvAssertRealNear(first.
x, second.
x, 1.0e-12,
"Resetting the Brownian RNG seed should reproduce the x displacement"));
591 PetscCall(
PicurvAssertRealNear(first.
y, second.
y, 1.0e-12,
"Resetting the Brownian RNG seed should reproduce the y displacement"));
592 PetscCall(
PicurvAssertRealNear(first.
z, second.
z, 1.0e-12,
"Resetting the Brownian RNG seed should reproduce the z displacement"));
596 PetscFunctionReturn(0);
605 PetscReal *positions = NULL;
606 PetscReal *velocities = NULL;
607 PetscReal *diffusivity = NULL;
608 Cmpnts *diffusivity_gradient = NULL;
609 PetscReal *psi = NULL;
610 PetscReal *weights = NULL;
611 PetscInt *cell_ids = NULL;
612 PetscInt *status = NULL;
613 PetscInt64 *pid = NULL;
615 PetscFunctionBeginUser;
620 PetscCall(DMSwarmGetField(user->
swarm,
"position", NULL, NULL, (
void **)&positions));
621 PetscCall(DMSwarmGetField(user->
swarm,
"velocity", NULL, NULL, (
void **)&velocities));
622 PetscCall(DMSwarmGetField(user->
swarm,
"Diffusivity", NULL, NULL, (
void **)&diffusivity));
623 PetscCall(DMSwarmGetField(user->
swarm,
"DiffusivityGradient", NULL, NULL, (
void **)&diffusivity_gradient));
624 PetscCall(DMSwarmGetField(user->
swarm,
"Psi", NULL, NULL, (
void **)&psi));
625 PetscCall(DMSwarmGetField(user->
swarm,
"weight", NULL, NULL, (
void **)&weights));
626 PetscCall(DMSwarmGetField(user->
swarm,
"DMSwarm_CellID", NULL, NULL, (
void **)&cell_ids));
627 PetscCall(DMSwarmGetField(user->
swarm,
"DMSwarm_location_status", NULL, NULL, (
void **)&status));
628 PetscCall(DMSwarmGetField(user->
swarm,
"DMSwarm_pid", NULL, NULL, (
void **)&pid));
630 positions[0] = 0.20; positions[1] = 0.30; positions[2] = 0.40;
631 velocities[0] = 0.40; velocities[1] = -0.20; velocities[2] = 0.10;
632 diffusivity[0] = 0.0;
633 diffusivity_gradient[0].
x = 0.10;
634 diffusivity_gradient[0].
y = 0.20;
635 diffusivity_gradient[0].
z = -0.10;
637 weights[0] = 0.5; weights[1] = 0.5; weights[2] = 0.5;
638 cell_ids[0] = 0; cell_ids[1] = 0; cell_ids[2] = 0;
642 PetscCall(DMSwarmRestoreField(user->
swarm,
"DMSwarm_pid", NULL, NULL, (
void **)&pid));
643 PetscCall(DMSwarmRestoreField(user->
swarm,
"DMSwarm_location_status", NULL, NULL, (
void **)&status));
644 PetscCall(DMSwarmRestoreField(user->
swarm,
"DMSwarm_CellID", NULL, NULL, (
void **)&cell_ids));
645 PetscCall(DMSwarmRestoreField(user->
swarm,
"weight", NULL, NULL, (
void **)&weights));
646 PetscCall(DMSwarmRestoreField(user->
swarm,
"Psi", NULL, NULL, (
void **)&psi));
647 PetscCall(DMSwarmRestoreField(user->
swarm,
"DiffusivityGradient", NULL, NULL, (
void **)&diffusivity_gradient));
648 PetscCall(DMSwarmRestoreField(user->
swarm,
"Diffusivity", NULL, NULL, (
void **)&diffusivity));
649 PetscCall(DMSwarmRestoreField(user->
swarm,
"velocity", NULL, NULL, (
void **)&velocities));
650 PetscCall(DMSwarmRestoreField(user->
swarm,
"position", NULL, NULL, (
void **)&positions));
654 PetscCall(DMSwarmGetField(user->
swarm,
"position", NULL, NULL, (
void **)&positions));
655 PetscCall(
PicurvAssertRealNear(0.325, positions[0], 1.0e-12,
"UpdateAllParticlePositions should advect x"));
656 PetscCall(
PicurvAssertRealNear(0.300, positions[1], 1.0e-12,
"UpdateAllParticlePositions should advect y"));
657 PetscCall(
PicurvAssertRealNear(0.400, positions[2], 1.0e-12,
"UpdateAllParticlePositions should advect z"));
658 PetscCall(DMSwarmRestoreField(user->
swarm,
"position", NULL, NULL, (
void **)&positions));
661 PetscFunctionReturn(0);
670 PetscReal *positions = NULL;
671 PetscReal *weights = NULL;
672 PetscInt *cell_ids = NULL;
673 PetscInt *status = NULL;
674 PetscInt64 *pid = NULL;
676 PetscFunctionBeginUser;
684 PetscCall(DMSwarmGetField(user->
swarm,
"position", NULL, NULL, (
void **)&positions));
685 PetscCall(DMSwarmGetField(user->
swarm,
"weight", NULL, NULL, (
void **)&weights));
686 PetscCall(DMSwarmGetField(user->
swarm,
"DMSwarm_CellID", NULL, NULL, (
void **)&cell_ids));
687 PetscCall(DMSwarmGetField(user->
swarm,
"DMSwarm_location_status", NULL, NULL, (
void **)&status));
688 PetscCall(DMSwarmGetField(user->
swarm,
"DMSwarm_pid", NULL, NULL, (
void **)&pid));
689 positions[0] = 0.375; positions[1] = 0.375; positions[2] = 0.375;
690 weights[0] = 0.5; weights[1] = 0.5; weights[2] = 0.5;
691 cell_ids[0] = 1; cell_ids[1] = 1; cell_ids[2] = 1;
694 PetscCall(DMSwarmRestoreField(user->
swarm,
"DMSwarm_pid", NULL, NULL, (
void **)&pid));
695 PetscCall(DMSwarmRestoreField(user->
swarm,
"DMSwarm_location_status", NULL, NULL, (
void **)&status));
696 PetscCall(DMSwarmRestoreField(user->
swarm,
"DMSwarm_CellID", NULL, NULL, (
void **)&cell_ids));
697 PetscCall(DMSwarmRestoreField(user->
swarm,
"weight", NULL, NULL, (
void **)&weights));
698 PetscCall(DMSwarmRestoreField(user->
swarm,
"position", NULL, NULL, (
void **)&positions));
702 PetscCall(DMSwarmGetField(user->
swarm,
"DMSwarm_CellID", NULL, NULL, (
void **)&cell_ids));
703 PetscCall(DMSwarmGetField(user->
swarm,
"DMSwarm_location_status", NULL, NULL, (
void **)&status));
704 PetscCall(
PicurvAssertIntEqual(1, cell_ids[0],
"prior-cell fast path should preserve the i cell id"));
705 PetscCall(
PicurvAssertIntEqual(1, cell_ids[1],
"prior-cell fast path should preserve the j cell id"));
706 PetscCall(
PicurvAssertIntEqual(1, cell_ids[2],
"prior-cell fast path should preserve the k cell id"));
715 PetscCall(DMSwarmRestoreField(user->
swarm,
"DMSwarm_location_status", NULL, NULL, (
void **)&status));
716 PetscCall(DMSwarmRestoreField(user->
swarm,
"DMSwarm_CellID", NULL, NULL, (
void **)&cell_ids));
719 PetscFunctionReturn(0);
728 PetscReal *positions = NULL;
729 PetscReal *weights = NULL;
730 PetscInt *cell_ids = NULL;
731 PetscInt *status = NULL;
732 PetscInt64 *pid = NULL;
734 PetscFunctionBeginUser;
742 PetscCall(DMSwarmGetField(user->
swarm,
"position", NULL, NULL, (
void **)&positions));
743 PetscCall(DMSwarmGetField(user->
swarm,
"weight", NULL, NULL, (
void **)&weights));
744 PetscCall(DMSwarmGetField(user->
swarm,
"DMSwarm_CellID", NULL, NULL, (
void **)&cell_ids));
745 PetscCall(DMSwarmGetField(user->
swarm,
"DMSwarm_location_status", NULL, NULL, (
void **)&status));
746 PetscCall(DMSwarmGetField(user->
swarm,
"DMSwarm_pid", NULL, NULL, (
void **)&pid));
747 positions[0] = 0.625; positions[1] = 0.625; positions[2] = 0.625;
748 weights[0] = 0.5; weights[1] = 0.5; weights[2] = 0.5;
749 cell_ids[0] = -1; cell_ids[1] = -1; cell_ids[2] = -1;
752 PetscCall(DMSwarmRestoreField(user->
swarm,
"DMSwarm_pid", NULL, NULL, (
void **)&pid));
753 PetscCall(DMSwarmRestoreField(user->
swarm,
"DMSwarm_location_status", NULL, NULL, (
void **)&status));
754 PetscCall(DMSwarmRestoreField(user->
swarm,
"DMSwarm_CellID", NULL, NULL, (
void **)&cell_ids));
755 PetscCall(DMSwarmRestoreField(user->
swarm,
"weight", NULL, NULL, (
void **)&weights));
756 PetscCall(DMSwarmRestoreField(user->
swarm,
"position", NULL, NULL, (
void **)&positions));
760 PetscCall(DMSwarmGetField(user->
swarm,
"DMSwarm_CellID", NULL, NULL, (
void **)&cell_ids));
761 PetscCall(DMSwarmGetField(user->
swarm,
"DMSwarm_location_status", NULL, NULL, (
void **)&status));
762 PetscCall(
PicurvAssertIntEqual(2, cell_ids[0],
"guess-path location should resolve the i cell id"));
763 PetscCall(
PicurvAssertIntEqual(2, cell_ids[1],
"guess-path location should resolve the j cell id"));
764 PetscCall(
PicurvAssertIntEqual(2, cell_ids[2],
"guess-path location should resolve the k cell id"));
773 PetscCall(DMSwarmRestoreField(user->
swarm,
"DMSwarm_location_status", NULL, NULL, (
void **)&status));
774 PetscCall(DMSwarmRestoreField(user->
swarm,
"DMSwarm_CellID", NULL, NULL, (
void **)&cell_ids));
777 PetscFunctionReturn(0);
789 PetscFunctionBeginUser;
790 PetscCall(PetscMemzero(&particle,
sizeof(particle)));
799 particle.
cell[0] = 1;
800 particle.
cell[1] = 1;
801 particle.
cell[2] = 1;
802 particle.
loc.
x = 0.375;
803 particle.
loc.
y = 0.375;
804 particle.
loc.
z = 0.375;
817 PetscFunctionReturn(0);
825 Cmpnts wall_velocity = {0.0, 0.0, 0.0};
826 Cmpnts reference_velocity = {2.0, 4.0, 6.0};
827 Cmpnts boundary_velocity = {0.0, 0.0, 0.0};
828 Cmpnts free_slip_reference = {2.0, 3.0, 4.0};
830 PetscFunctionBeginUser;
831 noslip(NULL, 2.0, 1.0, wall_velocity, reference_velocity, &boundary_velocity, 1.0, 0.0, 0.0);
836 freeslip(NULL, 2.0, 1.0, wall_velocity, free_slip_reference, &boundary_velocity, 1.0, 0.0, 0.0);
837 PetscCall(
PicurvAssertRealNear(1.0, boundary_velocity.
x, 1.0e-12,
"free-slip interpolated normal component"));
838 PetscCall(
PicurvAssertRealNear(3.0, boundary_velocity.
y, 1.0e-12,
"free-slip tangential y preserved"));
839 PetscCall(
PicurvAssertRealNear(4.0, boundary_velocity.
z, 1.0e-12,
"free-slip tangential z preserved"));
840 PetscFunctionReturn(0);
848 const PetscReal expected_smooth_e = PetscExpReal(0.41 * 5.5);
849 PetscReal e_coeff = 0.0;
850 PetscReal utau = 0.0;
851 PetscReal residual = 0.0;
853 PetscFunctionBeginUser;
854 e_coeff =
E_coeff(0.1, 0.0, 1.0e-3);
855 PetscCall(
PicurvAssertRealNear(expected_smooth_e, e_coeff, 1.0e-10,
"smooth-wall E coefficient"));
858 PetscCall(
PicurvAssertBool((PetscBool)(utau > 0.0),
"friction velocity should remain positive"));
859 residual =
f_hydset(1.0e-3, 1.0, 1.0e-2, utau, 0.0);
860 PetscCall(
PicurvAssertBool((PetscBool)(PetscAbsReal(residual) < 1.0e-5),
"Newton solve residual should be small"));
864 PetscFunctionReturn(0);
872 const PetscReal kinematic_viscosity = 1.0e-3;
873 const PetscReal wall_distance = 2.0e-2;
874 const PetscReal target_velocity = 1.0;
875 const PetscReal roughness_length = 1.0e-4;
876 PetscReal utau_loglaw = 0.0;
877 PetscReal utau_werner = 0.0;
878 PetscReal utau_cabot = 0.0;
879 PetscReal wall_shear_velocity = 0.0;
880 PetscReal wall_shear_normal = 0.0;
882 PetscFunctionBeginUser;
883 utau_loglaw =
find_utau_loglaw(target_velocity, wall_distance, roughness_length);
885 "simple log-law inversion should reconstruct the target velocity"));
887 utau_werner =
find_utau_Werner(kinematic_viscosity, target_velocity, wall_distance, 0.1);
888 PetscCall(
PicurvAssertBool((PetscBool)(utau_werner > 0.0),
"Werner-Wengle friction velocity should remain positive"));
890 "Werner-Wengle inversion should reconstruct the target velocity"));
892 find_utau_Cabot(kinematic_viscosity, target_velocity, wall_distance, 0.1, 0.0, 0.0,
893 &utau_cabot, &wall_shear_velocity, &wall_shear_normal);
894 PetscCall(
PicurvAssertBool((PetscBool)(utau_cabot > 0.0),
"Cabot friction velocity should remain positive"));
895 PetscCall(
PicurvAssertRealNear(target_velocity,
u_Cabot(kinematic_viscosity, wall_distance, utau_cabot, 0.0, wall_shear_velocity), 1.0e-6,
896 "Cabot inversion should reconstruct the target velocity when pressure gradient is zero"));
898 "zero normal pressure gradient should keep Cabot normal wall shear at zero"));
899 PetscFunctionReturn(0);
908 Cmpnts wall_velocity = {0.0, 0.0, 0.0};
909 Cmpnts reference_velocity = {0.0, 1.0, 0.0};
910 Cmpnts boundary_velocity = {0.0, 0.0, 0.0};
911 PetscReal friction_velocity = 0.0;
913 PetscFunctionBeginUser;
915 simCtx->
ren = 1000.0;
917 wall_function(user, 2.0e-2, 1.0e-2, wall_velocity, reference_velocity, &boundary_velocity, &friction_velocity, 1.0, 0.0, 0.0);
918 PetscCall(
PicurvAssertRealNear(0.0, boundary_velocity.
x, 1.0e-12,
"Werner wall function should preserve zero normal velocity"));
919 PetscCall(
PicurvAssertBool((PetscBool)(boundary_velocity.
y > 0.0 && boundary_velocity.
y < 1.0),
"Werner wall function should damp tangential velocity"));
920 PetscCall(
PicurvAssertBool((PetscBool)(friction_velocity > 0.0),
"Werner wall function should compute positive friction velocity"));
922 wall_function_loglaw(user, 1.0e-4, 2.0e-2, 1.0e-2, wall_velocity, reference_velocity, &boundary_velocity, &friction_velocity, 1.0, 0.0, 0.0);
923 PetscCall(
PicurvAssertRealNear(0.0, boundary_velocity.
x, 1.0e-12,
"log-law wall function should preserve zero normal velocity"));
924 PetscCall(
PicurvAssertBool((PetscBool)(boundary_velocity.
y > 0.0 && boundary_velocity.
y <= 1.0),
"log-law wall function should keep tangential velocity bounded"));
925 PetscCall(
PicurvAssertBool((PetscBool)(friction_velocity > 0.0),
"log-law wall function should compute positive friction velocity"));
927 wall_function_Cabot(user, 1.0e-4, 2.0e-2, 1.0e-2, wall_velocity, reference_velocity, &boundary_velocity, &friction_velocity,
928 1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 10);
929 PetscCall(
PicurvAssertRealNear(0.0, boundary_velocity.
x, 1.0e-12,
"Cabot wall function should preserve zero normal velocity"));
930 PetscCall(
PicurvAssertBool((PetscBool)(boundary_velocity.
y > 0.0 && boundary_velocity.
y <= 1.0),
"Cabot wall function should keep tangential velocity bounded"));
931 PetscCall(
PicurvAssertBool((PetscBool)(friction_velocity > 0.0),
"Cabot wall function should compute positive friction velocity"));
934 PetscFunctionReturn(0);
945 PetscFunctionBeginUser;
949 PetscFunctionReturn(0);
960 PetscFunctionBeginUser;
962 PetscCall(DMCreateGlobalVector(user->
da, &user->
CS));
972 PetscFunctionReturn(0);
983 PetscInt mx = 0, my = 0, mz = 0;
985 PetscFunctionBeginUser;
988 PetscCall(DMDAGetInfo(user->
da, NULL, &mx, &my, &mz, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
996 PetscCall(DMGetCoordinateDM(user->
da, &coord_dm));
998 "minimal fixture should derive fda from the coordinate-DM path"));
1001 PetscFunctionReturn(0);
1012 void *field_ptr = NULL;
1014 PetscFunctionBeginUser;
1018 PetscCall(DMSwarmGetField(user->
swarm,
"position", &bs, NULL, &field_ptr));
1020 PetscCall(
PicurvAssertBool((PetscBool)(field_ptr != NULL),
"position field should be retrievable"));
1021 PetscCall(DMSwarmRestoreField(user->
swarm,
"position", &bs, NULL, &field_ptr));
1023 PetscCall(DMSwarmGetField(user->
swarm,
"velocity", &bs, NULL, &field_ptr));
1025 PetscCall(DMSwarmRestoreField(user->
swarm,
"velocity", &bs, NULL, &field_ptr));
1027 PetscCall(DMSwarmGetField(user->
swarm,
"DMSwarm_CellID", &bs, NULL, &field_ptr));
1029 PetscCall(DMSwarmRestoreField(user->
swarm,
"DMSwarm_CellID", &bs, NULL, &field_ptr));
1031 PetscCall(DMSwarmGetField(user->
swarm,
"weight", &bs, NULL, &field_ptr));
1033 PetscCall(DMSwarmRestoreField(user->
swarm,
"weight", &bs, NULL, &field_ptr));
1035 PetscCall(DMSwarmGetField(user->
swarm,
"Diffusivity", &bs, NULL, &field_ptr));
1037 PetscCall(DMSwarmRestoreField(user->
swarm,
"Diffusivity", &bs, NULL, &field_ptr));
1039 PetscCall(DMSwarmGetField(user->
swarm,
"DiffusivityGradient", &bs, NULL, &field_ptr));
1040 PetscCall(
PicurvAssertIntEqual(3, bs,
"solver swarm should register particle diffusivity gradients"));
1041 PetscCall(DMSwarmRestoreField(user->
swarm,
"DiffusivityGradient", &bs, NULL, &field_ptr));
1043 PetscCall(DMSwarmGetField(user->
swarm,
"Psi", &bs, NULL, &field_ptr));
1045 PetscCall(DMSwarmRestoreField(user->
swarm,
"Psi", &bs, NULL, &field_ptr));
1047 PetscCall(DMSwarmGetField(user->
swarm,
"DMSwarm_location_status", &bs, NULL, &field_ptr));
1048 PetscCall(
PicurvAssertIntEqual(1, bs,
"solver swarm should register particle location status"));
1049 PetscCall(DMSwarmRestoreField(user->
swarm,
"DMSwarm_location_status", &bs, NULL, &field_ptr));
1052 PetscFunctionReturn(0);
1063 PetscFunctionBeginUser;
1066 PetscCall(VecSet(user->
Ucont, 11.0));
1067 PetscCall(VecSet(user->
Ucont_o, 7.0));
1068 PetscCall(VecSet(user->
Ucont_rm1, 3.0));
1069 PetscCall(VecSet(user->
Ucat, 5.0));
1070 PetscCall(VecSet(user->
Ucat_o, -1.0));
1071 PetscCall(VecSet(user->
P, 9.0));
1072 PetscCall(VecSet(user->
P_o, -2.0));
1084 PetscFunctionReturn(0);
1094 PetscInt xs = -1, ys = -1, zs = -1;
1095 PetscInt xm = -1, ym = -1, zm = -1;
1097 PetscFunctionBeginUser;
1112 PetscFunctionReturn(0);
1123 PetscFunctionBeginUser;
1135 PetscFunctionReturn(0);
1143 PetscReal seconds = 0.0;
1145 PetscFunctionBeginUser;
1152 PetscFunctionReturn(0);
1160 PetscReal ewma_fast = 0.0;
1161 PetscReal ewma_slow = 0.0;
1162 PetscReal conservative_fast = 0.0;
1163 PetscReal conservative_slow = 0.0;
1164 PetscReal required_headroom = 0.0;
1166 PetscFunctionBeginUser;
1175 PetscCall(
PicurvAssertRealNear(6.0, conservative_fast, 1.0e-12,
"conservative estimate tracks latest moderate step"));
1176 PetscCall(
PicurvAssertRealNear(12.0, conservative_slow, 1.0e-12,
"conservative estimate tracks newest slow step"));
1177 PetscCall(
PicurvAssertRealNear(24.0, required_headroom, 1.0e-12,
"required headroom scales with conservative estimate"));
1178 PetscFunctionReturn(0);
1186 PetscBool should_trigger = PETSC_FALSE;
1187 PetscReal required_headroom = 0.0;
1189 PetscFunctionBeginUser;
1191 PetscCall(
PicurvAssertBool((PetscBool)!should_trigger,
"guard should not trigger before warmup completes"));
1194 PetscCall(
PicurvAssertBool((PetscBool)!should_trigger,
"guard should not trigger when remaining walltime exceeds required headroom"));
1195 PetscCall(
PicurvAssertRealNear(28.0, required_headroom, 1.0e-12,
"required headroom after warmup"));
1198 PetscCall(
PicurvAssertBool(should_trigger,
"guard should trigger when remaining walltime reaches required headroom"));
1199 PetscCall(
PicurvAssertRealNear(28.0, required_headroom, 1.0e-12,
"required headroom remains unchanged at trigger threshold"));
1200 PetscFunctionReturn(0);
1208 PetscErrorCode ierr;
1246 ierr = PetscInitialize(&argc, &argv, NULL,
"PICurv runtime-kernel tests");
1251 ierr =
PicurvRunTests(
"unit-runtime", cases,
sizeof(cases) /
sizeof(cases[0]));
1257 ierr = PetscFinalize();
PetscErrorCode Validate_DrivenFlowConfiguration(UserCtx *user)
(Private) Validates all consistency rules for a driven flow (channel/pipe) setup.
Header file for Particle Motion and migration related functions.
PetscErrorCode CheckAndRemoveOutOfBoundsParticles(UserCtx *user, PetscInt *removedCountLocal, PetscInt *removedCountGlobal, const BoundingBox *bboxlist)
Checks for particles outside the physical domain boundaries and removes them using DMSwarmRemovePoint...
PetscErrorCode UpdateAllParticlePositions(UserCtx *user)
Loops over all local particles in the DMSwarm, updating their positions based on velocity and the glo...
PetscErrorCode CalculateParticleCountPerCell(UserCtx *user)
Counts particles in each cell of the DMDA 'da' and stores the result in user->ParticleCount.
PetscErrorCode CalculateBrownianDisplacement(UserCtx *user, PetscReal diff_eff, Cmpnts *displacement)
Calculates the stochastic displacement vector (Brownian motion) for a single particle.
PetscErrorCode LocateAllParticlesInGrid(UserCtx *user, BoundingBox *bboxlist)
Orchestrates the complete particle location and migration process for one timestep.
PetscErrorCode UpdateParticlePosition(UserCtx *user, Particle *particle)
Updates a particle's position based on its velocity and the timestep dt (stored in user->dt).
PetscErrorCode ResetAllParticleStatuses(UserCtx *user)
Marks all local particles as NEEDS_LOCATION for the next settlement pass.
PetscErrorCode CheckAndRemoveLostParticles(UserCtx *user, PetscInt *removedCountLocal, PetscInt *removedCountGlobal)
Removes particles that have been definitively flagged as LOST by the location algorithm.
Header file for Particle related physics modules.
PetscErrorCode UpdateParticleField(const char *fieldName, PetscReal dt, PetscReal *psi_io, PetscReal diffusivity, PetscReal mean_val, PetscReal cell_vol, PetscReal C_model)
Updates a single particle's field based on its state and physics model.
Header file for Particle Swarm management functions.
PetscErrorCode UpdateParticleWeights(PetscReal *d, Particle *particle)
Updates a particle's interpolation weights based on distances to cell faces.
PetscErrorCode DistributeParticles(PetscInt numParticles, PetscMPIInt rank, PetscMPIInt size, PetscInt *particlesPerProcess, PetscInt *remainder)
Distributes particles evenly across MPI processes, handling any remainders.
PetscBool IsParticleInsideBoundingBox(const BoundingBox *bbox, const Particle *particle)
Checks if a particle's location is within a specified bounding box.
PetscErrorCode ScatterAllParticleFieldsToEulerFields(UserCtx *user)
Scatters a predefined set of particle fields to their corresponding Eulerian fields.
PetscErrorCode SetInitialInteriorField(UserCtx *user, const char *fieldName)
Sets the initial values for the INTERIOR of a specified Eulerian field.
PetscErrorCode InterpolateAllFieldsToSwarm(UserCtx *user)
Interpolates all relevant fields from the DMDA to the DMSwarm.
PetscErrorCode ComputeSmagorinskyConstant(UserCtx *user)
Computes the dynamic Smagorinsky constant (Cs) for the LES model.
PetscReal RuntimeWalltimeGuardUpdateEWMA(PetscBool has_previous, PetscReal previous_ewma_seconds, PetscReal latest_step_seconds, PetscReal alpha)
Update an EWMA estimate for timestep wall-clock duration.
PetscReal RuntimeWalltimeGuardConservativeEstimate(PetscReal warmup_average_seconds, PetscReal ewma_seconds, PetscReal latest_step_seconds)
Return the conservative timestep estimate used by the walltime guard.
PetscErrorCode UpdateSolverHistoryVectors(UserCtx *user)
Copies the current time step's solution fields into history vectors (e.g., U(t_n) -> U_o,...
PetscReal RuntimeWalltimeGuardRequiredHeadroom(PetscReal min_seconds, PetscReal multiplier, PetscReal conservative_estimate_seconds)
Compute the required shutdown headroom from timestep estimate and floor.
PetscBool RuntimeWalltimeGuardShouldTrigger(PetscInt completed_steps, PetscInt warmup_steps, PetscReal remaining_seconds, PetscReal min_seconds, PetscReal multiplier, PetscReal warmup_average_seconds, PetscReal ewma_seconds, PetscReal latest_step_seconds, PetscReal *required_headroom_seconds_out)
Decide whether the runtime walltime guard should stop before another step.
PetscErrorCode GetOwnedCellRange(const DMDALocalInfo *info_nodes, PetscInt dim, PetscInt *xs_cell_global_out, PetscInt *xm_cell_local_out)
Determines the global starting index and number of CELLS owned by the current processor in a specifie...
PetscErrorCode ComputeAndStoreNeighborRanks(UserCtx *user)
Computes and stores the Cartesian neighbor ranks for the DMDA decomposition.
PetscBool RuntimeWalltimeGuardParsePositiveSeconds(const char *text, PetscReal *seconds_out)
Parse a positive floating-point seconds value from runtime metadata.
static PetscErrorCode TestUpdateParticleFieldIEMRelaxation(void)
Tests IEM relaxation updates for particle-carried fields.
static PetscErrorCode TestMinimalFixtureMirrorsProductionDMLayout(void)
Tests that the shared minimal fixture mirrors the production DA contract.
static PetscErrorCode TestUpdateParticlePositionDiffusivityGradientOnly(void)
Tests particle position updates driven only by diffusivity-gradient drift.
static PetscErrorCode TestSetInitialInteriorFieldConstantProfileOnZInlet(void)
Tests constant-profile interior initialization on a Z-direction inlet.
static PetscErrorCode TestLocateAllParticlesInGridGuessPathResolvesLocalParticle(void)
Tests the guess-then-verify orchestrator path for a local particle with an unknown prior cell.
static PetscErrorCode TestWallModelVelocityHelpers(void)
Tests closed-form and iterative wall-model velocity helpers against inverse reconstructions.
static PetscErrorCode TestSetInitialInteriorFieldIgnoresNonUcontRequest(void)
Tests that non-Ucont requests do not modify interior field initialization.
static PetscErrorCode TestLocateParticleOrFindMigrationTargetCountsReSearch(void)
Verifies that later settlement passes increment re-search metrics.
static PetscErrorCode TestRuntimeWalltimeGuardParsesPositiveSeconds(void)
Tests parsing of positive runtime walltime metadata values.
int main(int argc, char **argv)
Runs the unit-runtime PETSc test binary.
static PetscErrorCode TestInterpolateAllFieldsToSwarmConstantFields(void)
Tests direct interpolation from Eulerian fields to one localized swarm particle.
static PetscErrorCode SyncRuntimeFieldGhosts(UserCtx *user)
Synchronizes the minimal runtime fixture's global fields into their persistent local ghosts.
static PetscErrorCode TestDistributeParticlesRemainderHandling(void)
Tests particle distribution remainder handling across ranks.
static PetscErrorCode TestResetAllParticleStatusesLeavesLostParticlesUntouched(void)
Tests localized particle-status reset behavior for restart of the location workflow.
static PetscErrorCode SeedSingleParticle(UserCtx *user, PetscInt ci, PetscInt cj, PetscInt ck, PetscReal x, PetscReal y, PetscReal z, PetscReal wx, PetscReal wy, PetscReal wz, PetscInt status_value)
Seeds one localized swarm particle with the cell, position, weight, and status data used by runtime t...
static PetscErrorCode TestInterpolateAllFieldsToSwarmCornerAveragedConstantFields(void)
Tests the corner-averaged (legacy) interpolation path on constant fields.
static PetscErrorCode TestGetOwnedCellRangeSingleRankAccounting(void)
Tests owned-cell range accounting on a single MPI rank.
static PetscErrorCode TestUpdateParticlePositionWithoutBrownianContribution(void)
Tests particle position updates without Brownian forcing.
static PetscErrorCode TestComputeAndStoreNeighborRanksSingleRank(void)
Tests neighbor-rank discovery on a single MPI rank.
static PetscErrorCode TestCalculateParticleCountPerCellCountsGlobalCellIDs(void)
Tests particle counting by geometric cell IDs using the production +1 storage shift.
static PetscErrorCode TestValidateDrivenFlowConfigurationNoDrivenHandlers(void)
Tests driven-flow validation when no driven handlers are present.
static PetscErrorCode TestWallNoSlipAndFreeSlipHelpers(void)
Tests no-slip and free-slip wall helper kernels.
static PetscErrorCode TestMinimalFixtureRegistersProductionSwarmFields(void)
Tests that the shared swarm fixture registers the production field set.
static PetscErrorCode TestWallFunctionVectorWrappers(void)
Tests the vector wall-function wrappers on a tangential reference flow.
static PetscErrorCode TestLocateAllParticlesInGridPriorCellFastPath(void)
Tests the location orchestrator fast path when a particle already carries a valid prior cell.
static PetscErrorCode TestComputeSmagorinskyConstantConstantModel(void)
Tests the constant Smagorinsky model helper path.
static PetscErrorCode TestCalculateBrownianDisplacementDeterministicSeed(void)
Tests Brownian displacement generation against a duplicated seeded RNG stream.
static PetscErrorCode TestScatterAllParticleFieldsToEulerFieldsAveragesPsi(void)
Tests particle-to-grid scattering using known cell occupancy and scalar values.
static PetscErrorCode TestUpdateSolverHistoryVectorsShiftsStates(void)
Tests solver history-vector shifting between time levels.
static PetscErrorCode TestUpdateParticleWeightsComputesExpectedRatios(void)
Tests particle weight updates against expected ratios.
static PetscErrorCode TestCheckAndRemoveOutOfBoundsParticlesRemovesEscapedParticle(void)
Tests direct removal of particles that leave every rank bounding box.
static PetscErrorCode TestCheckAndRemoveLostParticlesRemovesLostEntries(void)
Tests direct removal of particles already marked LOST by the location workflow.
static PetscErrorCode TestRuntimeWalltimeGuardEstimatorHelpers(void)
Tests walltime-guard estimator helper calculations.
static PetscErrorCode TestIsParticleInsideBoundingBoxBasicCases(void)
Tests basic particle-inside-bounding-box classification cases.
static PetscErrorCode TestRuntimeWalltimeGuardTriggerDecision(void)
Tests runtime walltime-guard shutdown trigger decisions.
static PetscErrorCode TestUpdateAllParticlePositionsMovesSwarmEntries(void)
Tests swarm-wide particle position updates using the same transport path as the runtime loop.
static PetscErrorCode TestWallModelScalarHelpers(void)
Tests wall-model scalar helper kernels.
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 PicurvDestroyRuntimeContext(SimCtx **simCtx_ptr)
Finalizes and frees a runtime context built by PicurvBuildTinyRuntimeContext.
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 PicurvBuildTinyRuntimeContext(const char *bcs_contents, PetscBool enable_particles, SimCtx **simCtx_out, UserCtx **user_out, char *tmpdir, size_t tmpdir_len)
Builds a tiny runtime context through the real setup path for behavior-level tests.
PetscErrorCode PicurvAssertVecConstant(Vec vec, PetscScalar expected, PetscReal tol, const char *context)
Asserts that a PETSc vector is spatially constant within tolerance.
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.
PetscInt64 searchLocatedCount
PetscInt64 searchLostCount
BCFace identifiedInletBCFace
ParticleLocationStatus
Defines the state of a particle with respect to its location and migration status during the iterativ...
PetscInt64 traversalStepsSum
PetscInt FieldInitialization
Cmpnts max_coords
Maximum x, y, z coordinates of the bounding box.
Cmpnts diffusivitygradient
PetscInt64 searchPopulation
PetscInt currentSettlementPass
Cmpnts min_coords
Minimum x, y, z coordinates of the bounding box.
PetscInt64 bboxGuessFallbackCount
InterpolationMethod interpolationMethod
RankCellInfo * RankCellInfoMap
PetscInt64 bboxGuessSuccessCount
PetscInt64 maxParticlePassDepth
Cmpnts InitialConstantContra
SearchMetricsState searchMetrics
PetscRandom BrownianMotionRNG
PetscInt64 searchAttempts
PetscInt64 maxTraversalFailCount
Defines a 3D axis-aligned bounding box.
A 3D point or vector with PetscScalar components.
Defines a particle's core properties for Lagrangian tracking.
The master context for the entire simulation.
User-defined context containing data specific to a single computational grid level.
Header file for particle location functions using the walking search algorithm.
PetscErrorCode LocateParticleOrFindMigrationTarget(UserCtx *user, Particle *particle, ParticleLocationStatus *status_out)
Locates a particle's host cell or identifies its migration target using a robust walk search.
void wall_function_loglaw(UserCtx *user, double roughness_height, double distance_reference, double distance_boundary, Cmpnts velocity_wall, Cmpnts velocity_reference, Cmpnts *velocity_boundary, PetscReal *friction_velocity, double normal_x, double normal_y, double normal_z)
Applies log-law wall function with roughness correction.
double find_utau_loglaw(double velocity, double wall_distance, double roughness_length)
Solves for friction velocity using simple log-law (explicit formula)
double u_Werner(double kinematic_viscosity, double wall_distance, double friction_velocity)
Computes velocity using Werner-Wengle wall function.
void wall_function(UserCtx *user, double distance_reference, double distance_boundary, Cmpnts velocity_wall, Cmpnts velocity_reference, Cmpnts *velocity_boundary, PetscReal *friction_velocity, double normal_x, double normal_y, double normal_z)
Applies standard wall function with Werner-Wengle model.
double find_utau_hydset(double kinematic_viscosity, double known_velocity, double wall_distance, double initial_guess, double roughness_height)
Solves for friction velocity using Newton-Raphson iteration.
double u_Cabot(double kinematic_viscosity, double wall_distance, double friction_velocity, double pressure_gradient_tangent, double wall_shear_stress)
Computes velocity using Cabot wall function.
double u_loglaw(double wall_distance, double friction_velocity, double roughness_length)
Computes velocity using simple log-law (smooth wall with roughness offset)
void wall_function_Cabot(UserCtx *user, double roughness_height, double distance_reference, double distance_boundary, Cmpnts velocity_wall, Cmpnts velocity_reference, Cmpnts *velocity_boundary, PetscReal *friction_velocity, double normal_x, double normal_y, double normal_z, double pressure_gradient_x, double pressure_gradient_y, double pressure_gradient_z, int iteration_count)
Applies Cabot non-equilibrium wall function with pressure gradients.
void noslip(UserCtx *user, double distance_reference, double distance_boundary, Cmpnts velocity_wall, Cmpnts velocity_reference, Cmpnts *velocity_boundary, double normal_x, double normal_y, double normal_z)
Applies no-slip wall boundary condition with linear interpolation.
void find_utau_Cabot(double kinematic_viscosity, double velocity, double wall_distance, double initial_guess, double pressure_gradient_tangent, double pressure_gradient_normal, double *friction_velocity, double *wall_shear_velocity, double *wall_shear_normal)
Solves for friction velocity using Cabot wall function.
double integrate_1(double kinematic_viscosity, double wall_distance, double friction_velocity, int integration_mode)
Integrates eddy viscosity profile from wall to distance y.
double f_hydset(double kinematic_viscosity, double known_velocity, double wall_distance, double friction_velocity_guess, double roughness_height)
Residual function for friction velocity equation (log-law with roughness)
double E_coeff(double friction_velocity, double roughness_height, double kinematic_viscosity)
Computes roughness-modified log-law coefficient E.
double nu_t(double yplus)
Computes turbulent eddy viscosity ratio (ν_t / ν)
double find_utau_Werner(double kinematic_viscosity, double velocity, double wall_distance, double initial_guess)
Solves for friction velocity using Werner-Wengle wall function.
void freeslip(UserCtx *user, double distance_reference, double distance_boundary, Cmpnts velocity_wall, Cmpnts velocity_reference, Cmpnts *velocity_boundary, double normal_x, double normal_y, double normal_z)
Applies free-slip wall boundary condition.