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

C unit tests for grid bounding-box helper routines. More...

#include "test_support.h"
#include "grid.h"
#include <stdlib.h>
Include dependency graph for test_grid.c:

Go to the source code of this file.

Functions

static PetscErrorCode TestComputeLocalBoundingBoxUniformGrid (void)
 Tests local bounding-box construction on a uniform grid.
 
static PetscErrorCode TestGatherAndBroadcastBoundingBoxes (void)
 Tests collective gather and broadcast of rank-local bounding boxes.
 
int main (int argc, char **argv)
 Runs the unit-grid PETSc test binary.
 

Detailed Description

C unit tests for grid bounding-box helper routines.

Definition in file test_grid.c.

Function Documentation

◆ TestComputeLocalBoundingBoxUniformGrid()

static PetscErrorCode TestComputeLocalBoundingBoxUniformGrid ( void  )
static

Tests local bounding-box construction on a uniform grid.

Definition at line 15 of file test_grid.c.

16{
17 SimCtx *simCtx = NULL;
18 UserCtx *user = NULL;
19 BoundingBox bbox;
20 PetscReal expected_max = 0.0;
21
22 PetscFunctionBeginUser;
23 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
24 PetscCall(ComputeLocalBoundingBox(user, &bbox));
25 expected_max = ((PetscReal)(user->IM - 1) / (PetscReal)user->IM) + 1.0e-6;
26
27 PetscCall(PicurvAssertRealNear(-1.0e-6, bbox.min_coords.x, 1.0e-10, "bbox min x"));
28 PetscCall(PicurvAssertRealNear(-1.0e-6, bbox.min_coords.y, 1.0e-10, "bbox min y"));
29 PetscCall(PicurvAssertRealNear(-1.0e-6, bbox.min_coords.z, 1.0e-10, "bbox min z"));
30 PetscCall(PicurvAssertRealNear(expected_max, bbox.max_coords.x, 1.0e-10, "bbox max x"));
31 PetscCall(PicurvAssertRealNear(expected_max, bbox.max_coords.y, 1.0e-10, "bbox max y"));
32 PetscCall(PicurvAssertRealNear(expected_max, bbox.max_coords.z, 1.0e-10, "bbox max z"));
33
34 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
35 PetscFunctionReturn(0);
36}
PetscErrorCode ComputeLocalBoundingBox(UserCtx *user, BoundingBox *localBBox)
Computes the local bounding box of the grid on the current process.
Definition grid.c:679
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.
Cmpnts max_coords
Maximum x, y, z coordinates of the bounding box.
Definition variables.h:156
Cmpnts min_coords
Minimum x, y, z coordinates of the bounding box.
Definition variables.h:155
PetscScalar x
Definition variables.h:101
PetscScalar z
Definition variables.h:101
PetscScalar y
Definition variables.h:101
PetscInt IM
Definition variables.h:820
Defines a 3D axis-aligned bounding box.
Definition variables.h:154
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:

◆ TestGatherAndBroadcastBoundingBoxes()

static PetscErrorCode TestGatherAndBroadcastBoundingBoxes ( void  )
static

Tests collective gather and broadcast of rank-local bounding boxes.

Definition at line 41 of file test_grid.c.

42{
43 SimCtx *simCtx = NULL;
44 UserCtx *user = NULL;
45 BoundingBox local_bbox;
46 BoundingBox *boxes = NULL;
47 PetscMPIInt rank = 0, size = 1;
48
49 PetscFunctionBeginUser;
50 PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
51 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
52
53 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
54 PetscCall(ComputeLocalBoundingBox(user, &local_bbox));
55 PetscCall(GatherAllBoundingBoxes(user, &boxes));
56
57 if (rank == 0) {
58 PetscCall(PicurvAssertBool((PetscBool)(boxes != NULL), "root must receive gathered bounding boxes"));
59 } else {
60 PetscCall(PicurvAssertBool((PetscBool)(boxes == NULL), "non-root should have NULL pre-broadcast bbox pointer"));
61 }
62
63 PetscCall(BroadcastAllBoundingBoxes(user, &boxes));
64 PetscCall(PicurvAssertBool((PetscBool)(boxes != NULL), "all ranks must have bbox array after broadcast"));
65
66 for (PetscMPIInt r = 0; r < size; ++r) {
67 PetscCall(PicurvAssertBool((PetscBool)(boxes[r].min_coords.x <= boxes[r].max_coords.x), "bbox x range is valid"));
68 PetscCall(PicurvAssertBool((PetscBool)(boxes[r].min_coords.y <= boxes[r].max_coords.y), "bbox y range is valid"));
69 PetscCall(PicurvAssertBool((PetscBool)(boxes[r].min_coords.z <= boxes[r].max_coords.z), "bbox z range is valid"));
70 }
71
72 PetscCall(PicurvAssertRealNear(local_bbox.min_coords.x, boxes[rank].min_coords.x, 1.0e-10, "rank-local bbox min x"));
73 PetscCall(PicurvAssertRealNear(local_bbox.max_coords.x, boxes[rank].max_coords.x, 1.0e-10, "rank-local bbox max x"));
74
75 free(boxes);
76 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
77 PetscFunctionReturn(0);
78}
PetscErrorCode BroadcastAllBoundingBoxes(UserCtx *user, BoundingBox **bboxlist)
Broadcasts the bounding box information collected on rank 0 to all other ranks.
Definition grid.c:883
PetscErrorCode GatherAllBoundingBoxes(UserCtx *user, BoundingBox **allBBoxes)
Gathers local bounding boxes from all MPI processes to rank 0.
Definition grid.c:821
PetscErrorCode PicurvAssertBool(PetscBool value, const char *context)
Asserts that one boolean condition is true.
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-grid PETSc test binary.

Definition at line 83 of file test_grid.c.

84{
85 PetscErrorCode ierr;
86 const PicurvTestCase cases[] = {
87 {"compute-local-bbox-uniform-grid", TestComputeLocalBoundingBoxUniformGrid},
88 {"gather-and-broadcast-bboxes", TestGatherAndBroadcastBoundingBoxes},
89 };
90
91 ierr = PetscInitialize(&argc, &argv, NULL, "PICurv grid tests");
92 if (ierr) {
93 return (int)ierr;
94 }
95
96 ierr = PicurvRunTests("unit-grid", cases, sizeof(cases) / sizeof(cases[0]));
97 if (ierr) {
98 PetscFinalize();
99 return (int)ierr;
100 }
101
102 ierr = PetscFinalize();
103 return (int)ierr;
104}
static PetscErrorCode TestComputeLocalBoundingBoxUniformGrid(void)
Tests local bounding-box construction on a uniform grid.
Definition test_grid.c:15
static PetscErrorCode TestGatherAndBroadcastBoundingBoxes(void)
Tests collective gather and broadcast of rank-local bounding boxes.
Definition test_grid.c:41
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: