PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
Loading...
Searching...
No Matches
test_grid.c
Go to the documentation of this file.
1/**
2 * @file test_grid.c
3 * @brief C test module for PICurv.
4 */
5
6#include "test_support.h"
7
8#include "grid.h"
9
10#include <stdlib.h>
11/**
12 * @brief Test-local routine.
13 */
14
15static PetscErrorCode TestComputeLocalBoundingBoxUniformGrid(void)
16{
17 SimCtx *simCtx = NULL;
18 UserCtx *user = NULL;
19 BoundingBox bbox;
20
21 PetscFunctionBeginUser;
22 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
23 PetscCall(ComputeLocalBoundingBox(user, &bbox));
24
25 PetscCall(PicurvAssertRealNear(-1.0e-6, bbox.min_coords.x, 1.0e-10, "bbox min x"));
26 PetscCall(PicurvAssertRealNear(-1.0e-6, bbox.min_coords.y, 1.0e-10, "bbox min y"));
27 PetscCall(PicurvAssertRealNear(-1.0e-6, bbox.min_coords.z, 1.0e-10, "bbox min z"));
28 PetscCall(PicurvAssertRealNear(3.0 + 1.0e-6, bbox.max_coords.x, 1.0e-10, "bbox max x"));
29 PetscCall(PicurvAssertRealNear(3.0 + 1.0e-6, bbox.max_coords.y, 1.0e-10, "bbox max y"));
30 PetscCall(PicurvAssertRealNear(3.0 + 1.0e-6, bbox.max_coords.z, 1.0e-10, "bbox max z"));
31
32 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
33 PetscFunctionReturn(0);
34}
35/**
36 * @brief Test-local routine.
37 */
38
39static PetscErrorCode TestGatherAndBroadcastBoundingBoxes(void)
40{
41 SimCtx *simCtx = NULL;
42 UserCtx *user = NULL;
43 BoundingBox local_bbox;
44 BoundingBox *boxes = NULL;
45 PetscMPIInt rank = 0, size = 1;
46
47 PetscFunctionBeginUser;
48 PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
49 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
50
51 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
52 PetscCall(ComputeLocalBoundingBox(user, &local_bbox));
53 PetscCall(GatherAllBoundingBoxes(user, &boxes));
54
55 if (rank == 0) {
56 PetscCall(PicurvAssertBool((PetscBool)(boxes != NULL), "root must receive gathered bounding boxes"));
57 } else {
58 PetscCall(PicurvAssertBool((PetscBool)(boxes == NULL), "non-root should have NULL pre-broadcast bbox pointer"));
59 }
60
61 PetscCall(BroadcastAllBoundingBoxes(user, &boxes));
62 PetscCall(PicurvAssertBool((PetscBool)(boxes != NULL), "all ranks must have bbox array after broadcast"));
63
64 for (PetscMPIInt r = 0; r < size; ++r) {
65 PetscCall(PicurvAssertBool((PetscBool)(boxes[r].min_coords.x <= boxes[r].max_coords.x), "bbox x range is valid"));
66 PetscCall(PicurvAssertBool((PetscBool)(boxes[r].min_coords.y <= boxes[r].max_coords.y), "bbox y range is valid"));
67 PetscCall(PicurvAssertBool((PetscBool)(boxes[r].min_coords.z <= boxes[r].max_coords.z), "bbox z range is valid"));
68 }
69
70 PetscCall(PicurvAssertRealNear(local_bbox.min_coords.x, boxes[rank].min_coords.x, 1.0e-10, "rank-local bbox min x"));
71 PetscCall(PicurvAssertRealNear(local_bbox.max_coords.x, boxes[rank].max_coords.x, 1.0e-10, "rank-local bbox max x"));
72
73 free(boxes);
74 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
75 PetscFunctionReturn(0);
76}
77/**
78 * @brief Entry point for this unit-test binary.
79 */
80
81int main(int argc, char **argv)
82{
83 PetscErrorCode ierr;
84 const PicurvTestCase cases[] = {
85 {"compute-local-bbox-uniform-grid", TestComputeLocalBoundingBoxUniformGrid},
86 {"gather-and-broadcast-bboxes", TestGatherAndBroadcastBoundingBoxes},
87 };
88
89 ierr = PetscInitialize(&argc, &argv, NULL, "PICurv grid tests");
90 if (ierr) {
91 return (int)ierr;
92 }
93
94 ierr = PicurvRunTests("unit-grid", cases, sizeof(cases) / sizeof(cases[0]));
95 if (ierr) {
96 PetscFinalize();
97 return (int)ierr;
98 }
99
100 ierr = PetscFinalize();
101 return (int)ierr;
102}
Public interface for grid, solver, and metric setup routines.
PetscErrorCode BroadcastAllBoundingBoxes(UserCtx *user, BoundingBox **bboxlist)
Broadcasts the bounding box information collected on rank 0 to all other ranks.
Definition grid.c:879
PetscErrorCode ComputeLocalBoundingBox(UserCtx *user, BoundingBox *localBBox)
Computes the local bounding box of the grid on the current process.
Definition grid.c:675
PetscErrorCode GatherAllBoundingBoxes(UserCtx *user, BoundingBox **allBBoxes)
Gathers local bounding boxes from all MPI processes to rank 0.
Definition grid.c:817
static PetscErrorCode TestComputeLocalBoundingBoxUniformGrid(void)
Test-local routine.
Definition test_grid.c:15
static PetscErrorCode TestGatherAndBroadcastBoundingBoxes(void)
Test-local routine.
Definition test_grid.c:39
int main(int argc, char **argv)
Entry point for this unit-test binary.
Definition test_grid.c:81
PetscErrorCode PicurvCreateMinimalContexts(SimCtx **simCtx_out, UserCtx **user_out, PetscInt mx, PetscInt my, PetscInt mz)
Shared test-support routine.
PetscErrorCode PicurvAssertRealNear(PetscReal expected, PetscReal actual, PetscReal tol, const char *context)
Shared test-support routine.
PetscErrorCode PicurvDestroyMinimalContexts(SimCtx **simCtx_ptr, UserCtx **user_ptr)
Shared test-support routine.
PetscErrorCode PicurvRunTests(const char *suite_name, const PicurvTestCase *cases, size_t case_count)
Shared test-support routine.
PetscErrorCode PicurvAssertBool(PetscBool value, const char *context)
Shared test-support routine.
C test module for PICurv.
Named test case descriptor consumed by PicurvRunTests.
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
Defines a 3D axis-aligned bounding box.
Definition variables.h:154
The master context for the entire simulation.
Definition variables.h:591
User-defined context containing data specific to a single computational grid level.
Definition variables.h:738