PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
Loading...
Searching...
No Matches
test_runtime_kernels.c
Go to the documentation of this file.
1/**
2 * @file test_runtime_kernels.c
3 * @brief C unit tests for runtime, particle, wall, and walltime-guard helpers.
4 */
5
6#include "test_support.h"
7
8#include "BC_Handlers.h"
9#include "ParticleMotion.h"
10#include "ParticlePhysics.h"
11#include "ParticleSwarm.h"
12#include "interpolation.h"
13#include "initialcondition.h"
14#include "les.h"
15#include "runloop.h"
16#include "setup.h"
17#include "wallfunction.h"
18#include "walkingsearch.h"
19
20/**
21 * @brief Synchronizes the minimal runtime fixture's global fields into their persistent local ghosts.
22 */
23static PetscErrorCode SyncRuntimeFieldGhosts(UserCtx *user)
24{
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));
30 PetscCall(DMGlobalToLocalBegin(user->da, user->Diffusivity, INSERT_VALUES, user->lDiffusivity));
31 PetscCall(DMGlobalToLocalEnd(user->da, user->Diffusivity, INSERT_VALUES, user->lDiffusivity));
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));
36 PetscCall(DMGlobalToLocalBegin(user->fda, user->DiffusivityGradient, INSERT_VALUES, user->lDiffusivityGradient));
37 PetscCall(DMGlobalToLocalEnd(user->fda, user->DiffusivityGradient, INSERT_VALUES, user->lDiffusivityGradient));
38 PetscFunctionReturn(0);
39}
40
41/**
42 * @brief Seeds one localized swarm particle with the cell, position, weight, and status data used by runtime tests.
43 */
44static PetscErrorCode SeedSingleParticle(UserCtx *user,
45 PetscInt ci,
46 PetscInt cj,
47 PetscInt ck,
48 PetscReal x,
49 PetscReal y,
50 PetscReal z,
51 PetscReal wx,
52 PetscReal wy,
53 PetscReal wz,
54 PetscInt status_value)
55{
56 PetscReal *positions = NULL;
57 PetscReal *weights = NULL;
58 PetscInt *cell_ids = NULL;
59 PetscInt *status = NULL;
60
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));
66
67 positions[0] = x;
68 positions[1] = y;
69 positions[2] = z;
70 weights[0] = wx;
71 weights[1] = wy;
72 weights[2] = wz;
73 cell_ids[0] = ci;
74 cell_ids[1] = cj;
75 cell_ids[2] = ck;
76 status[0] = status_value;
77
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);
83}
84/**
85 * @brief Tests particle distribution remainder handling across ranks.
86 */
87
89{
90 PetscInt particles_per_rank = -1;
91 PetscInt remainder = -1;
92
93 PetscFunctionBeginUser;
94 PetscCall(DistributeParticles(10, 0, 3, &particles_per_rank, &remainder));
95 PetscCall(PicurvAssertIntEqual(4, particles_per_rank, "rank 0 should receive one remainder particle"));
96 PetscCall(PicurvAssertIntEqual(1, remainder, "remainder should be reported correctly"));
97
98 PetscCall(DistributeParticles(10, 2, 3, &particles_per_rank, &remainder));
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);
102}
103/**
104 * @brief Tests basic particle-inside-bounding-box classification cases.
105 */
106
108{
109 BoundingBox bbox;
110 Particle particle;
111
112 PetscFunctionBeginUser;
113 PetscCall(PetscMemzero(&bbox, sizeof(bbox)));
114 PetscCall(PetscMemzero(&particle, sizeof(particle)));
115
116 bbox.min_coords.x = 0.0;
117 bbox.min_coords.y = 0.0;
118 bbox.min_coords.z = 0.0;
119 bbox.max_coords.x = 1.0;
120 bbox.max_coords.y = 2.0;
121 bbox.max_coords.z = 3.0;
122
123 particle.loc.x = 0.25;
124 particle.loc.y = 1.0;
125 particle.loc.z = 2.5;
126 PetscCall(PicurvAssertBool(IsParticleInsideBoundingBox(&bbox, &particle), "particle should be inside bounding box"));
127
128 particle.loc.x = 1.5;
129 PetscCall(PicurvAssertBool((PetscBool)!IsParticleInsideBoundingBox(&bbox, &particle), "particle should be outside bounding box"));
130 PetscFunctionReturn(0);
131}
132/**
133 * @brief Tests particle weight updates against expected ratios.
134 */
135
137{
138 Particle particle;
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};
141
142 PetscFunctionBeginUser;
143 PetscCall(PetscMemzero(&particle, sizeof(particle)));
144 PetscCall(UpdateParticleWeights(distances, &particle));
145 PetscCall(PicurvAssertRealNear(0.25, particle.weights.x, 1.0e-12, "x interpolation weight"));
146 PetscCall(PicurvAssertRealNear(0.5, particle.weights.y, 1.0e-12, "y interpolation weight"));
147 PetscCall(PicurvAssertRealNear(0.2, particle.weights.z, 1.0e-12, "z interpolation weight"));
148
149 PetscCall(UpdateParticleWeights(clamped, &particle));
150 PetscCall(PicurvAssertRealNear(0.5, particle.weights.x, 1.0e-12, "clamped x weight remains centered"));
151 PetscCall(PicurvAssertRealNear(0.5, particle.weights.y, 1.0e-12, "clamped y weight remains centered"));
152 PetscCall(PicurvAssertRealNear(0.5, particle.weights.z, 1.0e-12, "clamped z weight remains centered"));
153 PetscFunctionReturn(0);
154}
155/**
156 * @brief Tests particle position updates without Brownian forcing.
157 */
158
160{
161 SimCtx *simCtx = NULL;
162 UserCtx *user = NULL;
163 Particle particle;
164
165 PetscFunctionBeginUser;
166 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
167 simCtx->dt = 0.25;
168
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;
176 particle.diffusivitygradient.x = 0.1;
177 particle.diffusivitygradient.y = 0.2;
178 particle.diffusivitygradient.z = -0.3;
179 particle.diffusivity = 0.0;
180
181 PetscCall(UpdateParticlePosition(user, &particle));
182 PetscCall(PicurvAssertRealNear(1.15, particle.loc.x, 1.0e-12, "x position update"));
183 PetscCall(PicurvAssertRealNear(-2.2, particle.loc.y, 1.0e-12, "y position update"));
184 PetscCall(PicurvAssertRealNear(3.425, particle.loc.z, 1.0e-12, "z position update"));
185
186 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
187 PetscFunctionReturn(0);
188}
189/**
190 * @brief Tests particle position updates driven only by diffusivity-gradient drift.
191 */
192
194{
195 SimCtx *simCtx = NULL;
196 UserCtx *user = NULL;
197 Particle particle;
198
199 PetscFunctionBeginUser;
200 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
201 simCtx->dt = 0.5;
202
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;
210 particle.diffusivitygradient.x = 0.2;
211 particle.diffusivitygradient.y = -0.1;
212 particle.diffusivitygradient.z = 0.3;
213 particle.diffusivity = 0.0;
214
215 PetscCall(UpdateParticlePosition(user, &particle));
216 PetscCall(PicurvAssertRealNear(0.35, particle.loc.x, 1.0e-12, "drift-only x position update"));
217 PetscCall(PicurvAssertRealNear(0.45, particle.loc.y, 1.0e-12, "drift-only y position update"));
218 PetscCall(PicurvAssertRealNear(0.90, particle.loc.z, 1.0e-12, "drift-only z position update"));
219
220 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
221 PetscFunctionReturn(0);
222}
223/**
224 * @brief Tests IEM relaxation updates for particle-carried fields.
225 */
226
227static PetscErrorCode TestUpdateParticleFieldIEMRelaxation(void)
228{
229 PetscReal psi = 1.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;
236
237 PetscFunctionBeginUser;
238 PetscCall(UpdateParticleField("Psi", dt, &psi, diffusivity, mean_val, cell_vol, c_model));
239 PetscCall(PicurvAssertRealNear(
240 mean_val + (1.0 - mean_val) * PetscExpReal(-(c_model * diffusivity / PetscPowReal(cell_vol, 0.6666667)) * dt),
241 psi,
242 1.0e-12,
243 "IEM update should match analytical relaxation"));
244
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);
248}
249/**
250 * @brief Tests that non-Ucont requests do not modify interior field initialization.
251 */
252
254{
255 SimCtx *simCtx = NULL;
256 UserCtx *user = NULL;
257
258 PetscFunctionBeginUser;
259 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
260 PetscCall(VecSet(user->Ucont, 7.0));
261
262 PetscCall(SetInitialInteriorField(user, "P"));
263 PetscCall(PicurvAssertVecConstant(user->Ucont, 7.0, 1.0e-12, "non-Ucont request should not modify Ucont"));
264
265 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
266 PetscFunctionReturn(0);
267}
268/**
269 * @brief Tests constant-profile interior initialization on a Z-direction inlet.
270 */
271
273{
274 SimCtx *simCtx = NULL;
275 UserCtx *user = NULL;
276 Cmpnts ***ucont = NULL;
277
278 PetscFunctionBeginUser;
279 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 5, 5, 5));
280 user->GridOrientation = 1;
282 simCtx->FieldInitialization = 1;
283 simCtx->InitialConstantContra.x = 0.0;
284 simCtx->InitialConstantContra.y = 0.0;
285 simCtx->InitialConstantContra.z = 2.0;
286 PetscCall(VecSet(user->Ucont, 0.0));
287
288 PetscCall(SetInitialInteriorField(user, "Ucont"));
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));
295
296 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
297 PetscFunctionReturn(0);
298}
299/**
300 * @brief Tests direct interpolation from Eulerian fields to one localized swarm particle.
301 */
302
304{
305 SimCtx *simCtx = NULL;
306 UserCtx *user = NULL;
307 Cmpnts ***grad = NULL;
308 PetscReal *velocity = NULL;
309 PetscReal *diffusivity = NULL;
310 PetscReal *diffusivity_gradient = NULL;
311
312 PetscFunctionBeginUser;
313 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
314 PetscCall(PicurvCreateSwarmPair(user, 1, "ske"));
315 PetscCall(VecSet(user->Ucat, 2.0));
316 PetscCall(VecSet(user->Diffusivity, 0.25));
317
318 PetscCall(DMDAVecGetArray(user->fda, user->DiffusivityGradient, &grad));
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;
325 }
326 }
327 }
328 PetscCall(DMDAVecRestoreArray(user->fda, user->DiffusivityGradient, &grad));
329 PetscCall(SyncRuntimeFieldGhosts(user));
330 PetscCall(SeedSingleParticle(user, 0, 0, 0, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, ACTIVE_AND_LOCATED));
331
332 PetscCall(InterpolateAllFieldsToSwarm(user));
333
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));
347
348 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
349 PetscFunctionReturn(0);
350}
351/**
352 * @brief Tests the corner-averaged (legacy) interpolation path on constant fields.
353 */
354
356{
357 SimCtx *simCtx = NULL;
358 UserCtx *user = NULL;
359 Cmpnts ***grad = NULL;
360 PetscReal *velocity = NULL;
361 PetscReal *diffusivity = NULL;
362 PetscReal *diffusivity_gradient = NULL;
363
364 PetscFunctionBeginUser;
365 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
367 PetscCall(PicurvCreateSwarmPair(user, 1, "ske"));
368 PetscCall(VecSet(user->Ucat, 2.0));
369 PetscCall(VecSet(user->Diffusivity, 0.25));
370
371 PetscCall(DMDAVecGetArray(user->fda, user->DiffusivityGradient, &grad));
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;
378 }
379 }
380 }
381 PetscCall(DMDAVecRestoreArray(user->fda, user->DiffusivityGradient, &grad));
382 PetscCall(SyncRuntimeFieldGhosts(user));
383 PetscCall(SeedSingleParticle(user, 0, 0, 0, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, ACTIVE_AND_LOCATED));
384
385 PetscCall(InterpolateAllFieldsToSwarm(user));
386
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));
400
401 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
402 PetscFunctionReturn(0);
403}
404/**
405 * @brief Tests particle-to-grid scattering using known cell occupancy and scalar values.
406 */
407
409{
410 SimCtx *simCtx = NULL;
411 UserCtx *user = NULL;
412 PetscInt *cell_ids = NULL;
413 PetscReal *psi = NULL;
414 PetscReal ***psi_grid = NULL;
415
416 PetscFunctionBeginUser;
417 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
418 PetscCall(PicurvCreateSwarmPair(user, 2, "ske"));
419
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));
424
425 PetscCall(DMSwarmGetField(user->swarm, "Psi", NULL, NULL, (void **)&psi));
426 psi[0] = 1.0;
427 psi[1] = 3.0;
428 PetscCall(DMSwarmRestoreField(user->swarm, "Psi", NULL, NULL, (void **)&psi));
429
431
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));
435
436 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
437 PetscFunctionReturn(0);
438}
439/**
440 * @brief Tests particle counting by geometric cell IDs using the production +1 storage shift.
441 */
442
444{
445 SimCtx *simCtx = NULL;
446 UserCtx *user = NULL;
447 PetscInt *cell_ids = NULL;
448 PetscReal ***counts = NULL;
449
450 PetscFunctionBeginUser;
451 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
452 PetscCall(PicurvCreateSwarmPair(user, 3, "ske"));
453
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));
459
460 PetscCall(CalculateParticleCountPerCell(user));
461
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));
466
467 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
468 PetscFunctionReturn(0);
469}
470/**
471 * @brief Tests localized particle-status reset behavior for restart of the location workflow.
472 */
473
475{
476 SimCtx *simCtx = NULL;
477 UserCtx *user = NULL;
478 PetscInt *status = NULL;
479
480 PetscFunctionBeginUser;
481 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
482 PetscCall(PicurvCreateSwarmPair(user, 3, "ske"));
483
484 PetscCall(DMSwarmGetField(user->swarm, "DMSwarm_location_status", NULL, NULL, (void **)&status));
485 status[0] = ACTIVE_AND_LOCATED;
486 status[1] = LOST;
487 status[2] = NEEDS_LOCATION;
488 PetscCall(DMSwarmRestoreField(user->swarm, "DMSwarm_location_status", NULL, NULL, (void **)&status));
489
490 PetscCall(ResetAllParticleStatuses(user));
491
492 PetscCall(DMSwarmGetField(user->swarm, "DMSwarm_location_status", NULL, NULL, (void **)&status));
493 PetscCall(PicurvAssertIntEqual(NEEDS_LOCATION, status[0], "ACTIVE_AND_LOCATED particles should be reset to NEEDS_LOCATION"));
494 PetscCall(PicurvAssertIntEqual(LOST, status[1], "LOST particles should remain LOST"));
495 PetscCall(PicurvAssertIntEqual(NEEDS_LOCATION, status[2], "NEEDS_LOCATION particles should remain unchanged"));
496 PetscCall(DMSwarmRestoreField(user->swarm, "DMSwarm_location_status", NULL, NULL, (void **)&status));
497
498 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
499 PetscFunctionReturn(0);
500}
501/**
502 * @brief Tests direct removal of particles that leave every rank bounding box.
503 */
504
506{
507 SimCtx *simCtx = NULL;
508 UserCtx *user = NULL;
509 PetscReal *positions = NULL;
510 PetscInt removed_local = 0;
511 PetscInt removed_global = 0;
512 PetscInt nlocal = 0;
513
514 PetscFunctionBeginUser;
515 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
516 PetscCall(PicurvCreateSwarmPair(user, 2, "ske"));
517
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));
522
523 PetscCall(CheckAndRemoveOutOfBoundsParticles(user, &removed_local, &removed_global, simCtx->bboxlist));
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"));
528
529 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
530 PetscFunctionReturn(0);
531}
532/**
533 * @brief Tests direct removal of particles already marked LOST by the location workflow.
534 */
535
537{
538 SimCtx *simCtx = NULL;
539 UserCtx *user = NULL;
540 PetscInt *status = NULL;
541 PetscInt removed_local = 0;
542 PetscInt removed_global = 0;
543 PetscInt nlocal = 0;
544
545 PetscFunctionBeginUser;
546 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
547 PetscCall(PicurvCreateSwarmPair(user, 3, "ske"));
548
549 PetscCall(DMSwarmGetField(user->swarm, "DMSwarm_location_status", NULL, NULL, (void **)&status));
550 status[0] = ACTIVE_AND_LOCATED;
551 status[1] = LOST;
552 status[2] = LOST;
553 PetscCall(DMSwarmRestoreField(user->swarm, "DMSwarm_location_status", NULL, NULL, (void **)&status));
554
555 PetscCall(CheckAndRemoveLostParticles(user, &removed_local, &removed_global));
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"));
559 PetscCall(PicurvAssertIntEqual(1, nlocal, "One non-LOST particle should remain"));
560
561 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
562 PetscFunctionReturn(0);
563}
564/**
565 * @brief Tests Brownian displacement generation against a duplicated seeded RNG stream.
566 */
567
569{
570 SimCtx *simCtx = NULL;
571 UserCtx *user = NULL;
572 char tmpdir[PETSC_MAX_PATH_LEN];
573 Cmpnts first;
574 Cmpnts second;
575
576 PetscFunctionBeginUser;
577 PetscCall(PicurvBuildTinyRuntimeContext(NULL, PETSC_FALSE, &simCtx, &user, tmpdir, sizeof(tmpdir)));
578 simCtx->dt = 0.25;
579 PetscCall(PicurvAssertBool((PetscBool)(simCtx->BrownianMotionRNG != NULL),
580 "runtime setup path should initialize the Brownian RNG"));
581
582 PetscCall(PetscRandomSetSeed(simCtx->BrownianMotionRNG, 12345));
583 PetscCall(PetscRandomSeed(simCtx->BrownianMotionRNG));
584
585 PetscCall(CalculateBrownianDisplacement(user, 0.5, &first));
586 PetscCall(PetscRandomSetSeed(simCtx->BrownianMotionRNG, 12345));
587 PetscCall(PetscRandomSeed(simCtx->BrownianMotionRNG));
588 PetscCall(CalculateBrownianDisplacement(user, 0.5, &second));
589
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"));
593
594 PetscCall(PicurvDestroyRuntimeContext(&simCtx));
595 PetscCall(PicurvRemoveTempDir(tmpdir));
596 PetscFunctionReturn(0);
597}
598/**
599 * @brief Tests swarm-wide particle position updates using the same transport path as the runtime loop.
600 */
602{
603 SimCtx *simCtx = NULL;
604 UserCtx *user = NULL;
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;
614
615 PetscFunctionBeginUser;
616 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
617 PetscCall(PicurvCreateSwarmPair(user, 1, "ske"));
618 simCtx->dt = 0.25;
619
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));
629
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;
636 psi[0] = 0.5;
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;
639 status[0] = ACTIVE_AND_LOCATED;
640 pid[0] = 7;
641
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));
651
652 PetscCall(UpdateAllParticlePositions(user));
653
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));
659
660 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
661 PetscFunctionReturn(0);
662}
663/**
664 * @brief Tests the location orchestrator fast path when a particle already carries a valid prior cell.
665 */
667{
668 SimCtx *simCtx = NULL;
669 UserCtx *user = NULL;
670 PetscReal *positions = NULL;
671 PetscReal *weights = NULL;
672 PetscInt *cell_ids = NULL;
673 PetscInt *status = NULL;
674 PetscInt64 *pid = NULL;
675
676 PetscFunctionBeginUser;
677 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
678 PetscCall(PetscMalloc1(simCtx->size, &user->RankCellInfoMap));
679 PetscCall(GetOwnedCellRange(&user->info, 0, &user->RankCellInfoMap[0].xs_cell, &user->RankCellInfoMap[0].xm_cell));
680 PetscCall(GetOwnedCellRange(&user->info, 1, &user->RankCellInfoMap[0].ys_cell, &user->RankCellInfoMap[0].ym_cell));
681 PetscCall(GetOwnedCellRange(&user->info, 2, &user->RankCellInfoMap[0].zs_cell, &user->RankCellInfoMap[0].zm_cell));
682 PetscCall(PicurvCreateSwarmPair(user, 1, "ske"));
683
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;
692 status[0] = NEEDS_LOCATION;
693 pid[0] = 11;
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));
699
700 PetscCall(LocateAllParticlesInGrid(user, simCtx->bboxlist));
701
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"));
707 PetscCall(PicurvAssertIntEqual(ACTIVE_AND_LOCATED, status[0], "prior-cell fast path should mark the particle ACTIVE_AND_LOCATED"));
708 PetscCall(PicurvAssertIntEqual(1, simCtx->searchMetrics.searchAttempts, "prior-cell fast path should record one search attempt"));
709 PetscCall(PicurvAssertIntEqual(1, (PetscInt)simCtx->searchMetrics.searchPopulation, "prior-cell fast path should record one input particle"));
710 PetscCall(PicurvAssertIntEqual(1, (PetscInt)simCtx->searchMetrics.searchLocatedCount, "prior-cell fast path should count one located particle"));
711 PetscCall(PicurvAssertIntEqual(0, (PetscInt)simCtx->searchMetrics.searchLostCount, "prior-cell fast path should not lose the particle"));
712 PetscCall(PicurvAssertIntEqual(0, (PetscInt)simCtx->searchMetrics.reSearchCount, "prior-cell fast path should not re-search on later passes"));
713 PetscCall(PicurvAssertBool((PetscBool)(simCtx->searchMetrics.traversalStepsSum > 0), "prior-cell fast path should accumulate traversal steps"));
714 PetscCall(PicurvAssertIntEqual(1, simCtx->searchMetrics.maxParticlePassDepth, "prior-cell fast path should report one settlement pass"));
715 PetscCall(DMSwarmRestoreField(user->swarm, "DMSwarm_location_status", NULL, NULL, (void **)&status));
716 PetscCall(DMSwarmRestoreField(user->swarm, "DMSwarm_CellID", NULL, NULL, (void **)&cell_ids));
717
718 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
719 PetscFunctionReturn(0);
720}
721/**
722 * @brief Tests the guess-then-verify orchestrator path for a local particle with an unknown prior cell.
723 */
725{
726 SimCtx *simCtx = NULL;
727 UserCtx *user = NULL;
728 PetscReal *positions = NULL;
729 PetscReal *weights = NULL;
730 PetscInt *cell_ids = NULL;
731 PetscInt *status = NULL;
732 PetscInt64 *pid = NULL;
733
734 PetscFunctionBeginUser;
735 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
736 PetscCall(PetscMalloc1(simCtx->size, &user->RankCellInfoMap));
737 PetscCall(GetOwnedCellRange(&user->info, 0, &user->RankCellInfoMap[0].xs_cell, &user->RankCellInfoMap[0].xm_cell));
738 PetscCall(GetOwnedCellRange(&user->info, 1, &user->RankCellInfoMap[0].ys_cell, &user->RankCellInfoMap[0].ym_cell));
739 PetscCall(GetOwnedCellRange(&user->info, 2, &user->RankCellInfoMap[0].zs_cell, &user->RankCellInfoMap[0].zm_cell));
740 PetscCall(PicurvCreateSwarmPair(user, 1, "ske"));
741
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;
750 status[0] = NEEDS_LOCATION;
751 pid[0] = 22;
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));
757
758 PetscCall(LocateAllParticlesInGrid(user, simCtx->bboxlist));
759
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"));
765 PetscCall(PicurvAssertIntEqual(ACTIVE_AND_LOCATED, status[0], "guess-path location should mark the particle ACTIVE_AND_LOCATED"));
766 PetscCall(PicurvAssertIntEqual(1, simCtx->searchMetrics.searchAttempts, "guess-path location should perform one robust search"));
767 PetscCall(PicurvAssertIntEqual(1, (PetscInt)simCtx->searchMetrics.searchPopulation, "guess-path location should count one input particle"));
768 PetscCall(PicurvAssertIntEqual(1, (PetscInt)simCtx->searchMetrics.searchLocatedCount, "guess-path location should count one located particle"));
769 PetscCall(PicurvAssertIntEqual(0, (PetscInt)simCtx->searchMetrics.searchLostCount, "guess-path location should not lose the particle"));
770 PetscCall(PicurvAssertIntEqual(0, (PetscInt)simCtx->searchMetrics.reSearchCount, "guess-path location should not count later-pass re-searches"));
771 PetscCall(PicurvAssertIntEqual(1, simCtx->searchMetrics.bboxGuessFallbackCount, "guess-path location should record one bbox fallback"));
772 PetscCall(PicurvAssertIntEqual(0, simCtx->searchMetrics.bboxGuessSuccessCount, "guess-path local resolution should not count as remote bbox success"));
773 PetscCall(DMSwarmRestoreField(user->swarm, "DMSwarm_location_status", NULL, NULL, (void **)&status));
774 PetscCall(DMSwarmRestoreField(user->swarm, "DMSwarm_CellID", NULL, NULL, (void **)&cell_ids));
775
776 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
777 PetscFunctionReturn(0);
778}
779/**
780 * @brief Verifies that later settlement passes increment re-search metrics.
781 */
783{
784 SimCtx *simCtx = NULL;
785 UserCtx *user = NULL;
786 Particle particle;
788
789 PetscFunctionBeginUser;
790 PetscCall(PetscMemzero(&particle, sizeof(particle)));
791 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
792 PetscCall(PetscMalloc1(simCtx->size, &user->RankCellInfoMap));
793 PetscCall(GetOwnedCellRange(&user->info, 0, &user->RankCellInfoMap[0].xs_cell, &user->RankCellInfoMap[0].xm_cell));
794 PetscCall(GetOwnedCellRange(&user->info, 1, &user->RankCellInfoMap[0].ys_cell, &user->RankCellInfoMap[0].ym_cell));
795 PetscCall(GetOwnedCellRange(&user->info, 2, &user->RankCellInfoMap[0].zs_cell, &user->RankCellInfoMap[0].zm_cell));
796
798 particle.PID = 33;
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;
805 particle.weights.x = 0.5;
806 particle.weights.y = 0.5;
807 particle.weights.z = 0.5;
808
809 PetscCall(LocateParticleOrFindMigrationTarget(user, &particle, &status));
810
811 PetscCall(PicurvAssertIntEqual(ACTIVE_AND_LOCATED, status, "direct re-search test should locate the particle"));
812 PetscCall(PicurvAssertIntEqual(1, (PetscInt)simCtx->searchMetrics.searchAttempts, "direct re-search test should record one robust walk"));
813 PetscCall(PicurvAssertIntEqual(1, (PetscInt)simCtx->searchMetrics.reSearchCount, "direct re-search test should increment re_search_count on later passes"));
814 PetscCall(PicurvAssertIntEqual(0, (PetscInt)simCtx->searchMetrics.maxTraversalFailCount, "direct re-search test should not hit MAX_TRAVERSAL"));
815
816 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
817 PetscFunctionReturn(0);
818}
819/**
820 * @brief Tests no-slip and free-slip wall helper kernels.
821 */
822
823static PetscErrorCode TestWallNoSlipAndFreeSlipHelpers(void)
824{
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};
829
830 PetscFunctionBeginUser;
831 noslip(NULL, 2.0, 1.0, wall_velocity, reference_velocity, &boundary_velocity, 1.0, 0.0, 0.0);
832 PetscCall(PicurvAssertRealNear(1.0, boundary_velocity.x, 1.0e-12, "no-slip interpolated x"));
833 PetscCall(PicurvAssertRealNear(2.0, boundary_velocity.y, 1.0e-12, "no-slip interpolated y"));
834 PetscCall(PicurvAssertRealNear(3.0, boundary_velocity.z, 1.0e-12, "no-slip interpolated z"));
835
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);
841}
842/**
843 * @brief Tests wall-model scalar helper kernels.
844 */
845
846static PetscErrorCode TestWallModelScalarHelpers(void)
847{
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;
852
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"));
856
857 utau = find_utau_hydset(1.0e-3, 1.0, 1.0e-2, 0.1, 0.0);
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"));
861
862 PetscCall(PicurvAssertRealNear(0.0, nu_t(0.0), 1.0e-12, "eddy viscosity ratio at wall"));
863 PetscCall(PicurvAssertBool((PetscBool)(integrate_1(1.0e-3, 1.0e-2, 0.1, 0) > 0.0), "integral helper should be positive"));
864 PetscFunctionReturn(0);
865}
866/**
867 * @brief Tests closed-form and iterative wall-model velocity helpers against inverse reconstructions.
868 */
869
870static PetscErrorCode TestWallModelVelocityHelpers(void)
871{
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;
881
882 PetscFunctionBeginUser;
883 utau_loglaw = find_utau_loglaw(target_velocity, wall_distance, roughness_length);
884 PetscCall(PicurvAssertRealNear(target_velocity, u_loglaw(wall_distance, utau_loglaw, roughness_length), 1.0e-12,
885 "simple log-law inversion should reconstruct the target velocity"));
886
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"));
889 PetscCall(PicurvAssertRealNear(target_velocity, u_Werner(kinematic_viscosity, wall_distance, utau_werner), 1.0e-6,
890 "Werner-Wengle inversion should reconstruct the target velocity"));
891
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"));
897 PetscCall(PicurvAssertRealNear(0.0, wall_shear_normal, 1.0e-10,
898 "zero normal pressure gradient should keep Cabot normal wall shear at zero"));
899 PetscFunctionReturn(0);
900}
901/**
902 * @brief Tests the vector wall-function wrappers on a tangential reference flow.
903 */
904static PetscErrorCode TestWallFunctionVectorWrappers(void)
905{
906 SimCtx *simCtx = NULL;
907 UserCtx *user = NULL;
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;
912
913 PetscFunctionBeginUser;
914 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
915 simCtx->ren = 1000.0;
916
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"));
921
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"));
926
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"));
932
933 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
934 PetscFunctionReturn(0);
935}
936/**
937 * @brief Tests driven-flow validation when no driven handlers are present.
938 */
939
941{
942 SimCtx *simCtx = NULL;
943 UserCtx *user = NULL;
944
945 PetscFunctionBeginUser;
946 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
947 PetscCall(Validate_DrivenFlowConfiguration(user));
948 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
949 PetscFunctionReturn(0);
950}
951/**
952 * @brief Tests the constant Smagorinsky model helper path.
953 */
954
956{
957 SimCtx *simCtx = NULL;
958 UserCtx *user = NULL;
959
960 PetscFunctionBeginUser;
961 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
962 PetscCall(DMCreateGlobalVector(user->da, &user->CS));
963 simCtx->step = 2;
964 simCtx->StartStep = 0;
965 simCtx->les = CONSTANT_SMAGORINSKY;
966 simCtx->Const_CS = 0.17;
967
968 PetscCall(ComputeSmagorinskyConstant(user));
969 PetscCall(PicurvAssertVecConstant(user->CS, 0.17, 1.0e-12, "constant Smagorinsky branch should fill CS"));
970
971 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
972 PetscFunctionReturn(0);
973}
974/**
975 * @brief Tests that the shared minimal fixture mirrors the production DA contract.
976 */
977
979{
980 SimCtx *simCtx = NULL;
981 UserCtx *user = NULL;
982 DM coord_dm = NULL;
983 PetscInt mx = 0, my = 0, mz = 0;
984
985 PetscFunctionBeginUser;
986 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 8, 6, 4));
987
988 PetscCall(DMDAGetInfo(user->da, NULL, &mx, &my, &mz, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
989 PetscCall(PicurvAssertIntEqual(user->IM + 1, mx, "minimal fixture should size da with IM+1 nodes"));
990 PetscCall(PicurvAssertIntEqual(user->JM + 1, my, "minimal fixture should size da with JM+1 nodes"));
991 PetscCall(PicurvAssertIntEqual(user->KM + 1, mz, "minimal fixture should size da with KM+1 nodes"));
992 PetscCall(PicurvAssertIntEqual(mx, user->info.mx, "user->info should be sourced from the production da"));
993 PetscCall(PicurvAssertIntEqual(my, user->info.my, "user->info my should match the da dimensions"));
994 PetscCall(PicurvAssertIntEqual(mz, user->info.mz, "user->info mz should match the da dimensions"));
995
996 PetscCall(DMGetCoordinateDM(user->da, &coord_dm));
997 PetscCall(PicurvAssertBool((PetscBool)(coord_dm == user->fda),
998 "minimal fixture should derive fda from the coordinate-DM path"));
999
1000 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
1001 PetscFunctionReturn(0);
1002}
1003/**
1004 * @brief Tests that the shared swarm fixture registers the production field set.
1005 */
1006
1008{
1009 SimCtx *simCtx = NULL;
1010 UserCtx *user = NULL;
1011 PetscInt bs = 0;
1012 void *field_ptr = NULL;
1013
1014 PetscFunctionBeginUser;
1015 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
1016 PetscCall(PicurvCreateSwarmPair(user, 2, "ske"));
1017
1018 PetscCall(DMSwarmGetField(user->swarm, "position", &bs, NULL, &field_ptr));
1019 PetscCall(PicurvAssertIntEqual(3, bs, "solver swarm should register particle position"));
1020 PetscCall(PicurvAssertBool((PetscBool)(field_ptr != NULL), "position field should be retrievable"));
1021 PetscCall(DMSwarmRestoreField(user->swarm, "position", &bs, NULL, &field_ptr));
1022
1023 PetscCall(DMSwarmGetField(user->swarm, "velocity", &bs, NULL, &field_ptr));
1024 PetscCall(PicurvAssertIntEqual(3, bs, "solver swarm should register particle velocity"));
1025 PetscCall(DMSwarmRestoreField(user->swarm, "velocity", &bs, NULL, &field_ptr));
1026
1027 PetscCall(DMSwarmGetField(user->swarm, "DMSwarm_CellID", &bs, NULL, &field_ptr));
1028 PetscCall(PicurvAssertIntEqual(3, bs, "solver swarm should register DMSwarm_CellID"));
1029 PetscCall(DMSwarmRestoreField(user->swarm, "DMSwarm_CellID", &bs, NULL, &field_ptr));
1030
1031 PetscCall(DMSwarmGetField(user->swarm, "weight", &bs, NULL, &field_ptr));
1032 PetscCall(PicurvAssertIntEqual(3, bs, "solver swarm should register particle weight"));
1033 PetscCall(DMSwarmRestoreField(user->swarm, "weight", &bs, NULL, &field_ptr));
1034
1035 PetscCall(DMSwarmGetField(user->swarm, "Diffusivity", &bs, NULL, &field_ptr));
1036 PetscCall(PicurvAssertIntEqual(1, bs, "solver swarm should register particle diffusivity"));
1037 PetscCall(DMSwarmRestoreField(user->swarm, "Diffusivity", &bs, NULL, &field_ptr));
1038
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));
1042
1043 PetscCall(DMSwarmGetField(user->swarm, "Psi", &bs, NULL, &field_ptr));
1044 PetscCall(PicurvAssertIntEqual(1, bs, "solver swarm should register particle scalar Psi"));
1045 PetscCall(DMSwarmRestoreField(user->swarm, "Psi", &bs, NULL, &field_ptr));
1046
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));
1050
1051 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
1052 PetscFunctionReturn(0);
1053}
1054/**
1055 * @brief Tests solver history-vector shifting between time levels.
1056 */
1057
1059{
1060 SimCtx *simCtx = NULL;
1061 UserCtx *user = NULL;
1062
1063 PetscFunctionBeginUser;
1064 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 5, 5, 5));
1065
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));
1073
1074 PetscCall(UpdateSolverHistoryVectors(user));
1075
1076 PetscCall(PicurvAssertVecConstant(user->Ucont_o, 11.0, 1.0e-12, "Ucont_o should receive current Ucont"));
1077 PetscCall(PicurvAssertVecConstant(user->Ucont_rm1, 7.0, 1.0e-12, "Ucont_rm1 should receive prior Ucont_o"));
1078 PetscCall(PicurvAssertVecConstant(user->Ucat_o, 5.0, 1.0e-12, "Ucat_o should receive current Ucat"));
1079 PetscCall(PicurvAssertVecConstant(user->P_o, 9.0, 1.0e-12, "P_o should receive current P"));
1080 PetscCall(PicurvAssertVecConstant(user->lUcont_o, 11.0, 1.0e-12, "lUcont_o ghost sync should match Ucont_o"));
1081 PetscCall(PicurvAssertVecConstant(user->lUcont_rm1, 7.0, 1.0e-12, "lUcont_rm1 ghost sync should match Ucont_rm1"));
1082
1083 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
1084 PetscFunctionReturn(0);
1085}
1086/**
1087 * @brief Tests owned-cell range accounting on a single MPI rank.
1088 */
1089
1091{
1092 SimCtx *simCtx = NULL;
1093 UserCtx *user = NULL;
1094 PetscInt xs = -1, ys = -1, zs = -1;
1095 PetscInt xm = -1, ym = -1, zm = -1;
1096
1097 PetscFunctionBeginUser;
1098 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 8, 6, 4));
1099
1100 PetscCall(GetOwnedCellRange(&user->info, 0, &xs, &xm));
1101 PetscCall(GetOwnedCellRange(&user->info, 1, &ys, &ym));
1102 PetscCall(GetOwnedCellRange(&user->info, 2, &zs, &zm));
1103
1104 PetscCall(PicurvAssertIntEqual(0, xs, "single-rank x cell-start index"));
1105 PetscCall(PicurvAssertIntEqual(0, ys, "single-rank y cell-start index"));
1106 PetscCall(PicurvAssertIntEqual(0, zs, "single-rank z cell-start index"));
1107 PetscCall(PicurvAssertIntEqual(user->info.mx - 2, xm, "single-rank x owned cell count"));
1108 PetscCall(PicurvAssertIntEqual(user->info.my - 2, ym, "single-rank y owned cell count"));
1109 PetscCall(PicurvAssertIntEqual(user->info.mz - 2, zm, "single-rank z owned cell count"));
1110
1111 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
1112 PetscFunctionReturn(0);
1113}
1114/**
1115 * @brief Tests neighbor-rank discovery on a single MPI rank.
1116 */
1117
1119{
1120 SimCtx *simCtx = NULL;
1121 UserCtx *user = NULL;
1122
1123 PetscFunctionBeginUser;
1124 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 6, 6, 6));
1125 PetscCall(ComputeAndStoreNeighborRanks(user));
1126
1127 PetscCall(PicurvAssertIntEqual(MPI_PROC_NULL, user->neighbors.rank_xm, "single-rank xm neighbor should be null"));
1128 PetscCall(PicurvAssertIntEqual(MPI_PROC_NULL, user->neighbors.rank_xp, "single-rank xp neighbor should be null"));
1129 PetscCall(PicurvAssertIntEqual(MPI_PROC_NULL, user->neighbors.rank_ym, "single-rank ym neighbor should be null"));
1130 PetscCall(PicurvAssertIntEqual(MPI_PROC_NULL, user->neighbors.rank_yp, "single-rank yp neighbor should be null"));
1131 PetscCall(PicurvAssertIntEqual(MPI_PROC_NULL, user->neighbors.rank_zm, "single-rank zm neighbor should be null"));
1132 PetscCall(PicurvAssertIntEqual(MPI_PROC_NULL, user->neighbors.rank_zp, "single-rank zp neighbor should be null"));
1133
1134 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
1135 PetscFunctionReturn(0);
1136}
1137/**
1138 * @brief Tests parsing of positive runtime walltime metadata values.
1139 */
1140
1142{
1143 PetscReal seconds = 0.0;
1144
1145 PetscFunctionBeginUser;
1146 PetscCall(PicurvAssertBool(RuntimeWalltimeGuardParsePositiveSeconds("300", &seconds), "plain positive seconds should parse"));
1147 PetscCall(PicurvAssertRealNear(300.0, seconds, 1.0e-12, "parsed walltime seconds"));
1148 PetscCall(PicurvAssertBool(RuntimeWalltimeGuardParsePositiveSeconds(" 42.5 ", &seconds), "whitespace-wrapped decimal seconds should parse"));
1149 PetscCall(PicurvAssertRealNear(42.5, seconds, 1.0e-12, "parsed decimal walltime seconds"));
1150 PetscCall(PicurvAssertBool((PetscBool)!RuntimeWalltimeGuardParsePositiveSeconds("nope", &seconds), "non-numeric metadata should fail parsing"));
1151 PetscCall(PicurvAssertBool((PetscBool)!RuntimeWalltimeGuardParsePositiveSeconds("-10", &seconds), "negative metadata should fail parsing"));
1152 PetscFunctionReturn(0);
1153}
1154/**
1155 * @brief Tests walltime-guard estimator helper calculations.
1156 */
1157
1159{
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;
1165
1166 PetscFunctionBeginUser;
1167 ewma_fast = RuntimeWalltimeGuardUpdateEWMA(PETSC_TRUE, 4.0, 6.0, 0.5);
1168 ewma_slow = RuntimeWalltimeGuardUpdateEWMA(PETSC_TRUE, ewma_fast, 12.0, 0.5);
1169 conservative_fast = RuntimeWalltimeGuardConservativeEstimate(5.0, ewma_fast, 6.0);
1170 conservative_slow = RuntimeWalltimeGuardConservativeEstimate(5.0, ewma_slow, 12.0);
1171 required_headroom = RuntimeWalltimeGuardRequiredHeadroom(8.0, 2.0, conservative_slow);
1172
1173 PetscCall(PicurvAssertRealNear(5.0, ewma_fast, 1.0e-12, "EWMA after moderate step"));
1174 PetscCall(PicurvAssertRealNear(8.5, ewma_slow, 1.0e-12, "EWMA after newer slow step"));
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);
1179}
1180/**
1181 * @brief Tests runtime walltime-guard shutdown trigger decisions.
1182 */
1183
1185{
1186 PetscBool should_trigger = PETSC_FALSE;
1187 PetscReal required_headroom = 0.0;
1188
1189 PetscFunctionBeginUser;
1190 should_trigger = RuntimeWalltimeGuardShouldTrigger(9, 10, 15.0, 5.0, 2.0, 6.0, 6.0, 6.0, &required_headroom);
1191 PetscCall(PicurvAssertBool((PetscBool)!should_trigger, "guard should not trigger before warmup completes"));
1192
1193 should_trigger = RuntimeWalltimeGuardShouldTrigger(10, 10, 40.0, 5.0, 2.0, 10.0, 12.0, 14.0, &required_headroom);
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"));
1196
1197 should_trigger = RuntimeWalltimeGuardShouldTrigger(10, 10, 28.0, 5.0, 2.0, 10.0, 12.0, 14.0, &required_headroom);
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);
1201}
1202/**
1203 * @brief Runs the unit-runtime PETSc test binary.
1204 */
1205
1206int main(int argc, char **argv)
1207{
1208 PetscErrorCode ierr;
1209 const PicurvTestCase cases[] = {
1210 {"distribute-particles-remainder-handling", TestDistributeParticlesRemainderHandling},
1211 {"is-particle-inside-bbox-basic-cases", TestIsParticleInsideBoundingBoxBasicCases},
1212 {"update-particle-weights-computes-expected-ratios", TestUpdateParticleWeightsComputesExpectedRatios},
1213 {"update-particle-position-without-brownian-contribution", TestUpdateParticlePositionWithoutBrownianContribution},
1214 {"update-particle-position-diffusivity-gradient-only", TestUpdateParticlePositionDiffusivityGradientOnly},
1215 {"update-particle-field-iem-relaxation", TestUpdateParticleFieldIEMRelaxation},
1216 {"set-initial-interior-field-ignores-non-ucont-request", TestSetInitialInteriorFieldIgnoresNonUcontRequest},
1217 {"set-initial-interior-field-constant-profile-on-z-inlet", TestSetInitialInteriorFieldConstantProfileOnZInlet},
1218 {"interpolate-all-fields-to-swarm-constant-fields", TestInterpolateAllFieldsToSwarmConstantFields},
1219 {"interpolate-all-fields-to-swarm-corner-averaged-constant-fields", TestInterpolateAllFieldsToSwarmCornerAveragedConstantFields},
1220 {"scatter-all-particle-fields-to-euler-fields-averages-psi", TestScatterAllParticleFieldsToEulerFieldsAveragesPsi},
1221 {"calculate-particle-count-per-cell-counts-global-cell-ids", TestCalculateParticleCountPerCellCountsGlobalCellIDs},
1222 {"reset-all-particle-statuses-leaves-lost-particles-untouched", TestResetAllParticleStatusesLeavesLostParticlesUntouched},
1223 {"check-and-remove-out-of-bounds-particles-removes-escaped-particle", TestCheckAndRemoveOutOfBoundsParticlesRemovesEscapedParticle},
1224 {"check-and-remove-lost-particles-removes-lost-entries", TestCheckAndRemoveLostParticlesRemovesLostEntries},
1225 {"calculate-brownian-displacement-deterministic-seed", TestCalculateBrownianDisplacementDeterministicSeed},
1226 {"update-all-particle-positions-moves-swarm-entries", TestUpdateAllParticlePositionsMovesSwarmEntries},
1227 {"locate-all-particles-in-grid-prior-cell-fast-path", TestLocateAllParticlesInGridPriorCellFastPath},
1228 {"locate-all-particles-in-grid-guess-path-resolves-local-particle", TestLocateAllParticlesInGridGuessPathResolvesLocalParticle},
1229 {"locate-particle-or-find-migration-target-counts-research", TestLocateParticleOrFindMigrationTargetCountsReSearch},
1230 {"wall-noslip-and-freeslip-helpers", TestWallNoSlipAndFreeSlipHelpers},
1231 {"wall-model-scalar-helpers", TestWallModelScalarHelpers},
1232 {"wall-model-velocity-helpers", TestWallModelVelocityHelpers},
1233 {"wall-function-vector-wrappers", TestWallFunctionVectorWrappers},
1234 {"validate-driven-flow-configuration-no-driven-handlers", TestValidateDrivenFlowConfigurationNoDrivenHandlers},
1235 {"compute-smagorinsky-constant-constant-model", TestComputeSmagorinskyConstantConstantModel},
1236 {"minimal-fixture-mirrors-production-dm-layout", TestMinimalFixtureMirrorsProductionDMLayout},
1237 {"minimal-fixture-registers-production-swarm-fields", TestMinimalFixtureRegistersProductionSwarmFields},
1238 {"update-solver-history-vectors-shifts-states", TestUpdateSolverHistoryVectorsShiftsStates},
1239 {"get-owned-cell-range-single-rank-accounting", TestGetOwnedCellRangeSingleRankAccounting},
1240 {"compute-and-store-neighbor-ranks-single-rank", TestComputeAndStoreNeighborRanksSingleRank},
1241 {"runtime-walltime-guard-parses-positive-seconds", TestRuntimeWalltimeGuardParsesPositiveSeconds},
1242 {"runtime-walltime-guard-estimator-helpers", TestRuntimeWalltimeGuardEstimatorHelpers},
1243 {"runtime-walltime-guard-trigger-decision", TestRuntimeWalltimeGuardTriggerDecision},
1244 };
1245
1246 ierr = PetscInitialize(&argc, &argv, NULL, "PICurv runtime-kernel tests");
1247 if (ierr) {
1248 return (int)ierr;
1249 }
1250
1251 ierr = PicurvRunTests("unit-runtime", cases, sizeof(cases) / sizeof(cases[0]));
1252 if (ierr) {
1253 PetscFinalize();
1254 return (int)ierr;
1255 }
1256
1257 ierr = PetscFinalize();
1258 return (int)ierr;
1259}
PetscErrorCode Validate_DrivenFlowConfiguration(UserCtx *user)
(Private) Validates all consistency rules for a driven flow (channel/pipe) setup.
Definition BC_Handlers.c:15
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.
Definition les.c:30
PetscReal RuntimeWalltimeGuardUpdateEWMA(PetscBool has_previous, PetscReal previous_ewma_seconds, PetscReal latest_step_seconds, PetscReal alpha)
Update an EWMA estimate for timestep wall-clock duration.
Definition runloop.c:150
PetscReal RuntimeWalltimeGuardConservativeEstimate(PetscReal warmup_average_seconds, PetscReal ewma_seconds, PetscReal latest_step_seconds)
Return the conservative timestep estimate used by the walltime guard.
Definition runloop.c:162
PetscErrorCode UpdateSolverHistoryVectors(UserCtx *user)
Copies the current time step's solution fields into history vectors (e.g., U(t_n) -> U_o,...
Definition runloop.c:319
PetscReal RuntimeWalltimeGuardRequiredHeadroom(PetscReal min_seconds, PetscReal multiplier, PetscReal conservative_estimate_seconds)
Compute the required shutdown headroom from timestep estimate and floor.
Definition runloop.c:173
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.
Definition runloop.c:184
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...
Definition setup.c:1883
PetscErrorCode ComputeAndStoreNeighborRanks(UserCtx *user)
Computes and stores the Cartesian neighbor ranks for the DMDA decomposition.
Definition setup.c:1980
PetscBool RuntimeWalltimeGuardParsePositiveSeconds(const char *text, PetscReal *seconds_out)
Parse a positive floating-point seconds value from runtime metadata.
Definition setup.c:18
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.
PetscMPIInt rank_zm
Definition variables.h:182
@ CONSTANT_SMAGORINSKY
Definition variables.h:490
Vec lDiffusivityGradient
Definition variables.h:841
Cmpnts vel
Definition variables.h:169
PetscInt ys_cell
Definition variables.h:187
PetscInt xs_cell
Definition variables.h:187
PetscMPIInt rank_yp
Definition variables.h:181
PetscInt64 searchLocatedCount
Definition variables.h:224
PetscInt64 searchLostCount
Definition variables.h:225
BCFace identifiedInletBCFace
Definition variables.h:831
PetscInt cell[3]
Definition variables.h:167
ParticleLocationStatus
Defines the state of a particle with respect to its location and migration status during the iterativ...
Definition variables.h:135
@ LOST
Definition variables.h:139
@ NEEDS_LOCATION
Definition variables.h:136
@ ACTIVE_AND_LOCATED
Definition variables.h:137
PetscMPIInt rank_ym
Definition variables.h:181
PetscMPIInt rank_xp
Definition variables.h:180
PetscInt KM
Definition variables.h:820
PetscInt64 traversalStepsSum
Definition variables.h:226
PetscInt FieldInitialization
Definition variables.h:696
PetscReal ren
Definition variables.h:691
Vec lUcont_rm1
Definition variables.h:845
PetscInt zm_cell
Definition variables.h:188
Cmpnts max_coords
Maximum x, y, z coordinates of the bounding box.
Definition variables.h:156
PetscInt zs_cell
Definition variables.h:187
Cmpnts diffusivitygradient
Definition variables.h:174
PetscInt64 searchPopulation
Definition variables.h:223
PetscReal dt
Definition variables.h:658
RankNeighbors neighbors
Definition variables.h:823
Vec lPsi
Definition variables.h:883
PetscInt currentSettlementPass
Definition variables.h:235
Vec DiffusivityGradient
Definition variables.h:841
Vec Ucont
Definition variables.h:837
PetscInt StartStep
Definition variables.h:653
Cmpnts min_coords
Minimum x, y, z coordinates of the bounding box.
Definition variables.h:155
PetscScalar x
Definition variables.h:101
Cmpnts loc
Definition variables.h:168
PetscInt64 reSearchCount
Definition variables.h:227
PetscInt xm_cell
Definition variables.h:188
Vec lUcont_o
Definition variables.h:844
PetscInt64 bboxGuessFallbackCount
Definition variables.h:233
InterpolationMethod interpolationMethod
Definition variables.h:744
RankCellInfo * RankCellInfoMap
Definition variables.h:881
PetscInt ym_cell
Definition variables.h:188
Vec Ucat_o
Definition variables.h:844
PetscInt64 bboxGuessSuccessCount
Definition variables.h:232
BoundingBox * bboxlist
Definition variables.h:742
PetscMPIInt rank_xm
Definition variables.h:180
PetscInt64 maxParticlePassDepth
Definition variables.h:234
PetscScalar z
Definition variables.h:101
@ INTERP_CORNER_AVERAGED
Definition variables.h:523
Vec Ucat
Definition variables.h:837
Vec ParticleCount
Definition variables.h:882
PetscReal Const_CS
Definition variables.h:734
Vec Ucont_o
Definition variables.h:844
PetscInt JM
Definition variables.h:820
Cmpnts InitialConstantContra
Definition variables.h:697
PetscMPIInt rank_zp
Definition variables.h:182
Vec Ucont_rm1
Definition variables.h:845
SearchMetricsState searchMetrics
Definition variables.h:752
PetscReal diffusivity
Definition variables.h:173
Vec lUcont
Definition variables.h:837
PetscInt step
Definition variables.h:651
Vec Diffusivity
Definition variables.h:840
PetscRandom BrownianMotionRNG
Definition variables.h:753
PetscInt GridOrientation
Definition variables.h:824
DMDALocalInfo info
Definition variables.h:818
Vec lUcat
Definition variables.h:837
PetscScalar y
Definition variables.h:101
PetscMPIInt size
Definition variables.h:647
PetscInt IM
Definition variables.h:820
Cmpnts weights
Definition variables.h:170
@ NUM_FACES
Definition variables.h:145
PetscInt les
Definition variables.h:732
Vec lDiffusivity
Definition variables.h:840
PetscInt64 searchAttempts
Definition variables.h:222
PetscInt64 PID
Definition variables.h:166
PetscInt64 maxTraversalFailCount
Definition variables.h:229
Vec Psi
Definition variables.h:883
Vec P_o
Definition variables.h:844
@ BC_FACE_NEG_Z
Definition variables.h:247
Defines a 3D axis-aligned bounding box.
Definition variables.h:154
A 3D point or vector with PetscScalar components.
Definition variables.h:100
Defines a particle's core properties for Lagrangian tracking.
Definition variables.h:165
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
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.