PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
Loading...
Searching...
No Matches
test_periodic_dev.c
Go to the documentation of this file.
1/**
2 * @file test_periodic_dev.c
3 * @brief Non-gating periodic-boundary development harnesses.
4 */
5
6#include "test_support.h"
7
8#include "Boundaries.h"
9
10/**
11 * @brief Appends one key/value pair to a linked list of boundary-condition parameters.
12 */
13static PetscErrorCode AppendBCParam(BC_Param **head, const char *key, const char *value)
14{
15 BC_Param *node = NULL;
16 BC_Param **cursor = head;
17
18 PetscFunctionBeginUser;
19 PetscCheck(head != NULL, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "BC param list head cannot be NULL.");
20 PetscCall(PetscCalloc1(1, &node));
21 PetscCall(PetscStrallocpy(key, &node->key));
22 PetscCall(PetscStrallocpy(value, &node->value));
23 while (*cursor) {
24 cursor = &(*cursor)->next;
25 }
26 *cursor = node;
27 PetscFunctionReturn(0);
28}
29/**
30 * @brief Destroys one boundary-condition handler allocated by a periodic test.
31 */
32static PetscErrorCode DestroyBoundaryHandler(BoundaryCondition **bc_ptr)
33{
34 BoundaryCondition *bc = NULL;
35
36 PetscFunctionBeginUser;
37 if (!bc_ptr || !*bc_ptr) PetscFunctionReturn(0);
38
39 bc = *bc_ptr;
40 if (bc->Destroy) {
41 PetscCall(bc->Destroy(bc));
42 }
43 PetscCall(PetscFree(bc));
44 *bc_ptr = NULL;
45 PetscFunctionReturn(0);
46}
47/**
48 * @brief Marks the x faces as periodic for periodic-transfer harnesses.
49 */
59/**
60 * @brief Marks the y faces as periodic for periodic-transfer harnesses.
61 */
71/**
72 * @brief Tests periodic configuration rejects an unpaired geometric periodic face.
73 */
75{
76 SimCtx *simCtx = NULL;
77 UserCtx *user = NULL;
78 PetscErrorCode ierr_validate = 0;
79
80 PetscFunctionBeginUser;
81 PetscCall(PicurvCreateMinimalContextsWithPeriodicity(&simCtx, &user, 4, 4, 4, PETSC_TRUE, PETSC_FALSE, PETSC_FALSE));
84
85 PetscCall(PetscPushErrorHandler(PetscIgnoreErrorHandler, NULL));
86 ierr_validate = BoundarySystem_Validate(user);
87 PetscCall(PetscPopErrorHandler());
88
89 PetscCall(PicurvAssertBool((PetscBool)(ierr_validate != 0), "unpaired periodic faces should be rejected"));
90 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
91 PetscFunctionReturn(0);
92}
93/**
94 * @brief Tests constant translational periodic geometry is accepted and stored.
95 */
96static PetscErrorCode TestPeriodicGeometryStoresTranslation(void)
97{
98 SimCtx *simCtx = NULL;
99 UserCtx *user = NULL;
100 Vec gcoor = NULL;
101 const Cmpnts ***coor = NULL;
102 PetscReal expected_x = 0.0;
103
104 PetscFunctionBeginUser;
105 PetscCall(PicurvCreateMinimalContextsWithPeriodicity(&simCtx, &user, 4, 4, 4, PETSC_TRUE, PETSC_FALSE, PETSC_FALSE));
106 MarkXPeriodic(user);
107
108 PetscCall(DMGetCoordinates(user->da, &gcoor));
109 PetscCall(DMDAVecGetArrayRead(user->fda, gcoor, &coor));
110 expected_x = coor[0][0][user->info.mx - 2].x - coor[0][0][0].x;
111 PetscCall(DMDAVecRestoreArrayRead(user->fda, gcoor, &coor));
112
113 PetscCall(ValidatePeriodicGeometry(user));
114 PetscCall(PicurvAssertBool(user->periodic_translation_valid[0], "X-periodic translation should be marked valid"));
115 PetscCall(PicurvAssertBool((PetscBool)!user->periodic_translation_valid[1], "Y translation should remain invalid"));
116 PetscCall(PicurvAssertBool((PetscBool)!user->periodic_translation_valid[2], "Z translation should remain invalid"));
117 PetscCall(PicurvAssertRealNear(expected_x, user->periodic_translation[0].x, 1.0e-12, "stored X translation"));
118 PetscCall(PicurvAssertRealNear(0.0, user->periodic_translation[0].y, 1.0e-12, "stored X translation y component"));
119 PetscCall(PicurvAssertRealNear(0.0, user->periodic_translation[0].z, 1.0e-12, "stored X translation z component"));
120
121 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
122 PetscFunctionReturn(0);
123}
124/**
125 * @brief Tests mixed periodic directions are validated and stored independently.
126 */
128{
129 SimCtx *simCtx = NULL;
130 UserCtx *user = NULL;
131
132 PetscFunctionBeginUser;
133 PetscCall(PicurvCreateMinimalContextsWithPeriodicity(&simCtx, &user, 4, 6, 4, PETSC_TRUE, PETSC_TRUE, PETSC_FALSE));
134 MarkXPeriodic(user);
135 MarkYPeriodic(user);
136
137 PetscCall(ValidatePeriodicGeometry(user));
138 PetscCall(PicurvAssertBool(user->periodic_translation_valid[0], "mixed-periodic X translation should be valid"));
139 PetscCall(PicurvAssertBool(user->periodic_translation_valid[1], "mixed-periodic Y translation should be valid"));
140 PetscCall(PicurvAssertBool((PetscBool)!user->periodic_translation_valid[2], "mixed-periodic Z translation should remain invalid"));
141 PetscCall(PicurvAssertBool((PetscBool)(PetscAbsReal(user->periodic_translation[0].x) > 0.0),
142 "mixed-periodic X translation should have nonzero x component"));
143 PetscCall(PicurvAssertBool((PetscBool)(PetscAbsReal(user->periodic_translation[1].y) > 0.0),
144 "mixed-periodic Y translation should have nonzero y component"));
145
146 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
147 PetscFunctionReturn(0);
148}
149/**
150 * @brief Tests non-translational seam geometry fails before metric construction.
151 */
153{
154 SimCtx *simCtx = NULL;
155 UserCtx *user = NULL;
156 Vec gcoor = NULL;
157 Cmpnts ***coor = NULL;
158 PetscErrorCode ierr_validate = 0;
159
160 PetscFunctionBeginUser;
161 PetscCall(PicurvCreateMinimalContextsWithPeriodicity(&simCtx, &user, 4, 4, 4, PETSC_TRUE, PETSC_FALSE, PETSC_FALSE));
162 MarkXPeriodic(user);
163
164 PetscCall(DMGetCoordinates(user->da, &gcoor));
165 PetscCall(DMDAVecGetArray(user->fda, gcoor, &coor));
166 coor[1][1][user->info.mx - 2].y += 0.25;
167 PetscCall(DMDAVecRestoreArray(user->fda, gcoor, &coor));
168 PetscCall(UpdateLocalGhosts(user, "Coordinates"));
169
170 PetscCall(PetscPushErrorHandler(PetscIgnoreErrorHandler, NULL));
171 ierr_validate = ValidatePeriodicGeometry(user);
172 PetscCall(PetscPopErrorHandler());
173
174 PetscCall(PicurvAssertBool((PetscBool)(ierr_validate != 0), "varying seam translation should be rejected"));
175 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
176 PetscFunctionReturn(0);
177}
178
179/**
180 * @brief Tests face-center coordinates remain geometrically continuous across a periodic seam.
181 */
183{
184 SimCtx *simCtx = NULL;
185 UserCtx *user = NULL;
186 Cmpnts ***centx = NULL;
187 const Cmpnts ***lcentx = NULL;
188 const char *fields[] = {"Centx"};
189 PetscReal translation, spacing;
190
191 PetscFunctionBeginUser;
192 PetscCall(PicurvCreateMinimalContextsWithPeriodicity(&simCtx, &user, 6, 4, 4, PETSC_TRUE, PETSC_FALSE, PETSC_FALSE));
193 MarkXPeriodic(user);
194 PetscCall(ValidatePeriodicGeometry(user));
195 translation = user->periodic_translation[0].x;
196 spacing = translation / (PetscReal)(user->info.mx - 2);
197
198 PetscCall(DMDAVecGetArray(user->fda, user->Centx, &centx));
199 for (PetscInt k = user->info.zs; k < user->info.zs + user->info.zm; k++) {
200 for (PetscInt j = user->info.ys; j < user->info.ys + user->info.ym; j++) {
201 for (PetscInt i = user->info.xs; i < user->info.xs + user->info.xm; i++) {
202 centx[k][j][i] = (Cmpnts){spacing * i, 2.0 + j, 3.0 + k};
203 }
204 }
205 }
206 PetscCall(DMDAVecRestoreArray(user->fda, user->Centx, &centx));
207
208 PetscCall(SynchronizePeriodicFaceFields(user, 'i', 1, fields));
209 PetscCall(DMDAVecGetArrayRead(user->fda, user->lCentx, &lcentx));
210 PetscCall(PicurvAssertRealNear(-spacing, lcentx[2][2][-1].x, 1.0e-12,
211 "translated Centx negative adjacent ghost"));
212 PetscCall(PicurvAssertRealNear(0.0, lcentx[2][2][0].x, 1.0e-12,
213 "translated Centx negative endpoint"));
214 PetscCall(PicurvAssertRealNear(translation + spacing, lcentx[2][2][user->info.mx - 1].x, 1.0e-12,
215 "translated Centx positive endpoint"));
216 PetscCall(PicurvAssertRealNear(translation + 2.0 * spacing, lcentx[2][2][user->info.mx].x, 1.0e-12,
217 "translated Centx positive adjacent ghost"));
218 PetscCall(DMDAVecRestoreArrayRead(user->fda, user->lCentx, &lcentx));
219
220 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
221 PetscFunctionReturn(0);
222}
223
224/**
225 * @brief Tests the dedicated QUICK outer-ghost repair for cell-centered inputs.
226 */
227static PetscErrorCode TestPeriodicQuickStencilPreparation(void)
228{
229 SimCtx *simCtx = NULL;
230 UserCtx *user = NULL;
231 Cmpnts ***ucat = NULL;
232 PetscReal ***nvert = NULL;
233
234 PetscFunctionBeginUser;
235 PetscCall(PicurvCreateMinimalContextsWithPeriodicity(&simCtx, &user, 6, 4, 4, PETSC_TRUE, PETSC_FALSE, PETSC_FALSE));
236 MarkXPeriodic(user);
237
238 PetscCall(DMDAVecGetArray(user->fda, user->Ucat, &ucat));
239 PetscCall(DMDAVecGetArray(user->da, user->Nvert, &nvert));
240 for (PetscInt k = user->info.zs; k < user->info.zs + user->info.zm; k++) {
241 for (PetscInt j = user->info.ys; j < user->info.ys + user->info.ym; j++) {
242 for (PetscInt i = user->info.xs; i < user->info.xs + user->info.xm; i++) {
243 ucat[k][j][i] = (Cmpnts){100.0 + i, 200.0 + i, 300.0 + i};
244 nvert[k][j][i] = 400.0 + i;
245 }
246 }
247 }
248 PetscCall(DMDAVecRestoreArray(user->da, user->Nvert, &nvert));
249 PetscCall(DMDAVecRestoreArray(user->fda, user->Ucat, &ucat));
250 PetscCall(UpdateLocalGhosts(user, "Ucat"));
251 PetscCall(UpdateLocalGhosts(user, "Nvert"));
252 PetscCall(PreparePeriodicQuickStencilFields(user, user->lUcat, user->lNvert));
253
254 PetscCall(DMDAVecGetArrayRead(user->fda, user->lUcat, &ucat));
255 PetscCall(DMDAVecGetArrayRead(user->da, user->lNvert, &nvert));
256 PetscCall(PicurvAssertRealNear(100.0 + user->info.mx - 3, ucat[2][2][-1].x, 1.0e-12,
257 "QUICK Ucat negative outer ghost"));
258 PetscCall(PicurvAssertRealNear(102.0, ucat[2][2][user->info.mx].x, 1.0e-12,
259 "QUICK Ucat positive outer ghost"));
260 PetscCall(PicurvAssertRealNear(400.0 + user->info.mx - 3, nvert[2][2][-1], 1.0e-12,
261 "QUICK Nvert negative outer ghost"));
262 PetscCall(PicurvAssertRealNear(402.0, nvert[2][2][user->info.mx], 1.0e-12,
263 "QUICK Nvert positive outer ghost"));
264 PetscCall(DMDAVecRestoreArrayRead(user->da, user->lNvert, &nvert));
265 PetscCall(DMDAVecRestoreArrayRead(user->fda, user->lUcat, &ucat));
266
267 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
268 PetscFunctionReturn(0);
269}
270/**
271 * @brief Tests periodic geometric factory construction in the non-gating periodic harness.
272 */
273static PetscErrorCode TestPeriodicGeometricFactoryAssignment(void)
274{
275 BoundaryCondition *bc = NULL;
276
277 PetscFunctionBeginUser;
279 PetscCall(PicurvAssertBool((PetscBool)(bc != NULL), "periodic geometric factory should allocate a handler"));
280 PetscCall(PicurvAssertIntEqual(BC_PRIORITY_WALL, bc->priority, "periodic geometric handler priority"));
281 PetscCall(PicurvAssertBool((PetscBool)(bc->Apply == NULL), "periodic geometric handler should not expose Apply"));
282 PetscCall(PicurvAssertBool((PetscBool)(bc->Initialize == NULL), "periodic geometric handler should not expose Initialize"));
283 PetscCall(DestroyBoundaryHandler(&bc));
284 PetscFunctionReturn(0);
285}
286/**
287 * @brief Tests direct periodic face transfer on the staggered velocity field.
288 */
290{
291 SimCtx *simCtx = NULL;
292 UserCtx *user = NULL;
293 Cmpnts ***ucont = NULL;
294 PetscReal expected_neg_face = 0.0;
295 PetscReal expected_pos_face = 0.0;
296 PetscReal expected_neg_ghost = 0.0;
297 PetscReal expected_pos_ghost = 0.0;
298
299 PetscFunctionBeginUser;
300 PetscCall(PicurvCreateMinimalContextsWithPeriodicity(&simCtx, &user, 4, 4, 4, PETSC_TRUE, PETSC_FALSE, PETSC_FALSE));
301 MarkXPeriodic(user);
302
303 PetscCall(DMDAVecGetArray(user->fda, user->Ucont, &ucont));
304 for (PetscInt k = user->info.zs; k < user->info.zs + user->info.zm; ++k) {
305 for (PetscInt j = user->info.ys; j < user->info.ys + user->info.ym; ++j) {
306 for (PetscInt i = user->info.xs; i < user->info.xs + user->info.xm; ++i) {
307 ucont[k][j][i].x = (PetscReal)i;
308 ucont[k][j][i].y = 10.0 + (PetscReal)i;
309 ucont[k][j][i].z = 20.0 + (PetscReal)i;
310 }
311 }
312 }
313 PetscCall(DMDAVecRestoreArray(user->fda, user->Ucont, &ucont));
314 PetscCall(DMGlobalToLocalBegin(user->fda, user->Ucont, INSERT_VALUES, user->lUcont));
315 PetscCall(DMGlobalToLocalEnd(user->fda, user->Ucont, INSERT_VALUES, user->lUcont));
316 PetscCall(DMDAVecGetArrayRead(user->fda, user->lUcont, &ucont));
317 expected_neg_face = ucont[2][2][user->info.mx - 2].x;
318 expected_pos_face = ucont[2][2][1].x;
319 expected_neg_ghost = ucont[2][2][user->info.mx - 3].x;
320 expected_pos_ghost = ucont[2][2][2].x;
321 PetscCall(DMDAVecRestoreArrayRead(user->fda, user->lUcont, &ucont));
322
323 PetscCall(TransferPeriodicFaceField(user, "Ucont"));
324
325 PetscCall(DMDAVecGetArrayRead(user->fda, user->lUcont, &ucont));
326 PetscCall(PicurvAssertRealNear(expected_neg_face, ucont[2][2][0].x, 1.0e-12, "NEG_X periodic face should copy from the opposite interior face"));
327 PetscCall(PicurvAssertRealNear(expected_pos_face, ucont[2][2][user->info.mx - 1].x, 1.0e-12, "POS_X periodic face should copy from the leading interior face"));
328 PetscCall(PicurvAssertRealNear(expected_neg_ghost, ucont[2][2][-1].x, 1.0e-12, "NEG_X ghost face should copy the second-to-last interior face"));
329 PetscCall(PicurvAssertRealNear(expected_pos_ghost, ucont[2][2][user->info.mx].x, 1.0e-12, "POS_X ghost face should copy the second interior face"));
330 PetscCall(DMDAVecRestoreArrayRead(user->fda, user->lUcont, &ucont));
331
332 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
333 PetscFunctionReturn(0);
334}
335/**
336 * @brief Tests periodic metric transfer through the aggregate periodic-metrics helper.
337 */
339{
340 SimCtx *simCtx = NULL;
341 UserCtx *user = NULL;
342 PetscReal ***aj = NULL;
343 PetscReal expected_neg_face = 0.0;
344 PetscReal expected_pos_face = 0.0;
345 PetscReal expected_neg_ghost = 0.0;
346 PetscReal expected_pos_ghost = 0.0;
347
348 PetscFunctionBeginUser;
349 PetscCall(PicurvCreateMinimalContextsWithPeriodicity(&simCtx, &user, 4, 4, 4, PETSC_TRUE, PETSC_FALSE, PETSC_FALSE));
350 MarkXPeriodic(user);
351 PetscCall(ValidatePeriodicGeometry(user));
352
353 PetscCall(DMDAVecGetArray(user->da, user->Aj, &aj));
354 for (PetscInt k = user->info.zs; k < user->info.zs + user->info.zm; ++k) {
355 for (PetscInt j = user->info.ys; j < user->info.ys + user->info.ym; ++j) {
356 for (PetscInt i = user->info.xs; i < user->info.xs + user->info.xm; ++i) {
357 aj[k][j][i] = 100.0 + (PetscReal)i;
358 }
359 }
360 }
361 PetscCall(DMDAVecRestoreArray(user->da, user->Aj, &aj));
362 PetscCall(DMGlobalToLocalBegin(user->da, user->Aj, INSERT_VALUES, user->lAj));
363 PetscCall(DMGlobalToLocalEnd(user->da, user->Aj, INSERT_VALUES, user->lAj));
364 PetscCall(DMDAVecGetArrayRead(user->da, user->lAj, &aj));
365 expected_neg_face = aj[2][2][user->info.mx - 2];
366 expected_pos_face = aj[2][2][1];
367 expected_neg_ghost = expected_pos_face;
368 expected_pos_ghost = expected_neg_face;
369 PetscCall(DMDAVecRestoreArrayRead(user->da, user->lAj, &aj));
370
371 PetscCall(ApplyMetricsPeriodicBCs(user));
372
373 PetscCall(DMDAVecGetArrayRead(user->da, user->lAj, &aj));
374 PetscCall(PicurvAssertRealNear(expected_neg_face, aj[2][2][0], 1.0e-12, "NEG_X periodic metric face should copy from the opposite interior face"));
375 PetscCall(PicurvAssertRealNear(expected_pos_face, aj[2][2][user->info.mx - 1], 1.0e-12, "POS_X periodic metric face should copy from the leading interior face"));
376 PetscCall(PicurvAssertRealNear(expected_neg_ghost, aj[2][2][-1], 1.0e-12, "cell-centered Aj negative ghost should retain PETSc wraparound"));
377 PetscCall(PicurvAssertRealNear(expected_pos_ghost, aj[2][2][user->info.mx], 1.0e-12, "cell-centered Aj positive ghost should retain PETSc wraparound"));
378 PetscCall(DMDAVecRestoreArrayRead(user->da, user->lAj, &aj));
379
380 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
381 PetscFunctionReturn(0);
382}
383/**
384 * @brief Tests ordered cell-field synchronization across two periodic directions.
385 */
387{
388 SimCtx *simCtx = NULL;
389 UserCtx *user = NULL;
390 PetscReal ***p = NULL;
391 PetscReal ***cs = NULL;
392 PetscReal ***diffusivity = NULL;
393 Cmpnts ***ucat = NULL;
394 PetscInt mx, my;
395 const char *fields[] = {"Ucat", "P", "Phi", "Nvert", "Nu_t", "CS", "Diffusivity"};
396
397 PetscFunctionBeginUser;
398 PetscCall(PicurvCreateMinimalContextsWithPeriodicity(&simCtx, &user, 4, 4, 4, PETSC_TRUE, PETSC_TRUE, PETSC_FALSE));
399 MarkXPeriodic(user);
400 MarkYPeriodic(user);
401
402 PetscCall(DMCreateGlobalVector(user->da, &user->Nu_t));
403 PetscCall(DMCreateLocalVector(user->da, &user->lNu_t));
404 PetscCall(VecSet(user->Nu_t, 5.0));
405 PetscCall(VecSet(user->lNu_t, 0.0));
406 PetscCall(DMCreateGlobalVector(user->da, &user->CS));
407 PetscCall(DMCreateLocalVector(user->da, &user->lCs));
408 PetscCall(VecSet(user->Phi, 3.0));
409 PetscCall(VecSet(user->Nvert, 4.0));
410
411 PetscCall(DMDAVecGetArray(user->da, user->P, &p));
412 PetscCall(DMDAVecGetArray(user->da, user->CS, &cs));
413 PetscCall(DMDAVecGetArray(user->da, user->Diffusivity, &diffusivity));
414 PetscCall(DMDAVecGetArray(user->fda, user->Ucat, &ucat));
415 for (PetscInt k = user->info.zs; k < user->info.zs + user->info.zm; ++k) {
416 for (PetscInt j = user->info.ys; j < user->info.ys + user->info.ym; ++j) {
417 for (PetscInt i = user->info.xs; i < user->info.xs + user->info.xm; ++i) {
418 PetscReal value = (PetscReal)(i + 10 * j + 100 * k);
419 p[k][j][i] = value;
420 cs[k][j][i] = value + 4000.0;
421 diffusivity[k][j][i] = value + 5000.0;
422 ucat[k][j][i].x = value + 1000.0;
423 ucat[k][j][i].y = value + 2000.0;
424 ucat[k][j][i].z = value + 3000.0;
425 }
426 }
427 }
428 PetscCall(DMDAVecRestoreArray(user->da, user->P, &p));
429 PetscCall(DMDAVecRestoreArray(user->da, user->CS, &cs));
430 PetscCall(DMDAVecRestoreArray(user->da, user->Diffusivity, &diffusivity));
431 PetscCall(DMDAVecRestoreArray(user->fda, user->Ucat, &ucat));
432
433 PetscCall(SynchronizePeriodicCellFields(user, 7, fields));
434
435 mx = user->info.mx;
436 my = user->info.my;
437 PetscCall(DMDAVecGetArrayRead(user->da, user->P, &p));
438 PetscCall(DMDAVecGetArrayRead(user->da, user->CS, &cs));
439 PetscCall(DMDAVecGetArrayRead(user->da, user->Diffusivity, &diffusivity));
440 PetscCall(DMDAVecGetArrayRead(user->fda, user->Ucat, &ucat));
441 PetscCall(PicurvAssertRealNear((PetscReal)((mx - 2) + 20 + 200), p[2][2][0], 1.0e-12,
442 "x-periodic negative endpoint should copy the opposite interior value"));
443 PetscCall(PicurvAssertRealNear((PetscReal)(1 + 20 + 200), p[2][2][mx - 1], 1.0e-12,
444 "x-periodic positive endpoint should copy the leading interior value"));
445 PetscCall(PicurvAssertRealNear((PetscReal)(2 + 10 * (my - 2) + 200), p[2][0][2], 1.0e-12,
446 "y-periodic negative endpoint should copy the opposite interior value"));
447 PetscCall(PicurvAssertRealNear((PetscReal)((mx - 2) + 10 * (my - 2) + 200), p[2][0][0], 1.0e-12,
448 "ordered synchronization should propagate the opposite periodic corner"));
449 PetscCall(PicurvAssertRealNear((PetscReal)((mx - 2) + 10 * (my - 2) + 1200), ucat[2][0][0].x, 1.0e-12,
450 "vector fields should use the same ordered periodic-corner protocol"));
451 PetscCall(PicurvAssertRealNear((PetscReal)((mx - 2) + 10 * (my - 2) + 4200), cs[2][0][0], 1.0e-12,
452 "CS should use the ordered periodic-corner protocol"));
453 PetscCall(PicurvAssertRealNear((PetscReal)((mx - 2) + 10 * (my - 2) + 5200), diffusivity[2][0][0], 1.0e-12,
454 "Diffusivity should use the ordered periodic-corner protocol"));
455 PetscCall(DMDAVecRestoreArrayRead(user->da, user->P, &p));
456 PetscCall(DMDAVecRestoreArrayRead(user->da, user->CS, &cs));
457 PetscCall(DMDAVecRestoreArrayRead(user->da, user->Diffusivity, &diffusivity));
458 PetscCall(DMDAVecRestoreArrayRead(user->fda, user->Ucat, &ucat));
459 PetscCall(PicurvAssertVecConstant(user->Nu_t, 5.0, 1.0e-12, "Nu_t should be accepted by periodic cell synchronization"));
460
461 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
462 PetscFunctionReturn(0);
463}
464
465/**
466 * @brief Tests ordered persistent synchronization for an I-face scalar metric.
467 */
469{
470 SimCtx *simCtx = NULL;
471 UserCtx *user = NULL;
472 PetscReal ***iaj = NULL;
473 Cmpnts ***csi = NULL;
474 const char *fields[] = {"Csi", "IAj"};
475 PetscInt mx, my;
476
477 PetscFunctionBeginUser;
478 PetscCall(PicurvCreateMinimalContextsWithPeriodicity(&simCtx, &user, 4, 4, 4, PETSC_TRUE, PETSC_TRUE, PETSC_FALSE));
479 MarkXPeriodic(user);
480 MarkYPeriodic(user);
481 mx = user->info.mx;
482 my = user->info.my;
483
484 PetscCall(DMDAVecGetArray(user->da, user->IAj, &iaj));
485 PetscCall(DMDAVecGetArray(user->fda, user->Csi, &csi));
486 for (PetscInt k = user->info.zs; k < user->info.zs + user->info.zm; k++) {
487 for (PetscInt j = user->info.ys; j < user->info.ys + user->info.ym; j++) {
488 for (PetscInt i = user->info.xs; i < user->info.xs + user->info.xm; i++) {
489 PetscReal value = (PetscReal)(i + 10 * j + 100 * k);
490 iaj[k][j][i] = value;
491 csi[k][j][i] = (Cmpnts){value + 1000.0, value + 2000.0, value + 3000.0};
492 }
493 }
494 }
495 PetscCall(DMDAVecRestoreArray(user->da, user->IAj, &iaj));
496 PetscCall(DMDAVecRestoreArray(user->fda, user->Csi, &csi));
497
498 PetscCall(SynchronizePeriodicFaceFields(user, 'i', 2, fields));
499
500 PetscCall(DMDAVecGetArrayRead(user->da, user->IAj, &iaj));
501 PetscCall(DMDAVecGetArrayRead(user->fda, user->Csi, &csi));
502 PetscCall(PicurvAssertRealNear((PetscReal)(mx - 2 + 20 + 200), iaj[2][2][0], 1.0e-12,
503 "I-face negative seam should copy the opposite physical seam face"));
504 PetscCall(PicurvAssertRealNear((PetscReal)(1 + 20 + 200), iaj[2][2][mx - 1], 1.0e-12,
505 "I-face positive dummy should copy the leading physical face"));
506 PetscCall(PicurvAssertRealNear((PetscReal)(2 + 10 * (my - 2) + 200), iaj[2][0][2], 1.0e-12,
507 "I-face field should use cell-style synchronization tangentially"));
508 PetscCall(PicurvAssertRealNear((PetscReal)(mx - 2 + 10 * (my - 2) + 200), iaj[2][0][0], 1.0e-12,
509 "ordered face synchronization should propagate mixed-axis corners"));
510 PetscCall(PicurvAssertRealNear((PetscReal)(mx - 2 + 10 * (my - 2) + 1200), csi[2][0][0].x, 1.0e-12,
511 "vector I-face fields should use the same ordered synchronization"));
512 PetscCall(DMDAVecRestoreArrayRead(user->da, user->IAj, &iaj));
513 PetscCall(DMDAVecRestoreArrayRead(user->fda, user->Csi, &csi));
514
515 PetscCall(DMDAVecGetArrayRead(user->da, user->lIAj, &iaj));
516 PetscCall(PicurvAssertRealNear((PetscReal)(mx - 3 + 20 + 200), iaj[2][2][-1], 1.0e-12,
517 "I-face negative adjacent ghost should be shifted to the prior face"));
518 PetscCall(PicurvAssertRealNear((PetscReal)(2 + 20 + 200), iaj[2][2][mx], 1.0e-12,
519 "I-face positive adjacent ghost should be shifted to the next face"));
520 PetscCall(PicurvAssertRealNear((PetscReal)(2 + 10 + 200), iaj[2][-1][2], 1.0e-12,
521 "I-face tangential ghosts should retain cell-style wraparound"));
522 PetscCall(DMDAVecRestoreArrayRead(user->da, user->lIAj, &iaj));
523
524 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
525 PetscFunctionReturn(0);
526}
527
528/**
529 * @brief Tests ordered persistent synchronization for component-staggered Ucont.
530 */
532{
533 SimCtx *simCtx = NULL;
534 UserCtx *user = NULL;
535 Cmpnts ***ucont = NULL;
536 const char *fields[] = {"Ucont"};
537 PetscInt mx, my;
538
539 PetscFunctionBeginUser;
540 PetscCall(PicurvCreateMinimalContextsWithPeriodicity(&simCtx, &user, 4, 4, 4, PETSC_TRUE, PETSC_TRUE, PETSC_FALSE));
541 MarkXPeriodic(user);
542 MarkYPeriodic(user);
543 mx = user->info.mx;
544 my = user->info.my;
545
546 PetscCall(DMDAVecGetArray(user->fda, user->Ucont, &ucont));
547 for (PetscInt k = user->info.zs; k < user->info.zs + user->info.zm; k++) {
548 for (PetscInt j = user->info.ys; j < user->info.ys + user->info.ym; j++) {
549 for (PetscInt i = user->info.xs; i < user->info.xs + user->info.xm; i++) {
550 PetscReal value = (PetscReal)(i + 10 * j + 100 * k);
551 ucont[k][j][i] = (Cmpnts){value + 1000.0, value + 2000.0, value + 3000.0};
552 }
553 }
554 }
555 PetscCall(DMDAVecRestoreArray(user->fda, user->Ucont, &ucont));
556
557 PetscCall(SynchronizePeriodicStaggeredFields(user, 1, fields));
558
559 PetscCall(DMDAVecGetArrayRead(user->fda, user->Ucont, &ucont));
560 PetscCall(PicurvAssertRealNear((PetscReal)(mx - 2 + 20 + 1200), ucont[2][2][0].x, 1.0e-12,
561 "Ucont.x negative X seam should copy the opposite physical seam"));
562 PetscCall(PicurvAssertRealNear((PetscReal)(1 + 20 + 2200), ucont[2][2][mx - 1].y, 1.0e-12,
563 "Ucont.y positive X dummy should copy the leading physical value"));
564 PetscCall(PicurvAssertRealNear((PetscReal)(2 + 10 * (my - 2) + 3200), ucont[2][0][2].z, 1.0e-12,
565 "Ucont.z should synchronize tangentially in Y"));
566 PetscCall(PicurvAssertRealNear((PetscReal)(mx - 2 + 10 * (my - 2) + 2200), ucont[2][0][0].y, 1.0e-12,
567 "ordered staggered synchronization should propagate mixed-axis corners"));
568 PetscCall(DMDAVecRestoreArrayRead(user->fda, user->Ucont, &ucont));
569
570 PetscCall(DMDAVecGetArrayRead(user->fda, user->lUcont, &ucont));
571 PetscCall(PicurvAssertRealNear((PetscReal)(mx - 3 + 20 + 1200), ucont[2][2][-1].x, 1.0e-12,
572 "Ucont.x should use face-normal X ghost repair"));
573 PetscCall(PicurvAssertRealNear((PetscReal)(1 + 20 + 2200), ucont[2][2][-1].y, 1.0e-12,
574 "Ucont.y should retain cell-style wraparound tangentially in X"));
575 PetscCall(PicurvAssertRealNear((PetscReal)(2 + 20 + 1200), ucont[2][2][mx].x, 1.0e-12,
576 "Ucont.x positive adjacent ghost should use the next face"));
577 PetscCall(PicurvAssertRealNear((PetscReal)(mx - 2 + 20 + 2200), ucont[2][2][mx].y, 1.0e-12,
578 "Ucont.y positive X ghost should retain cell-style wraparound"));
579 PetscCall(DMDAVecRestoreArrayRead(user->fda, user->lUcont, &ucont));
580
581 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
582 PetscFunctionReturn(0);
583}
584
585/**
586 * @brief Tests mixed-boundary post-projection cell finalization and Ucont preservation.
587 */
589{
590 SimCtx *simCtx = NULL;
591 UserCtx *user = NULL;
592 Vec ucont_before = NULL;
593 Vec ucont_difference = NULL;
594 Cmpnts ***ucat = NULL;
595 Cmpnts ***ucont = NULL;
596 PetscReal ***p = NULL;
597 PetscReal difference_norm = 0.0;
598 PetscInt mx;
599
600 PetscFunctionBeginUser;
601 PetscCall(PicurvCreateMinimalContextsWithPeriodicity(&simCtx, &user, 4, 4, 4, PETSC_TRUE, PETSC_FALSE, PETSC_FALSE));
602 MarkXPeriodic(user);
603 PetscCall(VecSet(user->Bcs.Ubcs, 0.0));
604
605 PetscCall(DMDAVecGetArray(user->fda, user->Ucat, &ucat));
606 PetscCall(DMDAVecGetArray(user->fda, user->Ucont, &ucont));
607 PetscCall(DMDAVecGetArray(user->da, user->P, &p));
608 for (PetscInt k = user->info.zs; k < user->info.zs + user->info.zm; ++k) {
609 for (PetscInt j = user->info.ys; j < user->info.ys + user->info.ym; ++j) {
610 for (PetscInt i = user->info.xs; i < user->info.xs + user->info.xm; ++i) {
611 PetscReal value = (PetscReal)(i + 10 * j + 100 * k);
612 ucat[k][j][i].x = value + 1000.0;
613 ucat[k][j][i].y = value + 2000.0;
614 ucat[k][j][i].z = value + 3000.0;
615 ucont[k][j][i].x = value + 4000.0;
616 ucont[k][j][i].y = value + 5000.0;
617 ucont[k][j][i].z = value + 6000.0;
618 p[k][j][i] = value;
619 }
620 }
621 }
622 PetscCall(DMDAVecRestoreArray(user->fda, user->Ucat, &ucat));
623 PetscCall(DMDAVecRestoreArray(user->fda, user->Ucont, &ucont));
624 PetscCall(DMDAVecRestoreArray(user->da, user->P, &p));
625
626 PetscCall(VecDuplicate(user->Ucont, &ucont_before));
627 PetscCall(VecDuplicate(user->Ucont, &ucont_difference));
628 PetscCall(VecCopy(user->Ucont, ucont_before));
629
630 PetscCall(FinalizePostProjectionCellFields(user));
631
632 PetscCall(VecWAXPY(ucont_difference, -1.0, user->Ucont, ucont_before));
633 PetscCall(VecNorm(ucont_difference, NORM_INFINITY, &difference_norm));
634 PetscCall(PicurvAssertRealNear(0.0, difference_norm, 1.0e-12,
635 "post-projection cell finalization must preserve mixed-boundary Ucont"));
636
637 mx = user->info.mx;
638 PetscCall(DMDAVecGetArrayRead(user->fda, user->Ucat, &ucat));
639 PetscCall(DMDAVecGetArrayRead(user->da, user->P, &p));
640 PetscCall(PicurvAssertRealNear((PetscReal)((mx - 2) + 20 + 200 + 1000), ucat[2][2][0].x, 1.0e-12,
641 "mixed finalization should restore the negative x-periodic Ucat endpoint"));
642 PetscCall(PicurvAssertRealNear((PetscReal)(-(2 + 10 + 200 + 1000)), ucat[2][0][2].x, 1.0e-12,
643 "mixed finalization should extrapolate the non-periodic y dummy face"));
644 PetscCall(PicurvAssertRealNear((PetscReal)(-((mx - 2) + 10 + 200 + 1000)), ucat[2][0][0].x, 1.0e-12,
645 "periodic synchronization should win at a mixed x-periodic/y-physical edge"));
646 PetscCall(PicurvAssertRealNear((PetscReal)((mx - 2) + 20 + 200), p[2][2][0], 1.0e-12,
647 "mixed finalization should restore the negative x-periodic pressure endpoint"));
648 PetscCall(DMDAVecRestoreArrayRead(user->fda, user->Ucat, &ucat));
649 PetscCall(DMDAVecRestoreArrayRead(user->da, user->P, &p));
650
651 PetscCall(VecDestroy(&ucont_difference));
652 PetscCall(VecDestroy(&ucont_before));
653 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
654 PetscFunctionReturn(0);
655}
656/**
657 * @brief Tests periodic driven-flow controller initialization, sensing, and trim application.
658 */
660{
661 SimCtx *simCtx = NULL;
662 UserCtx *user = NULL;
663 BoundaryCondition *bc = NULL;
664 BCContext ctx;
665 PetscReal dummy_inflow = 0.0;
666 PetscReal dummy_outflow = 0.0;
667 Cmpnts ***ucont = NULL;
668 Cmpnts ***uch = NULL;
669
670 PetscFunctionBeginUser;
671 PetscCall(PetscMemzero(&ctx, sizeof(ctx)));
672 PetscCall(PicurvCreateMinimalContextsWithPeriodicity(&simCtx, &user, 6, 6, 6, PETSC_FALSE, PETSC_FALSE, PETSC_TRUE));
673 PetscCall(PicurvPopulateIdentityMetrics(user));
674 PetscCall(VecSet(user->Ucont, 1.0));
675 PetscCall(DMGlobalToLocalBegin(user->fda, user->Ucont, INSERT_VALUES, user->lUcont));
676 PetscCall(DMGlobalToLocalEnd(user->fda, user->Ucont, INSERT_VALUES, user->lUcont));
677
684 PetscCall(AppendBCParam(&user->boundary_faces[BC_FACE_NEG_Z].params, "target_flux", "30.0"));
685 PetscCall(AppendBCParam(&user->boundary_faces[BC_FACE_NEG_Z].params, "apply_trim", "true"));
686
687 ctx.user = user;
690 PetscCall(bc->Initialize(bc, &ctx));
691 PetscCall(bc->PreStep(bc, &ctx, &dummy_inflow, &dummy_outflow));
692 PetscCall(PicurvAssertRealNear(30.0, simCtx->targetVolumetricFlux, 1.0e-12, "periodic driven initialization should store the target flux"));
693 PetscCall(PicurvAssertRealNear(0.2, simCtx->bulkVelocityCorrection, 1.0e-12, "periodic driven controller should compute the bulk correction from the measured flux"));
694
695 PetscCall(bc->Apply(bc, &ctx));
696 PetscCall(DMDAVecGetArrayRead(user->fda, user->lUcont, &ucont));
697 PetscCall(DMDAVecGetArrayRead(user->fda, user->Bcs.Uch, &uch));
698 PetscCall(PicurvAssertRealNear(1.2, ucont[0][3][3].z, 1.0e-12, "periodic driven trim should update the boundary face flux"));
699 PetscCall(PicurvAssertRealNear(0.2, uch[0][3][3].z, 1.0e-12, "periodic driven trim should be recorded in Uch"));
700 PetscCall(DMDAVecRestoreArrayRead(user->fda, user->lUcont, &ucont));
701 PetscCall(DMDAVecRestoreArrayRead(user->fda, user->Bcs.Uch, &uch));
702
703 PetscCall(DestroyBoundaryHandler(&bc));
705 user->boundary_faces[BC_FACE_NEG_Z].params = NULL;
706 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
707 PetscFunctionReturn(0);
708}
709/**
710 * @brief Tests that the periodic driven handler rejects non-periodic faces during initialization.
711 */
713{
714 SimCtx *simCtx = NULL;
715 UserCtx *user = NULL;
716 BoundaryCondition *bc = NULL;
717 BCContext ctx;
718 PetscErrorCode ierr_init = 0;
719
720 PetscFunctionBeginUser;
721 PetscCall(PetscMemzero(&ctx, sizeof(ctx)));
722 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 6, 6, 6));
726 PetscCall(AppendBCParam(&user->boundary_faces[BC_FACE_NEG_Z].params, "target_flux", "5.0"));
727 ctx.user = user;
729
731 PetscCall(PetscPushErrorHandler(PetscIgnoreErrorHandler, NULL));
732 ierr_init = bc->Initialize(bc, &ctx);
733 PetscCall(PetscPopErrorHandler());
734 PetscCall(PicurvAssertBool((PetscBool)(ierr_init != 0), "periodic driven initialization should reject non-periodic faces"));
735
736 PetscCall(DestroyBoundaryHandler(&bc));
738 user->boundary_faces[BC_FACE_NEG_Z].params = NULL;
739 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
740 PetscFunctionReturn(0);
741}
742/**
743 * @brief Runs the non-gating periodic development PETSc test binary.
744 */
745int main(int argc, char **argv)
746{
747 PetscErrorCode ierr;
748 const PicurvTestCase cases[] = {
749 {"periodic-geometric-factory-assignment", TestPeriodicGeometricFactoryAssignment},
750 {"periodic-configuration-requires-paired-faces", TestPeriodicConfigurationRequiresPairedFaces},
751 {"periodic-geometry-stores-translation", TestPeriodicGeometryStoresTranslation},
752 {"periodic-geometry-stores-mixed-translations", TestPeriodicGeometryStoresMixedTranslations},
753 {"periodic-geometry-rejects-varying-translation", TestPeriodicGeometryRejectsVaryingTranslation},
754 {"periodic-face-center-coordinate-synchronization", TestPeriodicFaceCenterCoordinateSynchronization},
755 {"periodic-quick-stencil-preparation", TestPeriodicQuickStencilPreparation},
756 {"transfer-periodic-face-field-copies-x-faces", TestTransferPeriodicFaceFieldCopiesXFaces},
757 {"apply-metrics-periodic-bcs-synchronizes-aj", TestApplyMetricsPeriodicBCsSynchronizesAj},
758 {"synchronize-periodic-cell-fields-copies-mixed-axes", TestSynchronizePeriodicCellFieldsCopiesMixedAxes},
759 {"synchronize-periodic-face-fields-copies-mixed-axes", TestSynchronizePeriodicFaceFieldsCopiesMixedAxes},
760 {"synchronize-periodic-staggered-fields-copies-mixed-axes", TestSynchronizePeriodicStaggeredFieldsCopiesMixedAxes},
761 {"finalize-post-projection-cell-fields-mixed-boundaries", TestFinalizePostProjectionCellFieldsMixedBoundaries},
762 {"periodic-driven-constant-handler-behavior", TestPeriodicDrivenConstantHandlerBehavior},
763 {"periodic-driven-constant-rejects-non-periodic-face", TestPeriodicDrivenConstantRejectsNonPeriodicFace},
764 };
765
766 ierr = PetscInitialize(&argc, &argv, NULL, "PICurv periodic development tests");
767 if (ierr) {
768 return (int)ierr;
769 }
770
771 ierr = PicurvRunTests("unit-periodic-dev", cases, sizeof(cases) / sizeof(cases[0]));
772 if (ierr) {
773 PetscFinalize();
774 return (int)ierr;
775 }
776
777 ierr = PetscFinalize();
778 return (int)ierr;
779}
PetscErrorCode PreparePeriodicQuickStencilFields(UserCtx *user, Vec local_vector_field, Vec local_scalar_field)
Repairs the outer adjacent periodic ghosts used by QUICK cell stencils.
PetscErrorCode TransferPeriodicFaceField(UserCtx *user, const char *field_name)
(Primitive) Copies periodic data from the interior to the local ghost cell region for a single field.
PetscErrorCode ApplyMetricsPeriodicBCs(UserCtx *user)
(Orchestrator) Updates all metric-related fields in the local ghost cell regions for periodic boundar...
PetscErrorCode SynchronizePeriodicStaggeredFields(UserCtx *user, PetscInt num_fields, const char *field_names[])
Synchronizes persistent component-staggered vector fields.
PetscErrorCode SynchronizePeriodicFaceFields(UserCtx *user, char face_direction, PetscInt num_fields, const char *field_names[])
Synchronizes persistent fields belonging to one face family.
PetscErrorCode BoundarySystem_Validate(UserCtx *user)
(Public) Validates the consistency and compatibility of the parsed boundary condition system.
Definition Boundaries.c:830
PetscErrorCode BoundaryCondition_Create(BCHandlerType handler_type, BoundaryCondition **new_bc_ptr)
(Private) Creates and configures a specific BoundaryCondition handler object.
Definition Boundaries.c:744
PetscErrorCode FinalizePostProjectionCellFields(UserCtx *user)
Finalizes cell-centered fields after the projection step.
PetscErrorCode SynchronizePeriodicCellFields(UserCtx *user, PetscInt num_fields, const char *field_names[])
Synchronizes periodic endpoint cells for a list of cell-centered fields.
PetscErrorCode ValidatePeriodicGeometry(UserCtx *user)
Validates that configured geometric periodic seams match by translation.
Definition grid.c:380
void FreeBC_ParamList(BC_Param *head)
Frees an entire linked list of boundary-condition parameters.
Definition io.c:302
PetscErrorCode UpdateLocalGhosts(UserCtx *user, const char *fieldName)
Updates the local vector (including ghost points) from its corresponding global vector.
Definition setup.c:1707
The "virtual table" struct for a boundary condition handler object.
Definition variables.h:321
PetscErrorCode(* PreStep)(BoundaryCondition *self, BCContext *ctx, PetscReal *local_inflow, PetscReal *local_outflow)
Definition variables.h:326
PetscErrorCode(* Destroy)(BoundaryCondition *self)
Definition variables.h:330
PetscErrorCode(* Initialize)(BoundaryCondition *self, BCContext *ctx)
Definition variables.h:325
PetscErrorCode(* Apply)(BoundaryCondition *self, BCContext *ctx)
Definition variables.h:327
BCPriorityType priority
Definition variables.h:323
static void MarkXPeriodic(UserCtx *user)
Marks the x faces as periodic for periodic-transfer harnesses.
static PetscErrorCode TestPeriodicDrivenConstantRejectsNonPeriodicFace(void)
Tests that the periodic driven handler rejects non-periodic faces during initialization.
static PetscErrorCode AppendBCParam(BC_Param **head, const char *key, const char *value)
Appends one key/value pair to a linked list of boundary-condition parameters.
static PetscErrorCode TestPeriodicGeometricFactoryAssignment(void)
Tests periodic geometric factory construction in the non-gating periodic harness.
static PetscErrorCode TestPeriodicFaceCenterCoordinateSynchronization(void)
Tests face-center coordinates remain geometrically continuous across a periodic seam.
int main(int argc, char **argv)
Runs the non-gating periodic development PETSc test binary.
static PetscErrorCode TestPeriodicGeometryStoresTranslation(void)
Tests constant translational periodic geometry is accepted and stored.
static PetscErrorCode TestSynchronizePeriodicFaceFieldsCopiesMixedAxes(void)
Tests ordered persistent synchronization for an I-face scalar metric.
static PetscErrorCode TestFinalizePostProjectionCellFieldsMixedBoundaries(void)
Tests mixed-boundary post-projection cell finalization and Ucont preservation.
static PetscErrorCode TestPeriodicConfigurationRequiresPairedFaces(void)
Tests periodic configuration rejects an unpaired geometric periodic face.
static PetscErrorCode TestPeriodicQuickStencilPreparation(void)
Tests the dedicated QUICK outer-ghost repair for cell-centered inputs.
static PetscErrorCode TestSynchronizePeriodicCellFieldsCopiesMixedAxes(void)
Tests ordered cell-field synchronization across two periodic directions.
static PetscErrorCode TestPeriodicGeometryRejectsVaryingTranslation(void)
Tests non-translational seam geometry fails before metric construction.
static PetscErrorCode TestSynchronizePeriodicStaggeredFieldsCopiesMixedAxes(void)
Tests ordered persistent synchronization for component-staggered Ucont.
static PetscErrorCode DestroyBoundaryHandler(BoundaryCondition **bc_ptr)
Destroys one boundary-condition handler allocated by a periodic test.
static PetscErrorCode TestPeriodicGeometryStoresMixedTranslations(void)
Tests mixed periodic directions are validated and stored independently.
static PetscErrorCode TestTransferPeriodicFaceFieldCopiesXFaces(void)
Tests direct periodic face transfer on the staggered velocity field.
static PetscErrorCode TestApplyMetricsPeriodicBCsSynchronizesAj(void)
Tests periodic metric transfer through the aggregate periodic-metrics helper.
static void MarkYPeriodic(UserCtx *user)
Marks the y faces as periodic for periodic-transfer harnesses.
static PetscErrorCode TestPeriodicDrivenConstantHandlerBehavior(void)
Tests periodic driven-flow controller initialization, sensing, and trim application.
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 PicurvCreateMinimalContextsWithPeriodicity(SimCtx **simCtx_out, UserCtx **user_out, PetscInt mx, PetscInt my, PetscInt mz, PetscBool x_periodic, PetscBool y_periodic, PetscBool z_periodic)
Builds minimal SimCtx and UserCtx fixtures for C unit tests with configurable periodicity.
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 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 PicurvPopulateIdentityMetrics(UserCtx *user)
Populates identity metric vectors on the minimal grid fixture.
PetscErrorCode PicurvAssertBool(PetscBool value, const char *context)
Asserts that one boolean condition is true.
Shared declarations for the PICurv C test fixture and assertion layer.
Named test case descriptor consumed by PicurvRunTests.
@ PERIODIC
Definition variables.h:260
@ WALL
Definition variables.h:254
BoundaryFaceConfig boundary_faces[6]
Definition variables.h:856
PetscReal targetVolumetricFlux
Definition variables.h:740
Vec lNvert
Definition variables.h:864
Vec Phi
Definition variables.h:864
struct BC_Param_s * next
Definition variables.h:307
Vec Csi
Definition variables.h:887
@ BC_HANDLER_PERIODIC_GEOMETRIC
Definition variables.h:284
@ BC_HANDLER_PERIODIC_DRIVEN_CONSTANT_FLUX
Definition variables.h:286
BCHandlerType handler_type
Definition variables.h:337
Vec lIAj
Definition variables.h:890
PetscReal bulkVelocityCorrection
Definition variables.h:741
BCFace face_id
Definition variables.h:313
Vec lCs
Definition variables.h:895
Vec Ucont
Definition variables.h:864
Vec Ubcs
Physical Cartesian velocity at boundary faces. Full 3D array but only boundary-face entries are meani...
Definition variables.h:121
PetscScalar x
Definition variables.h:101
Vec Centx
Definition variables.h:888
BCS Bcs
Definition variables.h:859
UserCtx * user
Definition variables.h:312
Vec lNu_t
Definition variables.h:895
Vec Nu_t
Definition variables.h:895
BC_Param * params
Definition variables.h:338
PetscScalar z
Definition variables.h:101
Vec Ucat
Definition variables.h:864
Vec IAj
Definition variables.h:890
Vec lCentx
Definition variables.h:889
Vec lUcont
Definition variables.h:864
Vec Diffusivity
Definition variables.h:867
Vec lAj
Definition variables.h:887
DMDALocalInfo info
Definition variables.h:843
@ BC_PRIORITY_WALL
Definition variables.h:295
Vec lUcat
Definition variables.h:864
PetscScalar y
Definition variables.h:101
Cmpnts periodic_translation[3]
Definition variables.h:852
Vec Nvert
Definition variables.h:864
PetscBool periodic_translation_valid[3]
Definition variables.h:853
BCType mathematical_type
Definition variables.h:336
Vec Uch
Characteristic velocity for boundary conditions.
Definition variables.h:122
@ BC_FACE_NEG_X
Definition variables.h:245
@ BC_FACE_POS_Z
Definition variables.h:247
@ BC_FACE_POS_Y
Definition variables.h:246
@ BC_FACE_NEG_Z
Definition variables.h:247
@ BC_FACE_POS_X
Definition variables.h:245
@ BC_FACE_NEG_Y
Definition variables.h:246
Provides execution context for a boundary condition handler.
Definition variables.h:311
A node in a linked list for storing key-value parameters from the bcs.dat file.
Definition variables.h:304
A 3D point or vector with PetscScalar components.
Definition variables.h:100
The master context for the entire simulation.
Definition variables.h:653
User-defined context containing data specific to a single computational grid level.
Definition variables.h:836
A generic C-style linked list node for integers.
Definition variables.h:407