PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
Loading...
Searching...
No Matches
Functions
test_periodic_dev.c File Reference

Non-gating periodic-boundary development harnesses. More...

#include "test_support.h"
#include "Boundaries.h"
Include dependency graph for test_periodic_dev.c:

Go to the source code of this file.

Functions

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 DestroyBoundaryHandler (BoundaryCondition **bc_ptr)
 Destroys one boundary-condition handler allocated by a periodic test.
 
static void MarkXPeriodic (UserCtx *user)
 Marks the x faces as periodic for periodic-transfer harnesses.
 
static PetscErrorCode TestPeriodicGeometricFactoryAssignment (void)
 Tests periodic geometric factory construction in the non-gating periodic harness.
 
static PetscErrorCode TestTransferPeriodicFaceFieldCopiesXFaces (void)
 Tests direct periodic face transfer on the staggered velocity field.
 
static PetscErrorCode TestApplyMetricsPeriodicBCsCopiesAjFaces (void)
 Tests periodic metric transfer through the aggregate periodic-metrics helper.
 
static PetscErrorCode TestPeriodicDrivenConstantHandlerBehavior (void)
 Tests periodic driven-flow controller initialization, sensing, and trim application.
 
static PetscErrorCode TestPeriodicDrivenConstantRejectsNonPeriodicFace (void)
 Tests that the periodic driven handler rejects non-periodic faces during initialization.
 
int main (int argc, char **argv)
 Runs the non-gating periodic development PETSc test binary.
 

Detailed Description

Non-gating periodic-boundary development harnesses.

Definition in file test_periodic_dev.c.

Function Documentation

◆ AppendBCParam()

static PetscErrorCode AppendBCParam ( BC_Param **  head,
const char *  key,
const char *  value 
)
static

Appends one key/value pair to a linked list of boundary-condition parameters.

Definition at line 13 of file test_periodic_dev.c.

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}
struct BC_Param_s * next
Definition variables.h:266
A node in a linked list for storing key-value parameters from the bcs.dat file.
Definition variables.h:263
A generic C-style linked list node for integers.
Definition variables.h:366
Here is the caller graph for this function:

◆ DestroyBoundaryHandler()

static PetscErrorCode DestroyBoundaryHandler ( BoundaryCondition **  bc_ptr)
static

Destroys one boundary-condition handler allocated by a periodic test.

Definition at line 32 of file test_periodic_dev.c.

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}
The "virtual table" struct for a boundary condition handler object.
Definition variables.h:280
PetscErrorCode(* Destroy)(BoundaryCondition *self)
Definition variables.h:289
Here is the caller graph for this function:

◆ MarkXPeriodic()

static void MarkXPeriodic ( UserCtx user)
static

Marks the x faces as periodic for periodic-transfer harnesses.

Definition at line 50 of file test_periodic_dev.c.

Here is the caller graph for this function:

◆ TestPeriodicGeometricFactoryAssignment()

static PetscErrorCode TestPeriodicGeometricFactoryAssignment ( void  )
static

Tests periodic geometric factory construction in the non-gating periodic harness.

Definition at line 62 of file test_periodic_dev.c.

63{
64 BoundaryCondition *bc = NULL;
65
66 PetscFunctionBeginUser;
68 PetscCall(PicurvAssertBool((PetscBool)(bc != NULL), "periodic geometric factory should allocate a handler"));
69 PetscCall(PicurvAssertIntEqual(BC_PRIORITY_WALL, bc->priority, "periodic geometric handler priority"));
70 PetscCall(PicurvAssertBool((PetscBool)(bc->Apply == NULL), "periodic geometric handler should not expose Apply"));
71 PetscCall(PicurvAssertBool((PetscBool)(bc->Initialize == NULL), "periodic geometric handler should not expose Initialize"));
72 PetscCall(DestroyBoundaryHandler(&bc));
73 PetscFunctionReturn(0);
74}
PetscErrorCode BoundaryCondition_Create(BCHandlerType handler_type, BoundaryCondition **new_bc_ptr)
(Private) Creates and configures a specific BoundaryCondition handler object.
Definition Boundaries.c:744
PetscErrorCode(* Initialize)(BoundaryCondition *self, BCContext *ctx)
Definition variables.h:284
PetscErrorCode(* Apply)(BoundaryCondition *self, BCContext *ctx)
Definition variables.h:286
BCPriorityType priority
Definition variables.h:282
static PetscErrorCode DestroyBoundaryHandler(BoundaryCondition **bc_ptr)
Destroys one boundary-condition handler allocated by a periodic test.
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.
@ BC_PRIORITY_WALL
Definition variables.h:254
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestTransferPeriodicFaceFieldCopiesXFaces()

static PetscErrorCode TestTransferPeriodicFaceFieldCopiesXFaces ( void  )
static

Tests direct periodic face transfer on the staggered velocity field.

Definition at line 78 of file test_periodic_dev.c.

79{
80 SimCtx *simCtx = NULL;
81 UserCtx *user = NULL;
82 Cmpnts ***ucont = NULL;
83 PetscReal expected_neg_face = 0.0;
84 PetscReal expected_pos_face = 0.0;
85 PetscReal expected_neg_ghost = 0.0;
86 PetscReal expected_pos_ghost = 0.0;
87
88 PetscFunctionBeginUser;
89 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
90 MarkXPeriodic(user);
91
92 PetscCall(DMDAVecGetArray(user->fda, user->Ucont, &ucont));
93 for (PetscInt k = user->info.zs; k < user->info.zs + user->info.zm; ++k) {
94 for (PetscInt j = user->info.ys; j < user->info.ys + user->info.ym; ++j) {
95 for (PetscInt i = user->info.xs; i < user->info.xs + user->info.xm; ++i) {
96 ucont[k][j][i].x = (PetscReal)i;
97 ucont[k][j][i].y = 10.0 + (PetscReal)i;
98 ucont[k][j][i].z = 20.0 + (PetscReal)i;
99 }
100 }
101 }
102 PetscCall(DMDAVecRestoreArray(user->fda, user->Ucont, &ucont));
103 PetscCall(DMGlobalToLocalBegin(user->fda, user->Ucont, INSERT_VALUES, user->lUcont));
104 PetscCall(DMGlobalToLocalEnd(user->fda, user->Ucont, INSERT_VALUES, user->lUcont));
105 PetscCall(DMDAVecGetArrayRead(user->fda, user->lUcont, &ucont));
106 expected_neg_face = ucont[2][2][user->info.mx - 2].x;
107 expected_pos_face = ucont[2][2][1].x;
108 expected_neg_ghost = ucont[2][2][user->info.mx - 3].x;
109 expected_pos_ghost = ucont[2][2][2].x;
110 PetscCall(DMDAVecRestoreArrayRead(user->fda, user->lUcont, &ucont));
111
112 PetscCall(TransferPeriodicFaceField(user, "Ucont"));
113
114 PetscCall(DMDAVecGetArrayRead(user->fda, user->lUcont, &ucont));
115 PetscCall(PicurvAssertRealNear(expected_neg_face, ucont[2][2][0].x, 1.0e-12, "NEG_X periodic face should copy from the opposite interior face"));
116 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"));
117 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"));
118 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"));
119 PetscCall(DMDAVecRestoreArrayRead(user->fda, user->lUcont, &ucont));
120
121 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
122 PetscFunctionReturn(0);
123}
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.
static void MarkXPeriodic(UserCtx *user)
Marks the x faces as periodic for periodic-transfer harnesses.
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.
Vec Ucont
Definition variables.h:790
PetscScalar x
Definition variables.h:101
PetscScalar z
Definition variables.h:101
Vec lUcont
Definition variables.h:790
DMDALocalInfo info
Definition variables.h:771
PetscScalar y
Definition variables.h:101
A 3D point or vector with PetscScalar components.
Definition variables.h:100
The master context for the entire simulation.
Definition variables.h:602
User-defined context containing data specific to a single computational grid level.
Definition variables.h:764
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestApplyMetricsPeriodicBCsCopiesAjFaces()

static PetscErrorCode TestApplyMetricsPeriodicBCsCopiesAjFaces ( void  )
static

Tests periodic metric transfer through the aggregate periodic-metrics helper.

Definition at line 127 of file test_periodic_dev.c.

128{
129 SimCtx *simCtx = NULL;
130 UserCtx *user = NULL;
131 PetscReal ***aj = NULL;
132 PetscReal expected_neg_face = 0.0;
133 PetscReal expected_pos_face = 0.0;
134 PetscReal expected_neg_ghost = 0.0;
135 PetscReal expected_pos_ghost = 0.0;
136
137 PetscFunctionBeginUser;
138 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
139 MarkXPeriodic(user);
140
141 PetscCall(DMDAVecGetArray(user->da, user->Aj, &aj));
142 for (PetscInt k = user->info.zs; k < user->info.zs + user->info.zm; ++k) {
143 for (PetscInt j = user->info.ys; j < user->info.ys + user->info.ym; ++j) {
144 for (PetscInt i = user->info.xs; i < user->info.xs + user->info.xm; ++i) {
145 aj[k][j][i] = 100.0 + (PetscReal)i;
146 }
147 }
148 }
149 PetscCall(DMDAVecRestoreArray(user->da, user->Aj, &aj));
150 PetscCall(DMGlobalToLocalBegin(user->da, user->Aj, INSERT_VALUES, user->lAj));
151 PetscCall(DMGlobalToLocalEnd(user->da, user->Aj, INSERT_VALUES, user->lAj));
152 PetscCall(DMDAVecGetArrayRead(user->da, user->lAj, &aj));
153 expected_neg_face = aj[2][2][user->info.mx - 2];
154 expected_pos_face = aj[2][2][1];
155 expected_neg_ghost = aj[2][2][user->info.mx - 3];
156 expected_pos_ghost = aj[2][2][2];
157 PetscCall(DMDAVecRestoreArrayRead(user->da, user->lAj, &aj));
158
159 PetscCall(ApplyMetricsPeriodicBCs(user));
160
161 PetscCall(DMDAVecGetArrayRead(user->da, user->lAj, &aj));
162 PetscCall(PicurvAssertRealNear(expected_neg_face, aj[2][2][0], 1.0e-12, "NEG_X periodic metric face should copy from the opposite interior face"));
163 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"));
164 PetscCall(PicurvAssertRealNear(expected_neg_ghost, aj[2][2][-1], 1.0e-12, "NEG_X periodic metric ghost should copy the second-to-last interior face"));
165 PetscCall(PicurvAssertRealNear(expected_pos_ghost, aj[2][2][user->info.mx], 1.0e-12, "POS_X periodic metric ghost should copy the second interior face"));
166 PetscCall(DMDAVecRestoreArrayRead(user->da, user->lAj, &aj));
167
168 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
169 PetscFunctionReturn(0);
170}
PetscErrorCode ApplyMetricsPeriodicBCs(UserCtx *user)
(Orchestrator) Updates all metric-related fields in the local ghost cell regions for periodic boundar...
Vec lAj
Definition variables.h:811
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestPeriodicDrivenConstantHandlerBehavior()

static PetscErrorCode TestPeriodicDrivenConstantHandlerBehavior ( void  )
static

Tests periodic driven-flow controller initialization, sensing, and trim application.

Definition at line 174 of file test_periodic_dev.c.

175{
176 SimCtx *simCtx = NULL;
177 UserCtx *user = NULL;
178 BoundaryCondition *bc = NULL;
179 BCContext ctx;
180 PetscReal dummy_inflow = 0.0;
181 PetscReal dummy_outflow = 0.0;
182 Cmpnts ***ucont = NULL;
183 Cmpnts ***uch = NULL;
184
185 PetscFunctionBeginUser;
186 PetscCall(PetscMemzero(&ctx, sizeof(ctx)));
187 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 6, 6, 6));
188 PetscCall(PicurvPopulateIdentityMetrics(user));
189 PetscCall(VecSet(user->Ucont, 1.0));
190 PetscCall(DMGlobalToLocalBegin(user->fda, user->Ucont, INSERT_VALUES, user->lUcont));
191 PetscCall(DMGlobalToLocalEnd(user->fda, user->Ucont, INSERT_VALUES, user->lUcont));
192
199 PetscCall(AppendBCParam(&user->boundary_faces[BC_FACE_NEG_Z].params, "target_flux", "30.0"));
200 PetscCall(AppendBCParam(&user->boundary_faces[BC_FACE_NEG_Z].params, "apply_trim", "true"));
201
202 ctx.user = user;
205 PetscCall(bc->Initialize(bc, &ctx));
206 PetscCall(bc->PreStep(bc, &ctx, &dummy_inflow, &dummy_outflow));
207 PetscCall(PicurvAssertRealNear(30.0, simCtx->targetVolumetricFlux, 1.0e-12, "periodic driven initialization should store the target flux"));
208 PetscCall(PicurvAssertRealNear(0.2, simCtx->bulkVelocityCorrection, 1.0e-12, "periodic driven controller should compute the bulk correction from the measured flux"));
209
210 PetscCall(bc->Apply(bc, &ctx));
211 PetscCall(DMDAVecGetArrayRead(user->fda, user->lUcont, &ucont));
212 PetscCall(DMDAVecGetArrayRead(user->fda, user->Bcs.Uch, &uch));
213 PetscCall(PicurvAssertRealNear(1.2, ucont[0][3][3].z, 1.0e-12, "periodic driven trim should update the boundary face flux"));
214 PetscCall(PicurvAssertRealNear(0.2, uch[0][3][3].z, 1.0e-12, "periodic driven trim should be recorded in Uch"));
215 PetscCall(DMDAVecRestoreArrayRead(user->fda, user->lUcont, &ucont));
216 PetscCall(DMDAVecRestoreArrayRead(user->fda, user->Bcs.Uch, &uch));
217
218 PetscCall(DestroyBoundaryHandler(&bc));
220 user->boundary_faces[BC_FACE_NEG_Z].params = NULL;
221 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
222 PetscFunctionReturn(0);
223}
void FreeBC_ParamList(BC_Param *head)
Frees an entire linked list of boundary-condition parameters.
Definition io.c:302
PetscErrorCode(* PreStep)(BoundaryCondition *self, BCContext *ctx, PetscReal *local_inflow, PetscReal *local_outflow)
Definition variables.h:285
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.
PetscErrorCode PicurvPopulateIdentityMetrics(UserCtx *user)
Populates identity metric vectors on the minimal grid fixture.
PetscReal targetVolumetricFlux
Definition variables.h:679
@ BC_HANDLER_PERIODIC_DRIVEN_CONSTANT_FLUX
Definition variables.h:245
PetscReal bulkVelocityCorrection
Definition variables.h:680
BCFace face_id
Definition variables.h:272
BCS Bcs
Definition variables.h:785
UserCtx * user
Definition variables.h:271
BC_Param * params
Definition variables.h:297
Vec Uch
Characteristic velocity for boundary conditions.
Definition variables.h:122
@ BC_FACE_POS_Z
Definition variables.h:206
@ BC_FACE_NEG_Z
Definition variables.h:206
Provides execution context for a boundary condition handler.
Definition variables.h:270
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestPeriodicDrivenConstantRejectsNonPeriodicFace()

static PetscErrorCode TestPeriodicDrivenConstantRejectsNonPeriodicFace ( void  )
static

Tests that the periodic driven handler rejects non-periodic faces during initialization.

Definition at line 227 of file test_periodic_dev.c.

228{
229 SimCtx *simCtx = NULL;
230 UserCtx *user = NULL;
231 BoundaryCondition *bc = NULL;
232 BCContext ctx;
233 PetscErrorCode ierr_init = 0;
234
235 PetscFunctionBeginUser;
236 PetscCall(PetscMemzero(&ctx, sizeof(ctx)));
237 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 6, 6, 6));
241 PetscCall(AppendBCParam(&user->boundary_faces[BC_FACE_NEG_Z].params, "target_flux", "5.0"));
242 ctx.user = user;
244
246 PetscCall(PetscPushErrorHandler(PetscIgnoreErrorHandler, NULL));
247 ierr_init = bc->Initialize(bc, &ctx);
248 PetscCall(PetscPopErrorHandler());
249 PetscCall(PicurvAssertBool((PetscBool)(ierr_init != 0), "periodic driven initialization should reject non-periodic faces"));
250
251 PetscCall(DestroyBoundaryHandler(&bc));
253 user->boundary_faces[BC_FACE_NEG_Z].params = NULL;
254 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
255 PetscFunctionReturn(0);
256}
@ WALL
Definition variables.h:213
Here is the call graph for this function:
Here is the caller graph for this function:

◆ main()

int main ( int  argc,
char **  argv 
)

Runs the non-gating periodic development PETSc test binary.

Definition at line 260 of file test_periodic_dev.c.

261{
262 PetscErrorCode ierr;
263 const PicurvTestCase cases[] = {
264 {"periodic-geometric-factory-assignment", TestPeriodicGeometricFactoryAssignment},
265 {"transfer-periodic-face-field-copies-x-faces", TestTransferPeriodicFaceFieldCopiesXFaces},
266 {"apply-metrics-periodic-bcs-copies-aj-faces", TestApplyMetricsPeriodicBCsCopiesAjFaces},
267 {"periodic-driven-constant-handler-behavior", TestPeriodicDrivenConstantHandlerBehavior},
268 {"periodic-driven-constant-rejects-non-periodic-face", TestPeriodicDrivenConstantRejectsNonPeriodicFace},
269 };
270
271 ierr = PetscInitialize(&argc, &argv, NULL, "PICurv periodic development tests");
272 if (ierr) {
273 return (int)ierr;
274 }
275
276 ierr = PicurvRunTests("unit-periodic-dev", cases, sizeof(cases) / sizeof(cases[0]));
277 if (ierr) {
278 PetscFinalize();
279 return (int)ierr;
280 }
281
282 ierr = PetscFinalize();
283 return (int)ierr;
284}
static PetscErrorCode TestPeriodicDrivenConstantRejectsNonPeriodicFace(void)
Tests that the periodic driven handler rejects non-periodic faces during initialization.
static PetscErrorCode TestPeriodicGeometricFactoryAssignment(void)
Tests periodic geometric factory construction in the non-gating periodic harness.
static PetscErrorCode TestApplyMetricsPeriodicBCsCopiesAjFaces(void)
Tests periodic metric transfer through the aggregate periodic-metrics helper.
static PetscErrorCode TestTransferPeriodicFaceFieldCopiesXFaces(void)
Tests direct periodic face transfer on the staggered velocity field.
static PetscErrorCode TestPeriodicDrivenConstantHandlerBehavior(void)
Tests periodic driven-flow controller initialization, sensing, and trim application.
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.
Named test case descriptor consumed by PicurvRunTests.
Here is the call graph for this function: