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

C unit tests for boundary factories, inlet ownership, and face helpers. More...

#include "test_support.h"
#include "Boundaries.h"
Include dependency graph for test_boundaries.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 boundary test.
 
static PetscInt InteriorSampleIndex (PetscInt nodes)
 Chooses a stable interior node index for face-sample assertions on tiny test grids.
 
static PetscReal GetFaceNormalComponent (Cmpnts value, BCFace face)
 Extracts the velocity component aligned with the supplied boundary-face normal.
 
static PetscReal GetFaceOrientationSign (BCFace face)
 Returns the sign convention used for face-normal flux expectations.
 
static const char * GetInletParamKey (BCFace face)
 Maps an inlet face to the matching configuration key used by the handler parser.
 
static PetscReal GetFaceInteriorPointCount (const UserCtx *user, BCFace face)
 Computes the number of face-interior sample points for a given boundary face.
 
static PetscErrorCode GetRepresentativeFaceSlots (const UserCtx *user, BCFace face, PetscInt *ucont_k, PetscInt *ucont_j, PetscInt *ucont_i, PetscInt *ubcs_k, PetscInt *ubcs_j, PetscInt *ubcs_i)
 Selects representative Ucont and Ubcs slots for face-matrix assertions.
 
static PetscErrorCode GetParabolicSampleSlots (const UserCtx *user, BCFace face, PetscInt *center_k, PetscInt *center_j, PetscInt *center_i, PetscInt *off_k, PetscInt *off_j, PetscInt *off_i, PetscInt *wall_k, PetscInt *wall_j, PetscInt *wall_i, PetscInt *ubcs_center_k, PetscInt *ubcs_center_j, PetscInt *ubcs_center_i)
 Selects center, off-center, and wall-adjacent slots for parabolic-face checks.
 
static void ResetBoundaryFaceConfig (BoundaryFaceConfig *cfg)
 Resets one boundary-face configuration entry to a neutral test-local baseline.
 
static PetscErrorCode TestCanRankServiceFaceMatchesInletWhenDefined (void)
 Tests that face-service detection matches a defined inlet face.
 
static PetscErrorCode TestCanRankServiceInletFaceRequiresDefinition (void)
 Tests that inlet-face service requires an inlet face to be defined.
 
static PetscErrorCode TestBoundaryConditionFactoryAssignments (void)
 Tests the boundary-condition factory assignments for representative handlers.
 
static PetscErrorCode TestBoundaryConditionFactoryImplementedHandlerMatrix (void)
 Tests the implemented-handler matrix exposed by the factory.
 
static PetscErrorCode TestBoundaryConditionFactoryRejectsUnsupportedHandler (void)
 Tests that unsupported handler types are rejected by the factory.
 
static PetscErrorCode TestGetDeterministicFaceGridLocationFaceMatrix (void)
 Tests deterministic inlet-face grid-location helpers across all faces.
 
static PetscErrorCode TestGetRandomCellAndLogicalCoordsOnInletFaceMatrix (void)
 Tests random inlet-face cell selection across all supported faces.
 
static PetscErrorCode TestInletConstantVelocityHandlerBehavior (void)
 Tests constant inlet handler initialization and face writes on a Z inlet.
 
static PetscErrorCode TestWallNoSlipHandlerFaceMatrix (void)
 Tests wall no-slip application across the full face matrix.
 
static PetscErrorCode TestInletConstantVelocityHandlerFaceMatrix (void)
 Tests constant-inlet initialization, flux accounting, and face writes across all faces.
 
static PetscErrorCode TestInletParabolicProfileHandlerBehavior (void)
 Tests parabolic inlet handler shape on a tiny Z-face.
 
static PetscErrorCode TestInletParabolicProfileHandlerFaceMatrix (void)
 Tests parabolic-inlet centerline, wall, and flux behavior across all faces.
 
static PetscErrorCode TestOutletConservationHandlerBehavior (void)
 Tests outlet conservation handler correction and post-step flux accounting.
 
static PetscErrorCode TestOutletConservationHandlerFaceMatrix (void)
 Tests outlet-conservation correction and flux accounting across all outlet faces.
 
int main (int argc, char **argv)
 Runs the unit-boundaries PETSc test binary.
 

Detailed Description

C unit tests for boundary factories, inlet ownership, and face helpers.

Definition in file test_boundaries.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_boundaries.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:307
A node in a linked list for storing key-value parameters from the bcs.dat file.
Definition variables.h:304
A generic C-style linked list node for integers.
Definition variables.h:407
Here is the caller graph for this function:

◆ DestroyBoundaryHandler()

static PetscErrorCode DestroyBoundaryHandler ( BoundaryCondition **  bc_ptr)
static

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

Definition at line 33 of file test_boundaries.c.

34{
35 BoundaryCondition *bc = NULL;
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:321
PetscErrorCode(* Destroy)(BoundaryCondition *self)
Definition variables.h:330
Here is the caller graph for this function:

◆ InteriorSampleIndex()

static PetscInt InteriorSampleIndex ( PetscInt  nodes)
static

Chooses a stable interior node index for face-sample assertions on tiny test grids.

Definition at line 51 of file test_boundaries.c.

52{
53 PetscInt index = 3;
54
55 if (index >= nodes - 1) {
56 index = nodes - 2;
57 }
58 if (index < 1) {
59 index = 1;
60 }
61 return index;
62}
Here is the caller graph for this function:

◆ GetFaceNormalComponent()

static PetscReal GetFaceNormalComponent ( Cmpnts  value,
BCFace  face 
)
static

Extracts the velocity component aligned with the supplied boundary-face normal.

Definition at line 67 of file test_boundaries.c.

68{
69 switch (face) {
70 case BC_FACE_NEG_X:
71 case BC_FACE_POS_X:
72 return value.x;
73 case BC_FACE_NEG_Y:
74 case BC_FACE_POS_Y:
75 return value.y;
76 case BC_FACE_NEG_Z:
77 case BC_FACE_POS_Z:
78 return value.z;
79 }
80 return 0.0;
81}
PetscScalar x
Definition variables.h:101
PetscScalar z
Definition variables.h:101
PetscScalar y
Definition variables.h:101
@ 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
Here is the caller graph for this function:

◆ GetFaceOrientationSign()

static PetscReal GetFaceOrientationSign ( BCFace  face)
static

Returns the sign convention used for face-normal flux expectations.

Definition at line 86 of file test_boundaries.c.

87{
88 switch (face) {
89 case BC_FACE_NEG_X:
90 case BC_FACE_NEG_Y:
91 case BC_FACE_NEG_Z:
92 return 1.0;
93 case BC_FACE_POS_X:
94 case BC_FACE_POS_Y:
95 case BC_FACE_POS_Z:
96 return -1.0;
97 }
98 return 0.0;
99}
Here is the caller graph for this function:

◆ GetInletParamKey()

static const char * GetInletParamKey ( BCFace  face)
static

Maps an inlet face to the matching configuration key used by the handler parser.

Definition at line 104 of file test_boundaries.c.

105{
106 switch (face) {
107 case BC_FACE_NEG_X:
108 case BC_FACE_POS_X:
109 return "vx";
110 case BC_FACE_NEG_Y:
111 case BC_FACE_POS_Y:
112 return "vy";
113 case BC_FACE_NEG_Z:
114 case BC_FACE_POS_Z:
115 return "vz";
116 }
117 return "";
118}
Here is the caller graph for this function:

◆ GetFaceInteriorPointCount()

static PetscReal GetFaceInteriorPointCount ( const UserCtx user,
BCFace  face 
)
static

Computes the number of face-interior sample points for a given boundary face.

Definition at line 123 of file test_boundaries.c.

124{
125 switch (face) {
126 case BC_FACE_NEG_X:
127 case BC_FACE_POS_X:
128 return (PetscReal)(user->info.my - 2) * (PetscReal)(user->info.mz - 2);
129 case BC_FACE_NEG_Y:
130 case BC_FACE_POS_Y:
131 return (PetscReal)(user->info.mx - 2) * (PetscReal)(user->info.mz - 2);
132 case BC_FACE_NEG_Z:
133 case BC_FACE_POS_Z:
134 return (PetscReal)(user->info.mx - 2) * (PetscReal)(user->info.my - 2);
135 }
136 return 0.0;
137}
DMDALocalInfo info
Definition variables.h:818
Here is the caller graph for this function:

◆ GetRepresentativeFaceSlots()

static PetscErrorCode GetRepresentativeFaceSlots ( const UserCtx user,
BCFace  face,
PetscInt *  ucont_k,
PetscInt *  ucont_j,
PetscInt *  ucont_i,
PetscInt *  ubcs_k,
PetscInt *  ubcs_j,
PetscInt *  ubcs_i 
)
static

Selects representative Ucont and Ubcs slots for face-matrix assertions.

Definition at line 142 of file test_boundaries.c.

150{
151 const PetscInt sample_i = InteriorSampleIndex(user->info.mx);
152 const PetscInt sample_j = InteriorSampleIndex(user->info.my);
153 const PetscInt sample_k = InteriorSampleIndex(user->info.mz);
154
155 PetscFunctionBeginUser;
156 switch (face) {
157 case BC_FACE_NEG_X:
158 *ucont_k = sample_k; *ucont_j = sample_j; *ucont_i = 0;
159 *ubcs_k = sample_k; *ubcs_j = sample_j; *ubcs_i = 0;
160 break;
161 case BC_FACE_POS_X:
162 *ucont_k = sample_k; *ucont_j = sample_j; *ucont_i = user->info.mx - 2;
163 *ubcs_k = sample_k; *ubcs_j = sample_j; *ubcs_i = user->info.mx - 1;
164 break;
165 case BC_FACE_NEG_Y:
166 *ucont_k = sample_k; *ucont_j = 0; *ucont_i = sample_i;
167 *ubcs_k = sample_k; *ubcs_j = 0; *ubcs_i = sample_i;
168 break;
169 case BC_FACE_POS_Y:
170 *ucont_k = sample_k; *ucont_j = user->info.my - 2; *ucont_i = sample_i;
171 *ubcs_k = sample_k; *ubcs_j = user->info.my - 1; *ubcs_i = sample_i;
172 break;
173 case BC_FACE_NEG_Z:
174 *ucont_k = 0; *ucont_j = sample_j; *ucont_i = sample_i;
175 *ubcs_k = 0; *ubcs_j = sample_j; *ubcs_i = sample_i;
176 break;
177 case BC_FACE_POS_Z:
178 *ucont_k = user->info.mz - 2; *ucont_j = sample_j; *ucont_i = sample_i;
179 *ubcs_k = user->info.mz - 1; *ubcs_j = sample_j; *ubcs_i = sample_i;
180 break;
181 }
182 PetscFunctionReturn(0);
183}
static PetscInt InteriorSampleIndex(PetscInt nodes)
Chooses a stable interior node index for face-sample assertions on tiny test grids.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetParabolicSampleSlots()

static PetscErrorCode GetParabolicSampleSlots ( const UserCtx user,
BCFace  face,
PetscInt *  center_k,
PetscInt *  center_j,
PetscInt *  center_i,
PetscInt *  off_k,
PetscInt *  off_j,
PetscInt *  off_i,
PetscInt *  wall_k,
PetscInt *  wall_j,
PetscInt *  wall_i,
PetscInt *  ubcs_center_k,
PetscInt *  ubcs_center_j,
PetscInt *  ubcs_center_i 
)
static

Selects center, off-center, and wall-adjacent slots for parabolic-face checks.

Definition at line 188 of file test_boundaries.c.

202{
203 PetscFunctionBeginUser;
204 switch (face) {
205 case BC_FACE_NEG_X:
206 *center_k = 3; *center_j = 3; *center_i = 0;
207 *off_k = 2; *off_j = 2; *off_i = 0;
208 *wall_k = 1; *wall_j = 1; *wall_i = 0;
209 *ubcs_center_k = 3; *ubcs_center_j = 3; *ubcs_center_i = 0;
210 break;
211 case BC_FACE_POS_X:
212 *center_k = 3; *center_j = 3; *center_i = user->info.mx - 2;
213 *off_k = 2; *off_j = 2; *off_i = user->info.mx - 2;
214 *wall_k = 1; *wall_j = 1; *wall_i = user->info.mx - 2;
215 *ubcs_center_k = 3; *ubcs_center_j = 3; *ubcs_center_i = user->info.mx - 1;
216 break;
217 case BC_FACE_NEG_Y:
218 *center_k = 3; *center_j = 0; *center_i = 3;
219 *off_k = 2; *off_j = 0; *off_i = 2;
220 *wall_k = 1; *wall_j = 0; *wall_i = 1;
221 *ubcs_center_k = 3; *ubcs_center_j = 0; *ubcs_center_i = 3;
222 break;
223 case BC_FACE_POS_Y:
224 *center_k = 3; *center_j = user->info.my - 2; *center_i = 3;
225 *off_k = 2; *off_j = user->info.my - 2; *off_i = 2;
226 *wall_k = 1; *wall_j = user->info.my - 2; *wall_i = 1;
227 *ubcs_center_k = 3; *ubcs_center_j = user->info.my - 1; *ubcs_center_i = 3;
228 break;
229 case BC_FACE_NEG_Z:
230 *center_k = 0; *center_j = 3; *center_i = 3;
231 *off_k = 0; *off_j = 2; *off_i = 2;
232 *wall_k = 0; *wall_j = 1; *wall_i = 1;
233 *ubcs_center_k = 0; *ubcs_center_j = 3; *ubcs_center_i = 3;
234 break;
235 case BC_FACE_POS_Z:
236 *center_k = user->info.mz - 2; *center_j = 3; *center_i = 3;
237 *off_k = user->info.mz - 2; *off_j = 2; *off_i = 2;
238 *wall_k = user->info.mz - 2; *wall_j = 1; *wall_i = 1;
239 *ubcs_center_k = user->info.mz - 1; *ubcs_center_j = 3; *ubcs_center_i = 3;
240 break;
241 }
242 PetscFunctionReturn(0);
243}
Here is the caller graph for this function:

◆ ResetBoundaryFaceConfig()

static void ResetBoundaryFaceConfig ( BoundaryFaceConfig cfg)
static

Resets one boundary-face configuration entry to a neutral test-local baseline.

Definition at line 248 of file test_boundaries.c.

249{
252 cfg->face_id = BC_FACE_NEG_X;
253 if (cfg->params) {
255 cfg->params = NULL;
256 }
257}
void FreeBC_ParamList(BC_Param *head)
Frees an entire linked list of boundary-condition parameters.
Definition io.c:302
@ INTERFACE
Definition variables.h:253
@ BC_HANDLER_UNDEFINED
Definition variables.h:272
BCHandlerType handler_type
Definition variables.h:337
BC_Param * params
Definition variables.h:338
BCType mathematical_type
Definition variables.h:336
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestCanRankServiceFaceMatchesInletWhenDefined()

static PetscErrorCode TestCanRankServiceFaceMatchesInletWhenDefined ( void  )
static

Tests that face-service detection matches a defined inlet face.

Definition at line 262 of file test_boundaries.c.

263{
264 SimCtx *simCtx = NULL;
265 UserCtx *user = NULL;
266 PetscBool generic_face = PETSC_FALSE;
267 PetscBool inlet_face = PETSC_FALSE;
268 PetscBool any_face_serviceable = PETSC_FALSE;
269
270 PetscFunctionBeginUser;
271 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
272
273 for (BCFace face = BC_FACE_NEG_X; face <= BC_FACE_POS_Z; face++) {
274 PetscCall(CanRankServiceFace(&user->info, user->IM, user->JM, user->KM, face, &generic_face));
275 if (generic_face) {
276 any_face_serviceable = PETSC_TRUE;
277 }
278
279 user->inletFaceDefined = PETSC_TRUE;
280 user->identifiedInletBCFace = face;
281 PetscCall(CanRankServiceInletFace(user, &user->info, user->IM, user->JM, user->KM, &inlet_face));
282 PetscCall(PicurvAssertBool((PetscBool)(inlet_face == generic_face), "inlet face service check should match generic face check"));
283 }
284 PetscCall(PicurvAssertBool(any_face_serviceable, "at least one global face should be serviceable in a non-degenerate domain"));
285
286 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
287 PetscFunctionReturn(0);
288}
PetscErrorCode CanRankServiceInletFace(UserCtx *user, const DMDALocalInfo *info, PetscInt IM_nodes_global, PetscInt JM_nodes_global, PetscInt KM_nodes_global, PetscBool *can_service_inlet_out)
Determines if the current MPI rank owns any part of the globally defined inlet face,...
Definition Boundaries.c:11
PetscErrorCode CanRankServiceFace(const DMDALocalInfo *info, PetscInt IM_nodes_global, PetscInt JM_nodes_global, PetscInt KM_nodes_global, BCFace face_id, PetscBool *can_service_out)
Determines if the current MPI rank owns any part of a specified global face.
Definition Boundaries.c:126
PetscErrorCode PicurvCreateMinimalContexts(SimCtx **simCtx_out, UserCtx **user_out, PetscInt mx, PetscInt my, PetscInt mz)
Builds minimal SimCtx and UserCtx fixtures for C unit tests.
PetscErrorCode PicurvDestroyMinimalContexts(SimCtx **simCtx_ptr, UserCtx **user_ptr)
Destroys minimal SimCtx/UserCtx fixtures and all owned PETSc objects.
PetscErrorCode PicurvAssertBool(PetscBool value, const char *context)
Asserts that one boolean condition is true.
PetscBool inletFaceDefined
Definition variables.h:830
BCFace identifiedInletBCFace
Definition variables.h:831
PetscInt KM
Definition variables.h:820
PetscInt JM
Definition variables.h:820
PetscInt IM
Definition variables.h:820
BCFace
Identifies the six logical faces of a structured computational block.
Definition variables.h:244
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
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestCanRankServiceInletFaceRequiresDefinition()

static PetscErrorCode TestCanRankServiceInletFaceRequiresDefinition ( void  )
static

Tests that inlet-face service requires an inlet face to be defined.

Definition at line 293 of file test_boundaries.c.

294{
295 SimCtx *simCtx = NULL;
296 UserCtx *user = NULL;
297 PetscBool can_service = PETSC_FALSE;
298
299 PetscFunctionBeginUser;
300 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
301
302 user->inletFaceDefined = PETSC_FALSE;
303 for (BCFace face = BC_FACE_NEG_X; face <= BC_FACE_POS_Z; face++) {
304 user->identifiedInletBCFace = face;
305 PetscCall(CanRankServiceInletFace(user, &user->info, user->IM, user->JM, user->KM, &can_service));
306 PetscCall(PicurvAssertBool((PetscBool)!can_service, "undefined inlet face should not be serviceable"));
307 }
308
309 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
310 PetscFunctionReturn(0);
311}
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestBoundaryConditionFactoryAssignments()

static PetscErrorCode TestBoundaryConditionFactoryAssignments ( void  )
static

Tests the boundary-condition factory assignments for representative handlers.

Definition at line 316 of file test_boundaries.c.

317{
318 BoundaryCondition *wall = NULL;
319 BoundaryCondition *inlet = NULL;
320
321 PetscFunctionBeginUser;
323 PetscCall(PicurvAssertBool((PetscBool)(wall != NULL), "factory should return wall handler"));
324 PetscCall(PicurvAssertIntEqual(BC_PRIORITY_WALL, wall->priority, "wall handler priority"));
325 PetscCall(PicurvAssertBool((PetscBool)(wall->Apply != NULL), "wall handler should expose Apply"));
326 PetscCall(PicurvAssertBool((PetscBool)(wall->Destroy == NULL), "wall handler should not require destroy hook"));
327 PetscCall(DestroyBoundaryHandler(&wall));
328
330 PetscCall(PicurvAssertBool((PetscBool)(inlet != NULL), "factory should return inlet handler"));
331 PetscCall(PicurvAssertIntEqual(BC_PRIORITY_INLET, inlet->priority, "inlet handler priority"));
332 PetscCall(PicurvAssertBool((PetscBool)(inlet->Initialize != NULL), "inlet handler should expose Initialize"));
333 PetscCall(PicurvAssertBool((PetscBool)(inlet->Apply != NULL), "inlet handler should expose Apply"));
334 PetscCall(PicurvAssertBool((PetscBool)(inlet->Destroy != NULL), "inlet handler should expose its destroy hook"));
335 PetscCall(DestroyBoundaryHandler(&inlet));
336
337 PetscFunctionReturn(0);
338}
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:325
PetscErrorCode(* Apply)(BoundaryCondition *self, BCContext *ctx)
Definition variables.h:327
BCPriorityType priority
Definition variables.h:323
static PetscErrorCode DestroyBoundaryHandler(BoundaryCondition **bc_ptr)
Destroys one boundary-condition handler allocated by a boundary test.
PetscErrorCode PicurvAssertIntEqual(PetscInt expected, PetscInt actual, const char *context)
Asserts that two integer values are equal.
@ BC_HANDLER_INLET_CONSTANT_VELOCITY
Definition variables.h:276
@ BC_HANDLER_WALL_NOSLIP
Definition variables.h:273
@ BC_PRIORITY_WALL
Definition variables.h:295
@ BC_PRIORITY_INLET
Definition variables.h:293
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestBoundaryConditionFactoryImplementedHandlerMatrix()

static PetscErrorCode TestBoundaryConditionFactoryImplementedHandlerMatrix ( void  )
static

Tests the implemented-handler matrix exposed by the factory.

Definition at line 343 of file test_boundaries.c.

344{
345 struct HandlerExpectation {
346 BCHandlerType handler;
347 BCPriorityType priority;
348 PetscBool expect_apply;
349 PetscBool expect_initialize;
350 };
351 const struct HandlerExpectation expectations[] = {
352 {BC_HANDLER_WALL_NOSLIP, BC_PRIORITY_WALL, PETSC_TRUE, PETSC_FALSE},
353 {BC_HANDLER_INLET_CONSTANT_VELOCITY, BC_PRIORITY_INLET, PETSC_TRUE, PETSC_TRUE},
354 {BC_HANDLER_INLET_PARABOLIC, BC_PRIORITY_INLET, PETSC_TRUE, PETSC_TRUE},
355 {BC_HANDLER_OUTLET_CONSERVATION, BC_PRIORITY_OUTLET, PETSC_TRUE, PETSC_FALSE},
356 };
357
358 PetscFunctionBeginUser;
359 for (size_t i = 0; i < sizeof(expectations) / sizeof(expectations[0]); ++i) {
360 BoundaryCondition *bc = NULL;
361 PetscCall(BoundaryCondition_Create(expectations[i].handler, &bc));
362 PetscCall(PicurvAssertBool((PetscBool)(bc != NULL), "factory should allocate a handler object"));
363 PetscCall(PicurvAssertIntEqual(expectations[i].priority, bc->priority, "handler priority should match expectation"));
364 PetscCall(PicurvAssertBool((PetscBool)((bc->Apply != NULL) == expectations[i].expect_apply), "Apply hook expectation mismatch"));
365 PetscCall(PicurvAssertBool((PetscBool)((bc->Initialize != NULL) == expectations[i].expect_initialize), "Initialize hook expectation mismatch"));
366 PetscCall(DestroyBoundaryHandler(&bc));
367 }
368 PetscFunctionReturn(0);
369}
BCHandlerType
Defines the specific computational "strategy" for a boundary handler.
Definition variables.h:271
@ BC_HANDLER_INLET_PARABOLIC
Definition variables.h:277
@ BC_HANDLER_OUTLET_CONSERVATION
Definition variables.h:282
BCPriorityType
Definition variables.h:291
@ BC_PRIORITY_OUTLET
Definition variables.h:296
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestBoundaryConditionFactoryRejectsUnsupportedHandler()

static PetscErrorCode TestBoundaryConditionFactoryRejectsUnsupportedHandler ( void  )
static

Tests that unsupported handler types are rejected by the factory.

Definition at line 374 of file test_boundaries.c.

375{
376 BoundaryCondition *bc = NULL;
377 PetscErrorCode ierr_create;
378
379 PetscFunctionBeginUser;
380 PetscCall(PetscPushErrorHandler(PetscIgnoreErrorHandler, NULL));
382 PetscCall(PetscPopErrorHandler());
383
384 PetscCall(PicurvAssertBool((PetscBool)(ierr_create != 0), "unsupported handler should return a non-zero error code"));
385 if (bc) {
386 PetscCall(PetscFree(bc));
387 }
388 PetscFunctionReturn(0);
389}
@ BC_HANDLER_OUTLET_PRESSURE
Definition variables.h:283
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestGetDeterministicFaceGridLocationFaceMatrix()

static PetscErrorCode TestGetDeterministicFaceGridLocationFaceMatrix ( void  )
static

Tests deterministic inlet-face grid-location helpers across all faces.

Definition at line 394 of file test_boundaries.c.

395{
396 SimCtx *simCtx = NULL;
397 UserCtx *user = NULL;
398
399 PetscFunctionBeginUser;
400 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 8, 8, 8));
401 simCtx->np = 128;
402
403 for (BCFace face = BC_FACE_NEG_X; face <= BC_FACE_POS_Z; ++face) {
404 PetscInt ci = -1, cj = -1, ck = -1;
405 PetscReal xi = -1.0, eta = -1.0, zta = -1.0;
406 PetscBool placed = PETSC_FALSE;
407
408 user->identifiedInletBCFace = face;
410 user, &user->info,
411 user->info.xs, user->info.ys, user->info.zs,
412 user->info.mx, user->info.my, user->info.mz,
413 0,
414 &ci, &cj, &ck, &xi, &eta, &zta, &placed));
415
416 PetscCall(PicurvAssertBool(placed, "single-rank deterministic face placement should succeed"));
417 PetscCall(PicurvAssertBool((PetscBool)(ci >= user->info.xs && ci < user->info.xs + user->info.xm), "ci must map to owned node window"));
418 PetscCall(PicurvAssertBool((PetscBool)(cj >= user->info.ys && cj < user->info.ys + user->info.ym), "cj must map to owned node window"));
419 PetscCall(PicurvAssertBool((PetscBool)(ck >= user->info.zs && ck < user->info.zs + user->info.zm), "ck must map to owned node window"));
420 PetscCall(PicurvAssertBool((PetscBool)(xi >= 0.0 && xi < 1.0), "xi should be in [0,1)"));
421 PetscCall(PicurvAssertBool((PetscBool)(eta >= 0.0 && eta < 1.0), "eta should be in [0,1)"));
422 PetscCall(PicurvAssertBool((PetscBool)(zta >= 0.0 && zta < 1.0), "zta should be in [0,1)"));
423
424 if (face == BC_FACE_NEG_X || face == BC_FACE_POS_X) {
425 PetscCall(PicurvAssertRealNear(0.5, xi, 1.0e-10, "deterministic x-face placement should sit halfway into boundary-adjacent cell"));
426 } else if (face == BC_FACE_NEG_Y || face == BC_FACE_POS_Y) {
427 PetscCall(PicurvAssertRealNear(0.5, eta, 1.0e-10, "deterministic y-face placement should sit halfway into boundary-adjacent cell"));
428 } else {
429 PetscCall(PicurvAssertRealNear(0.5, zta, 1.0e-10, "deterministic z-face placement should sit halfway into boundary-adjacent cell"));
430 }
431 }
432
433 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
434 PetscFunctionReturn(0);
435}
PetscErrorCode GetDeterministicFaceGridLocation(UserCtx *user, const DMDALocalInfo *info, PetscInt xs_gnode_rank, PetscInt ys_gnode_rank, PetscInt zs_gnode_rank, PetscInt IM_cells_global, PetscInt JM_cells_global, PetscInt KM_cells_global, PetscInt64 particle_global_id, PetscInt *ci_metric_lnode_out, PetscInt *cj_metric_lnode_out, PetscInt *ck_metric_lnode_out, PetscReal *xi_metric_logic_out, PetscReal *eta_metric_logic_out, PetscReal *zta_metric_logic_out, PetscBool *placement_successful_out)
Places particles in a deterministic grid/raster pattern on a specified domain face.
Definition Boundaries.c:212
PetscErrorCode PicurvAssertRealNear(PetscReal expected, PetscReal actual, PetscReal tol, const char *context)
Asserts that two real values agree within tolerance.
PetscInt np
Definition variables.h:739
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestGetRandomCellAndLogicalCoordsOnInletFaceMatrix()

static PetscErrorCode TestGetRandomCellAndLogicalCoordsOnInletFaceMatrix ( void  )
static

Tests random inlet-face cell selection across all supported faces.

Definition at line 440 of file test_boundaries.c.

441{
442 SimCtx *simCtx = NULL;
443 UserCtx *user = NULL;
444 PetscRandom rand_i = NULL, rand_j = NULL, rand_k = NULL;
445
446 PetscFunctionBeginUser;
447 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 8, 8, 8));
448
449 PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rand_i));
450 PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rand_j));
451 PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rand_k));
452 PetscCall(PetscRandomSetInterval(rand_i, 0.0, 1.0));
453 PetscCall(PetscRandomSetInterval(rand_j, 0.0, 1.0));
454 PetscCall(PetscRandomSetInterval(rand_k, 0.0, 1.0));
455 PetscCall(PetscRandomSetFromOptions(rand_i));
456 PetscCall(PetscRandomSetFromOptions(rand_j));
457 PetscCall(PetscRandomSetFromOptions(rand_k));
458
459 for (BCFace face = BC_FACE_NEG_X; face <= BC_FACE_POS_Z; ++face) {
460 PetscInt ci = -1, cj = -1, ck = -1;
461 PetscReal xi = -1.0, eta = -1.0, zta = -1.0;
462 const PetscReal boundary_eps = 5.0e-4;
463
464 user->identifiedInletBCFace = face;
466 user, &user->info,
467 user->info.xs, user->info.ys, user->info.zs,
468 user->info.mx, user->info.my, user->info.mz,
469 &rand_i, &rand_j, &rand_k,
470 &ci, &cj, &ck,
471 &xi, &eta, &zta));
472
473 PetscCall(PicurvAssertBool((PetscBool)(ci >= user->info.xs && ci < user->info.xs + user->info.xm), "ci must map to owned node window"));
474 PetscCall(PicurvAssertBool((PetscBool)(cj >= user->info.ys && cj < user->info.ys + user->info.ym), "cj must map to owned node window"));
475 PetscCall(PicurvAssertBool((PetscBool)(ck >= user->info.zs && ck < user->info.zs + user->info.zm), "ck must map to owned node window"));
476 PetscCall(PicurvAssertBool((PetscBool)(xi >= 0.0 && xi <= 1.0), "xi should be in [0,1]"));
477 PetscCall(PicurvAssertBool((PetscBool)(eta >= 0.0 && eta <= 1.0), "eta should be in [0,1]"));
478 PetscCall(PicurvAssertBool((PetscBool)(zta >= 0.0 && zta <= 1.0), "zta should be in [0,1]"));
479
480 switch (face) {
481 case BC_FACE_NEG_X:
482 PetscCall(PicurvAssertBool((PetscBool)(xi <= boundary_eps), "NEG_X should pin xi near 0"));
483 break;
484 case BC_FACE_POS_X:
485 PetscCall(PicurvAssertBool((PetscBool)(xi >= 1.0 - boundary_eps), "POS_X should pin xi near 1"));
486 break;
487 case BC_FACE_NEG_Y:
488 PetscCall(PicurvAssertBool((PetscBool)(eta <= boundary_eps), "NEG_Y should pin eta near 0"));
489 break;
490 case BC_FACE_POS_Y:
491 PetscCall(PicurvAssertBool((PetscBool)(eta >= 1.0 - boundary_eps), "POS_Y should pin eta near 1"));
492 break;
493 case BC_FACE_NEG_Z:
494 PetscCall(PicurvAssertBool((PetscBool)(zta <= boundary_eps), "NEG_Z should pin zta near 0"));
495 break;
496 case BC_FACE_POS_Z:
497 PetscCall(PicurvAssertBool((PetscBool)(zta >= 1.0 - boundary_eps), "POS_Z should pin zta near 1"));
498 break;
499 }
500 }
501
502 PetscCall(PetscRandomDestroy(&rand_i));
503 PetscCall(PetscRandomDestroy(&rand_j));
504 PetscCall(PetscRandomDestroy(&rand_k));
505 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
506 PetscFunctionReturn(0);
507}
PetscErrorCode GetRandomCellAndLogicalCoordsOnInletFace(UserCtx *user, const DMDALocalInfo *info, PetscInt xs_gnode_rank, PetscInt ys_gnode_rank, PetscInt zs_gnode_rank, PetscInt IM_nodes_global, PetscInt JM_nodes_global, PetscInt KM_nodes_global, PetscRandom *rand_logic_i_ptr, PetscRandom *rand_logic_j_ptr, PetscRandom *rand_logic_k_ptr, PetscInt *ci_metric_lnode_out, PetscInt *cj_metric_lnode_out, PetscInt *ck_metric_lnode_out, PetscReal *xi_metric_logic_out, PetscReal *eta_metric_logic_out, PetscReal *zta_metric_logic_out)
Assuming the current rank services the inlet face, this function selects a random cell (owned by this...
Definition Boundaries.c:399
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestInletConstantVelocityHandlerBehavior()

static PetscErrorCode TestInletConstantVelocityHandlerBehavior ( void  )
static

Tests constant inlet handler initialization and face writes on a Z inlet.

Definition at line 512 of file test_boundaries.c.

513{
514 SimCtx *simCtx = NULL;
515 UserCtx *user = NULL;
516 BoundaryCondition *bc = NULL;
517 BCContext ctx;
518 Cmpnts ***ucont = NULL;
519 Cmpnts ***ubcs = NULL;
520 PetscReal local_inflow = 0.0;
521 PetscReal local_outflow = 0.0;
522 PetscReal expected_flux = 0.0;
523
524 PetscFunctionBeginUser;
525 PetscCall(PetscMemzero(&ctx, sizeof(ctx)));
526 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 6, 6, 6));
527 PetscCall(VecSet(user->Ucont, 0.0));
528 PetscCall(VecSet(user->Bcs.Ubcs, 0.0));
529 PetscCall(PicurvPopulateIdentityMetrics(user));
530
534 PetscCall(AppendBCParam(&user->boundary_faces[BC_FACE_NEG_Z].params, "vz", "2.5"));
535
536 ctx.user = user;
539 PetscCall(bc->Initialize(bc, &ctx));
540 PetscCall(bc->PostStep(bc, &ctx, &local_inflow, &local_outflow));
541
542 PetscCall(DMDAVecGetArrayRead(user->fda, user->Ucont, &ucont));
543 PetscCall(DMDAVecGetArrayRead(user->fda, user->Bcs.Ubcs, &ubcs));
544 PetscCall(PicurvAssertRealNear(2.5, ucont[0][3][3].z, 1.0e-12, "constant inlet should set Ucont normal component"));
545 PetscCall(PicurvAssertRealNear(2.5, ubcs[0][3][3].z, 1.0e-12, "constant inlet should set Ubcs normal component"));
546 PetscCall(DMDAVecRestoreArrayRead(user->fda, user->Ucont, &ucont));
547 PetscCall(DMDAVecRestoreArrayRead(user->fda, user->Bcs.Ubcs, &ubcs));
548
549 expected_flux = (PetscReal)(user->info.mx - 2) * (PetscReal)(user->info.my - 2) * 2.5;
550 PetscCall(PicurvAssertRealNear(expected_flux, local_inflow, 1.0e-12, "constant inlet PostStep should sum the face flux"));
551 PetscCall(PicurvAssertRealNear(0.0, local_outflow, 1.0e-12, "constant inlet should not contribute to outflow accumulation"));
552
553 PetscCall(DestroyBoundaryHandler(&bc));
555 user->boundary_faces[BC_FACE_NEG_Z].params = NULL;
556 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
557 PetscFunctionReturn(0);
558}
PetscErrorCode(* PostStep)(BoundaryCondition *self, BCContext *ctx, PetscReal *local_inflow, PetscReal *local_outflow)
Definition variables.h:328
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.
@ INLET
Definition variables.h:258
BoundaryFaceConfig boundary_faces[6]
Definition variables.h:829
BCFace face_id
Definition variables.h:313
Vec Ucont
Definition variables.h:837
Vec Ubcs
Physical Cartesian velocity at boundary faces. Full 3D array but only boundary-face entries are meani...
Definition variables.h:121
BCS Bcs
Definition variables.h:832
UserCtx * user
Definition variables.h:312
Provides execution context for a boundary condition handler.
Definition variables.h:311
A 3D point or vector with PetscScalar components.
Definition variables.h:100
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestWallNoSlipHandlerFaceMatrix()

static PetscErrorCode TestWallNoSlipHandlerFaceMatrix ( void  )
static

Tests wall no-slip application across the full face matrix.

Definition at line 563 of file test_boundaries.c.

564{
565 SimCtx *simCtx = NULL;
566 UserCtx *user = NULL;
567 BoundaryCondition *bc = NULL;
568 BCContext ctx;
569 Cmpnts ***ucont = NULL;
570 Cmpnts ***ubcs = NULL;
571 const BCFace faces[] = {
575 };
576
577 PetscFunctionBeginUser;
578 PetscCall(PetscMemzero(&ctx, sizeof(ctx)));
579 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 6, 6, 6));
581
582 for (size_t n = 0; n < sizeof(faces) / sizeof(faces[0]); ++n) {
583 const BCFace face = faces[n];
584 PetscInt ucont_k, ucont_j, ucont_i, ubcs_k, ubcs_j, ubcs_i;
585
586 PetscCall(VecSet(user->Ucont, 5.0));
587 PetscCall(VecSet(user->Bcs.Ubcs, 7.0));
588 ctx.user = user;
589 ctx.face_id = face;
590 PetscCall(bc->Apply(bc, &ctx));
591
592 PetscCall(GetRepresentativeFaceSlots(user, face, &ucont_k, &ucont_j, &ucont_i, &ubcs_k, &ubcs_j, &ubcs_i));
593 PetscCall(DMDAVecGetArrayRead(user->fda, user->Ucont, &ucont));
594 PetscCall(DMDAVecGetArrayRead(user->fda, user->Bcs.Ubcs, &ubcs));
595 PetscCall(PicurvAssertRealNear(0.0, GetFaceNormalComponent(ucont[ucont_k][ucont_j][ucont_i], face), 1.0e-12, "wall no-slip should zero the normal flux component"));
596 PetscCall(PicurvAssertRealNear(0.0, ubcs[ubcs_k][ubcs_j][ubcs_i].x, 1.0e-12, "wall no-slip should zero Ubcs.x"));
597 PetscCall(PicurvAssertRealNear(0.0, ubcs[ubcs_k][ubcs_j][ubcs_i].y, 1.0e-12, "wall no-slip should zero Ubcs.y"));
598 PetscCall(PicurvAssertRealNear(0.0, ubcs[ubcs_k][ubcs_j][ubcs_i].z, 1.0e-12, "wall no-slip should zero Ubcs.z"));
599 PetscCall(DMDAVecRestoreArrayRead(user->fda, user->Ucont, &ucont));
600 PetscCall(DMDAVecRestoreArrayRead(user->fda, user->Bcs.Ubcs, &ubcs));
601 }
602
603 PetscCall(DestroyBoundaryHandler(&bc));
604 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
605 PetscFunctionReturn(0);
606}
static PetscReal GetFaceNormalComponent(Cmpnts value, BCFace face)
Extracts the velocity component aligned with the supplied boundary-face normal.
static PetscErrorCode GetRepresentativeFaceSlots(const UserCtx *user, BCFace face, PetscInt *ucont_k, PetscInt *ucont_j, PetscInt *ucont_i, PetscInt *ubcs_k, PetscInt *ubcs_j, PetscInt *ubcs_i)
Selects representative Ucont and Ubcs slots for face-matrix assertions.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestInletConstantVelocityHandlerFaceMatrix()

static PetscErrorCode TestInletConstantVelocityHandlerFaceMatrix ( void  )
static

Tests constant-inlet initialization, flux accounting, and face writes across all faces.

Definition at line 611 of file test_boundaries.c.

612{
613 struct ConstantInletCase {
614 BCFace face;
615 PetscReal velocity;
616 const char *value_text;
617 };
618 const struct ConstantInletCase cases[] = {
619 {BC_FACE_NEG_X, 1.5, "1.5"},
620 {BC_FACE_POS_X, 2.0, "2.0"},
621 {BC_FACE_NEG_Y, 2.5, "2.5"},
622 {BC_FACE_POS_Y, 3.0, "3.0"},
623 {BC_FACE_NEG_Z, 3.5, "3.5"},
624 {BC_FACE_POS_Z, 4.0, "4.0"},
625 };
626 SimCtx *simCtx = NULL;
627 UserCtx *user = NULL;
628 BCContext ctx;
629 Cmpnts ***ucont = NULL;
630 Cmpnts ***ubcs = NULL;
631
632 PetscFunctionBeginUser;
633 PetscCall(PetscMemzero(&ctx, sizeof(ctx)));
634 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 6, 6, 6));
635 PetscCall(PicurvPopulateIdentityMetrics(user));
636
637 for (size_t n = 0; n < sizeof(cases) / sizeof(cases[0]); ++n) {
638 const struct ConstantInletCase *test_case = &cases[n];
639 BoundaryCondition *bc = NULL;
640 PetscReal local_inflow = 0.0;
641 PetscReal local_outflow = 0.0;
642 PetscInt ucont_k, ucont_j, ucont_i, ubcs_k, ubcs_j, ubcs_i;
643 const PetscReal expected_value = GetFaceOrientationSign(test_case->face) * test_case->velocity;
644 const PetscReal expected_flux = GetFaceInteriorPointCount(user, test_case->face) * expected_value;
645
646 PetscCall(VecSet(user->Ucont, 0.0));
647 PetscCall(VecSet(user->Bcs.Ubcs, 0.0));
648 ResetBoundaryFaceConfig(&user->boundary_faces[test_case->face]);
649 user->boundary_faces[test_case->face].face_id = test_case->face;
650 user->boundary_faces[test_case->face].mathematical_type = INLET;
652 PetscCall(AppendBCParam(&user->boundary_faces[test_case->face].params,
653 GetInletParamKey(test_case->face),
654 test_case->value_text));
655
656 ctx.user = user;
657 ctx.face_id = test_case->face;
659 PetscCall(bc->PreStep(bc, &ctx, &local_inflow, &local_outflow));
660 PetscCall(PicurvAssertRealNear(0.0, local_inflow, 1.0e-12, "constant inlet PreStep should leave inflow unchanged"));
661 PetscCall(PicurvAssertRealNear(0.0, local_outflow, 1.0e-12, "constant inlet PreStep should leave outflow unchanged"));
662 PetscCall(bc->Initialize(bc, &ctx));
663 PetscCall(bc->PostStep(bc, &ctx, &local_inflow, &local_outflow));
664
665 PetscCall(GetRepresentativeFaceSlots(user, test_case->face, &ucont_k, &ucont_j, &ucont_i, &ubcs_k, &ubcs_j, &ubcs_i));
666 PetscCall(DMDAVecGetArrayRead(user->fda, user->Ucont, &ucont));
667 PetscCall(DMDAVecGetArrayRead(user->fda, user->Bcs.Ubcs, &ubcs));
668 PetscCall(PicurvAssertRealNear(expected_value, GetFaceNormalComponent(ucont[ucont_k][ucont_j][ucont_i], test_case->face), 1.0e-12, "constant inlet should set the face-normal Ucont component"));
669 PetscCall(PicurvAssertRealNear(expected_value, GetFaceNormalComponent(ubcs[ubcs_k][ubcs_j][ubcs_i], test_case->face), 1.0e-12, "constant inlet should set the face-normal Ubcs component"));
670 PetscCall(PicurvAssertRealNear((test_case->face == BC_FACE_NEG_X || test_case->face == BC_FACE_POS_X) ? expected_value : 0.0, ubcs[ubcs_k][ubcs_j][ubcs_i].x, 1.0e-12, "constant inlet should only set the expected Ubcs axis"));
671 PetscCall(PicurvAssertRealNear((test_case->face == BC_FACE_NEG_Y || test_case->face == BC_FACE_POS_Y) ? expected_value : 0.0, ubcs[ubcs_k][ubcs_j][ubcs_i].y, 1.0e-12, "constant inlet should only set the expected Ubcs axis"));
672 PetscCall(PicurvAssertRealNear((test_case->face == BC_FACE_NEG_Z || test_case->face == BC_FACE_POS_Z) ? expected_value : 0.0, ubcs[ubcs_k][ubcs_j][ubcs_i].z, 1.0e-12, "constant inlet should only set the expected Ubcs axis"));
673 PetscCall(DMDAVecRestoreArrayRead(user->fda, user->Ucont, &ucont));
674 PetscCall(DMDAVecRestoreArrayRead(user->fda, user->Bcs.Ubcs, &ubcs));
675 PetscCall(PicurvAssertRealNear(expected_flux, local_inflow, 1.0e-12, "constant inlet PostStep should integrate the face flux with orientation"));
676 PetscCall(PicurvAssertRealNear(0.0, local_outflow, 1.0e-12, "constant inlet should not add to outflow"));
677
678 PetscCall(DestroyBoundaryHandler(&bc));
679 ResetBoundaryFaceConfig(&user->boundary_faces[test_case->face]);
680 }
681
682 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
683 PetscFunctionReturn(0);
684}
PetscErrorCode(* PreStep)(BoundaryCondition *self, BCContext *ctx, PetscReal *local_inflow, PetscReal *local_outflow)
Definition variables.h:326
static PetscReal GetFaceInteriorPointCount(const UserCtx *user, BCFace face)
Computes the number of face-interior sample points for a given boundary face.
static PetscReal GetFaceOrientationSign(BCFace face)
Returns the sign convention used for face-normal flux expectations.
static const char * GetInletParamKey(BCFace face)
Maps an inlet face to the matching configuration key used by the handler parser.
static void ResetBoundaryFaceConfig(BoundaryFaceConfig *cfg)
Resets one boundary-face configuration entry to a neutral test-local baseline.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestInletParabolicProfileHandlerBehavior()

static PetscErrorCode TestInletParabolicProfileHandlerBehavior ( void  )
static

Tests parabolic inlet handler shape on a tiny Z-face.

Definition at line 689 of file test_boundaries.c.

690{
691 SimCtx *simCtx = NULL;
692 UserCtx *user = NULL;
693 BoundaryCondition *bc = NULL;
694 BCContext ctx;
695 Cmpnts ***ucont = NULL;
696
697 PetscFunctionBeginUser;
698 PetscCall(PetscMemzero(&ctx, sizeof(ctx)));
699 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 6, 6, 6));
700 PetscCall(PicurvPopulateIdentityMetrics(user));
701
705 PetscCall(AppendBCParam(&user->boundary_faces[BC_FACE_NEG_Z].params, "v_max", "4.0"));
706
707 ctx.user = user;
710 PetscCall(bc->Initialize(bc, &ctx));
711
712 PetscCall(DMDAVecGetArrayRead(user->fda, user->Ucont, &ucont));
713 PetscCall(PicurvAssertRealNear(4.0, ucont[0][3][3].z, 1.0e-12, "parabolic inlet should peak at the face centerline"));
714 PetscCall(PicurvAssertRealNear(0.0, ucont[0][1][1].z, 1.0e-12, "parabolic inlet should vanish at wall-adjacent locations"));
715 PetscCall(PicurvAssertBool((PetscBool)(ucont[0][3][3].z > ucont[0][2][2].z), "parabolic inlet centerline should exceed off-center velocity"));
716 PetscCall(DMDAVecRestoreArrayRead(user->fda, user->Ucont, &ucont));
717
718 PetscCall(DestroyBoundaryHandler(&bc));
720 user->boundary_faces[BC_FACE_NEG_Z].params = NULL;
721 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
722 PetscFunctionReturn(0);
723}
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestInletParabolicProfileHandlerFaceMatrix()

static PetscErrorCode TestInletParabolicProfileHandlerFaceMatrix ( void  )
static

Tests parabolic-inlet centerline, wall, and flux behavior across all faces.

Definition at line 728 of file test_boundaries.c.

729{
730 const BCFace faces[] = {
734 };
735 SimCtx *simCtx = NULL;
736 UserCtx *user = NULL;
737 BCContext ctx;
738 Cmpnts ***ucont = NULL;
739 Cmpnts ***ubcs = NULL;
740
741 PetscFunctionBeginUser;
742 PetscCall(PetscMemzero(&ctx, sizeof(ctx)));
743 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 6, 6, 6));
744 PetscCall(PicurvPopulateIdentityMetrics(user));
745
746 for (size_t n = 0; n < sizeof(faces) / sizeof(faces[0]); ++n) {
747 const BCFace face = faces[n];
748 BoundaryCondition *bc = NULL;
749 PetscReal local_inflow = 0.0;
750 PetscReal local_outflow = 0.0;
751 PetscInt center_k, center_j, center_i;
752 PetscInt off_k, off_j, off_i;
753 PetscInt wall_k, wall_j, wall_i;
754 PetscInt ubcs_center_k, ubcs_center_j, ubcs_center_i;
755 const PetscReal expected_center = GetFaceOrientationSign(face) * 4.0;
756
757 PetscCall(VecSet(user->Ucont, 0.0));
758 PetscCall(VecSet(user->Bcs.Ubcs, 0.0));
760 user->boundary_faces[face].face_id = face;
763 PetscCall(AppendBCParam(&user->boundary_faces[face].params, "v_max", "4.0"));
764
765 ctx.user = user;
766 ctx.face_id = face;
768 PetscCall(bc->PreStep(bc, &ctx, &local_inflow, &local_outflow));
769 PetscCall(PicurvAssertRealNear(0.0, local_inflow, 1.0e-12, "parabolic inlet PreStep should leave inflow unchanged"));
770 PetscCall(PicurvAssertRealNear(0.0, local_outflow, 1.0e-12, "parabolic inlet PreStep should leave outflow unchanged"));
771 PetscCall(bc->Initialize(bc, &ctx));
772 PetscCall(bc->PostStep(bc, &ctx, &local_inflow, &local_outflow));
773
774 PetscCall(GetParabolicSampleSlots(user, face,
775 &center_k, &center_j, &center_i,
776 &off_k, &off_j, &off_i,
777 &wall_k, &wall_j, &wall_i,
778 &ubcs_center_k, &ubcs_center_j, &ubcs_center_i));
779 PetscCall(DMDAVecGetArrayRead(user->fda, user->Ucont, &ucont));
780 PetscCall(DMDAVecGetArrayRead(user->fda, user->Bcs.Ubcs, &ubcs));
781 PetscCall(PicurvAssertRealNear(expected_center, GetFaceNormalComponent(ucont[center_k][center_j][center_i], face), 1.0e-12, "parabolic inlet should peak at the face centerline"));
782 PetscCall(PicurvAssertRealNear(expected_center, GetFaceNormalComponent(ubcs[ubcs_center_k][ubcs_center_j][ubcs_center_i], face), 1.0e-12, "parabolic inlet should write the centerline boundary velocity"));
783 PetscCall(PicurvAssertRealNear(0.0, GetFaceNormalComponent(ucont[wall_k][wall_j][wall_i], face), 1.0e-12, "parabolic inlet should vanish at the wall"));
784 PetscCall(PicurvAssertBool((PetscBool)(PetscAbsReal(GetFaceNormalComponent(ucont[center_k][center_j][center_i], face)) >
785 PetscAbsReal(GetFaceNormalComponent(ucont[off_k][off_j][off_i], face))),
786 "parabolic inlet centerline magnitude should exceed the off-center magnitude"));
787 PetscCall(DMDAVecRestoreArrayRead(user->fda, user->Ucont, &ucont));
788 PetscCall(DMDAVecRestoreArrayRead(user->fda, user->Bcs.Ubcs, &ubcs));
789 PetscCall(PicurvAssertBool((PetscBool)(PetscAbsReal(local_inflow) > 0.0), "parabolic inlet PostStep should accumulate a non-zero flux"));
790 PetscCall(PicurvAssertBool((PetscBool)((local_inflow > 0.0) == (GetFaceOrientationSign(face) > 0.0)),
791 "parabolic inlet PostStep should preserve the face orientation"));
792 PetscCall(PicurvAssertRealNear(0.0, local_outflow, 1.0e-12, "parabolic inlet should not add to outflow"));
793
794 PetscCall(DestroyBoundaryHandler(&bc));
796 }
797
798 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
799 PetscFunctionReturn(0);
800}
static PetscErrorCode GetParabolicSampleSlots(const UserCtx *user, BCFace face, PetscInt *center_k, PetscInt *center_j, PetscInt *center_i, PetscInt *off_k, PetscInt *off_j, PetscInt *off_i, PetscInt *wall_k, PetscInt *wall_j, PetscInt *wall_i, PetscInt *ubcs_center_k, PetscInt *ubcs_center_j, PetscInt *ubcs_center_i)
Selects center, off-center, and wall-adjacent slots for parabolic-face checks.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestOutletConservationHandlerBehavior()

static PetscErrorCode TestOutletConservationHandlerBehavior ( void  )
static

Tests outlet conservation handler correction and post-step flux accounting.

Definition at line 805 of file test_boundaries.c.

806{
807 SimCtx *simCtx = NULL;
808 UserCtx *user = NULL;
809 BoundaryCondition *bc = NULL;
810 BCContext ctx;
811 PetscReal global_inflow = 30.0;
812 PetscReal global_farfield_in = 0.0;
813 PetscReal global_farfield_out = 0.0;
814 PetscReal global_outflow = 25.0;
815 PetscReal local_inflow = 0.0;
816 PetscReal local_outflow = 0.0;
817 Cmpnts ***ucont = NULL;
818 Cmpnts ***ucat = NULL;
819
820 PetscFunctionBeginUser;
821 PetscCall(PetscMemzero(&ctx, sizeof(ctx)));
822 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 6, 6, 6));
823 PetscCall(PicurvPopulateIdentityMetrics(user));
824 PetscCall(VecSet(user->Ucat, 1.0));
825 PetscCall(DMGlobalToLocalBegin(user->fda, user->Ucat, INSERT_VALUES, user->lUcat));
826 PetscCall(DMGlobalToLocalEnd(user->fda, user->Ucat, INSERT_VALUES, user->lUcat));
827 simCtx->AreaOutSum = (PetscReal)(user->info.mx - 2) * (PetscReal)(user->info.my - 2);
828
832
833 ctx.user = user;
835 ctx.global_inflow_sum = &global_inflow;
836 ctx.global_farfield_inflow_sum = &global_farfield_in;
837 ctx.global_farfield_outflow_sum = &global_farfield_out;
838 ctx.global_outflow_sum = &global_outflow;
839
841 PetscCall(bc->Apply(bc, &ctx));
842 PetscCall(bc->PostStep(bc, &ctx, &local_inflow, &local_outflow));
843
844 PetscCall(DMDAVecGetArrayRead(user->fda, user->Ucont, &ucont));
845 PetscCall(DMDAVecGetArrayRead(user->fda, user->lUcat, &ucat));
846 PetscCall(PicurvAssertRealNear(1.2, ucont[5][3][3].z, 1.0e-12, "outlet correction should add the expected flux trim"));
847 PetscCall(PicurvAssertRealNear(1.0, ucat[5][3][3].z, 1.0e-12, "outlet handler should preserve the interior Ucat reference"));
848 PetscCall(DMDAVecRestoreArrayRead(user->fda, user->Ucont, &ucont));
849 PetscCall(DMDAVecRestoreArrayRead(user->fda, user->lUcat, &ucat));
850
851 PetscCall(PicurvAssertRealNear(30.0, local_outflow, 1.0e-12, "outlet PostStep should report corrected flux"));
852 PetscCall(PicurvAssertRealNear(0.0, local_inflow, 1.0e-12, "outlet handler should not accumulate inflow"));
853
854 PetscCall(DestroyBoundaryHandler(&bc));
855 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
856 PetscFunctionReturn(0);
857}
const PetscReal * global_outflow_sum
Definition variables.h:317
@ OUTLET
Definition variables.h:257
const PetscReal * global_farfield_outflow_sum
Definition variables.h:316
const PetscReal * global_inflow_sum
Definition variables.h:314
Vec Ucat
Definition variables.h:837
const PetscReal * global_farfield_inflow_sum
Definition variables.h:315
PetscReal AreaOutSum
Definition variables.h:726
Vec lUcat
Definition variables.h:837
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestOutletConservationHandlerFaceMatrix()

static PetscErrorCode TestOutletConservationHandlerFaceMatrix ( void  )
static

Tests outlet-conservation correction and flux accounting across all outlet faces.

Definition at line 862 of file test_boundaries.c.

863{
864 const BCFace faces[] = {
868 };
869 SimCtx *simCtx = NULL;
870 UserCtx *user = NULL;
871 BCContext ctx;
872 Cmpnts ***ucont = NULL;
873 Cmpnts ***ubcs = NULL;
874
875 PetscFunctionBeginUser;
876 PetscCall(PetscMemzero(&ctx, sizeof(ctx)));
877 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 6, 6, 6));
878 PetscCall(PicurvPopulateIdentityMetrics(user));
879 PetscCall(VecSet(user->Ucat, 1.0));
880 PetscCall(DMGlobalToLocalBegin(user->fda, user->Ucat, INSERT_VALUES, user->lUcat));
881 PetscCall(DMGlobalToLocalEnd(user->fda, user->Ucat, INSERT_VALUES, user->lUcat));
882
883 for (size_t n = 0; n < sizeof(faces) / sizeof(faces[0]); ++n) {
884 const BCFace face = faces[n];
885 BoundaryCondition *bc = NULL;
886 PetscReal global_inflow = 0.0;
887 PetscReal global_farfield_in = 0.0;
888 PetscReal global_farfield_out = 0.0;
889 PetscReal global_outflow = 0.0;
890 PetscReal local_inflow = 0.0;
891 PetscReal local_outflow = 0.0;
892 PetscReal area = GetFaceInteriorPointCount(user, face);
893 PetscInt ucont_k, ucont_j, ucont_i, ubcs_k, ubcs_j, ubcs_i;
894
895 PetscCall(VecSet(user->Ucont, 0.0));
896 PetscCall(VecSet(user->Bcs.Ubcs, 0.0));
898 user->boundary_faces[face].face_id = face;
901 user->simCtx->AreaOutSum = area;
902
903 ctx.user = user;
904 ctx.face_id = face;
905 ctx.global_inflow_sum = &global_inflow;
906 ctx.global_farfield_inflow_sum = &global_farfield_in;
907 ctx.global_farfield_outflow_sum = &global_farfield_out;
908 ctx.global_outflow_sum = &global_outflow;
909
911 PetscCall(bc->PreStep(bc, &ctx, &local_inflow, &local_outflow));
912 PetscCall(PicurvAssertRealNear(area, local_outflow, 1.0e-12, "outlet PreStep should measure the uncorrected face flux"));
913
914 global_inflow = area + 0.2 * area;
915 global_outflow = local_outflow;
916 local_outflow = 0.0;
917 PetscCall(bc->Apply(bc, &ctx));
918 PetscCall(bc->PostStep(bc, &ctx, &local_inflow, &local_outflow));
919
920 PetscCall(GetRepresentativeFaceSlots(user, face, &ucont_k, &ucont_j, &ucont_i, &ubcs_k, &ubcs_j, &ubcs_i));
921 PetscCall(DMDAVecGetArrayRead(user->fda, user->Ucont, &ucont));
922 PetscCall(DMDAVecGetArrayRead(user->fda, user->Bcs.Ubcs, &ubcs));
923 PetscCall(PicurvAssertRealNear(1.2, GetFaceNormalComponent(ucont[ucont_k][ucont_j][ucont_i], face), 1.0e-12, "outlet correction should add the configured flux trim"));
924 PetscCall(PicurvAssertRealNear(1.0, ubcs[ubcs_k][ubcs_j][ubcs_i].x, 1.0e-12, "outlet handler should copy Ucat into Ubcs.x"));
925 PetscCall(PicurvAssertRealNear(1.0, ubcs[ubcs_k][ubcs_j][ubcs_i].y, 1.0e-12, "outlet handler should copy Ucat into Ubcs.y"));
926 PetscCall(PicurvAssertRealNear(1.0, ubcs[ubcs_k][ubcs_j][ubcs_i].z, 1.0e-12, "outlet handler should copy Ucat into Ubcs.z"));
927 PetscCall(DMDAVecRestoreArrayRead(user->fda, user->Ucont, &ucont));
928 PetscCall(DMDAVecRestoreArrayRead(user->fda, user->Bcs.Ubcs, &ubcs));
929 PetscCall(PicurvAssertRealNear(0.0, local_inflow, 1.0e-12, "outlet handler should not add to inflow"));
930 PetscCall(PicurvAssertRealNear(1.2 * area, local_outflow, 1.0e-12, "outlet PostStep should report the corrected flux"));
931
932 PetscCall(DestroyBoundaryHandler(&bc));
934 }
935
936 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
937 PetscFunctionReturn(0);
938}
SimCtx * simCtx
Back-pointer to the master simulation context.
Definition variables.h:814
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 unit-boundaries PETSc test binary.

Definition at line 943 of file test_boundaries.c.

944{
945 PetscErrorCode ierr;
946 const PicurvTestCase cases[] = {
947 {"can-rank-service-face-matches-inlet-when-defined", TestCanRankServiceFaceMatchesInletWhenDefined},
948 {"can-rank-service-inlet-face-requires-definition", TestCanRankServiceInletFaceRequiresDefinition},
949 {"boundary-condition-factory-assignments", TestBoundaryConditionFactoryAssignments},
950 {"boundary-condition-factory-implemented-handler-matrix", TestBoundaryConditionFactoryImplementedHandlerMatrix},
951 {"boundary-condition-factory-rejects-unsupported-handler", TestBoundaryConditionFactoryRejectsUnsupportedHandler},
952 {"deterministic-face-grid-location-matrix", TestGetDeterministicFaceGridLocationFaceMatrix},
953 {"random-inlet-face-location-matrix", TestGetRandomCellAndLogicalCoordsOnInletFaceMatrix},
954 {"wall-no-slip-handler-face-matrix", TestWallNoSlipHandlerFaceMatrix},
955 {"inlet-constant-velocity-handler-behavior", TestInletConstantVelocityHandlerBehavior},
956 {"inlet-constant-velocity-handler-face-matrix", TestInletConstantVelocityHandlerFaceMatrix},
957 {"inlet-parabolic-profile-handler-behavior", TestInletParabolicProfileHandlerBehavior},
958 {"inlet-parabolic-profile-handler-face-matrix", TestInletParabolicProfileHandlerFaceMatrix},
959 {"outlet-conservation-handler-behavior", TestOutletConservationHandlerBehavior},
960 {"outlet-conservation-handler-face-matrix", TestOutletConservationHandlerFaceMatrix},
961 };
962
963 ierr = PetscInitialize(&argc, &argv, NULL, "PICurv boundary tests");
964 if (ierr) {
965 return (int)ierr;
966 }
967
968 ierr = PicurvRunTests("unit-boundaries", cases, sizeof(cases) / sizeof(cases[0]));
969 if (ierr) {
970 PetscFinalize();
971 return (int)ierr;
972 }
973
974 ierr = PetscFinalize();
975 return (int)ierr;
976}
static PetscErrorCode TestOutletConservationHandlerFaceMatrix(void)
Tests outlet-conservation correction and flux accounting across all outlet faces.
static PetscErrorCode TestCanRankServiceFaceMatchesInletWhenDefined(void)
Tests that face-service detection matches a defined inlet face.
static PetscErrorCode TestBoundaryConditionFactoryAssignments(void)
Tests the boundary-condition factory assignments for representative handlers.
static PetscErrorCode TestInletParabolicProfileHandlerBehavior(void)
Tests parabolic inlet handler shape on a tiny Z-face.
static PetscErrorCode TestGetDeterministicFaceGridLocationFaceMatrix(void)
Tests deterministic inlet-face grid-location helpers across all faces.
static PetscErrorCode TestInletConstantVelocityHandlerBehavior(void)
Tests constant inlet handler initialization and face writes on a Z inlet.
static PetscErrorCode TestBoundaryConditionFactoryRejectsUnsupportedHandler(void)
Tests that unsupported handler types are rejected by the factory.
static PetscErrorCode TestGetRandomCellAndLogicalCoordsOnInletFaceMatrix(void)
Tests random inlet-face cell selection across all supported faces.
static PetscErrorCode TestBoundaryConditionFactoryImplementedHandlerMatrix(void)
Tests the implemented-handler matrix exposed by the factory.
static PetscErrorCode TestCanRankServiceInletFaceRequiresDefinition(void)
Tests that inlet-face service requires an inlet face to be defined.
static PetscErrorCode TestInletParabolicProfileHandlerFaceMatrix(void)
Tests parabolic-inlet centerline, wall, and flux behavior across all faces.
static PetscErrorCode TestWallNoSlipHandlerFaceMatrix(void)
Tests wall no-slip application across the full face matrix.
static PetscErrorCode TestInletConstantVelocityHandlerFaceMatrix(void)
Tests constant-inlet initialization, flux accounting, and face writes across all faces.
static PetscErrorCode TestOutletConservationHandlerBehavior(void)
Tests outlet conservation handler correction and post-step flux accounting.
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: