PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
Loading...
Searching...
No Matches
test_particle_kernels.c
Go to the documentation of this file.
1/**
2 * @file test_particle_kernels.c
3 * @brief C test module for PICurv.
4 */
5
6#include "test_support.h"
7
8#include "walkingsearch.h"
9/**
10 * @brief Test-local routine.
11 */
12
13static PetscErrorCode TestCheckCellWithinLocalGrid(void)
14{
15 SimCtx *simCtx = NULL;
16 UserCtx *user = NULL;
17 PetscBool within = PETSC_FALSE;
18
19 PetscFunctionBeginUser;
20 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
21
22 PetscCall(CheckCellWithinLocalGrid(user, 1, 1, 1, &within));
23 PetscCall(PicurvAssertBool(within, "cell (1,1,1) should be within the serial local grid"));
24
25 PetscCall(CheckCellWithinLocalGrid(user, 3, 1, 1, &within));
26 PetscCall(PicurvAssertBool((PetscBool)!within, "cell (3,1,1) should fall outside the valid local cell range"));
27
28 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
29 PetscFunctionReturn(0);
30}
31/**
32 * @brief Test-local routine.
33 */
34
35static PetscErrorCode TestInitializeTraversalParameters(void)
36{
37 SimCtx *simCtx = NULL;
38 UserCtx *user = NULL;
39 Particle particle;
40 PetscInt idx = -1, idy = -1, idz = -1, steps = -1;
41
42 PetscFunctionBeginUser;
43 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
44 PetscCall(PetscMemzero(&particle, sizeof(particle)));
45
46 particle.PID = 42;
47 particle.cell[0] = 1;
48 particle.cell[1] = 1;
49 particle.cell[2] = 1;
50
51 PetscCall(InitializeTraversalParameters(user, &particle, &idx, &idy, &idz, &steps));
52 PetscCall(PicurvAssertIntEqual(1, idx, "InitializeTraversalParameters should preserve prior i"));
53 PetscCall(PicurvAssertIntEqual(1, idy, "InitializeTraversalParameters should preserve prior j"));
54 PetscCall(PicurvAssertIntEqual(1, idz, "InitializeTraversalParameters should preserve prior k"));
55 PetscCall(PicurvAssertIntEqual(0, steps, "InitializeTraversalParameters should reset traversal steps"));
56
57 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
58 PetscFunctionReturn(0);
59}
60/**
61 * @brief Test-local routine.
62 */
63
64static PetscErrorCode TestRetrieveCurrentCell(void)
65{
66 SimCtx *simCtx = NULL;
67 UserCtx *user = NULL;
68 Cell cell;
69
70 PetscFunctionBeginUser;
71 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 4, 4, 4));
72 PetscCall(PetscMemzero(&cell, sizeof(cell)));
73
74 PetscCall(RetrieveCurrentCell(user, 1, 1, 1, &cell));
75 PetscCall(PicurvAssertRealNear(1.0, cell.vertices[0].x, 1.0e-12, "vertex 0 x coordinate"));
76 PetscCall(PicurvAssertRealNear(1.0, cell.vertices[0].y, 1.0e-12, "vertex 0 y coordinate"));
77 PetscCall(PicurvAssertRealNear(1.0, cell.vertices[0].z, 1.0e-12, "vertex 0 z coordinate"));
78 PetscCall(PicurvAssertRealNear(2.0, cell.vertices[5].x, 1.0e-12, "vertex 5 x coordinate"));
79 PetscCall(PicurvAssertRealNear(2.0, cell.vertices[5].y, 1.0e-12, "vertex 5 y coordinate"));
80 PetscCall(PicurvAssertRealNear(2.0, cell.vertices[5].z, 1.0e-12, "vertex 5 z coordinate"));
81
82 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
83 PetscFunctionReturn(0);
84}
85/**
86 * @brief Entry point for this unit-test binary.
87 */
88
89int main(int argc, char **argv)
90{
91 PetscErrorCode ierr;
92 const PicurvTestCase cases[] = {
93 {"check-cell-within-local-grid", TestCheckCellWithinLocalGrid},
94 {"initialize-traversal-parameters", TestInitializeTraversalParameters},
95 {"retrieve-current-cell", TestRetrieveCurrentCell},
96 };
97
98 ierr = PetscInitialize(&argc, &argv, NULL, "PICurv particle kernel tests");
99 if (ierr) {
100 return (int)ierr;
101 }
102
103 ierr = PicurvRunTests("unit-particles", cases, sizeof(cases) / sizeof(cases[0]));
104 if (ierr) {
105 PetscFinalize();
106 return (int)ierr;
107 }
108
109 ierr = PetscFinalize();
110 return (int)ierr;
111}
static PetscErrorCode TestCheckCellWithinLocalGrid(void)
Test-local routine.
int main(int argc, char **argv)
Entry point for this unit-test binary.
static PetscErrorCode TestInitializeTraversalParameters(void)
Test-local routine.
static PetscErrorCode TestRetrieveCurrentCell(void)
Test-local routine.
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 PicurvAssertIntEqual(PetscInt expected, PetscInt actual, const char *context)
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.
PetscInt cell[3]
Definition variables.h:167
PetscScalar x
Definition variables.h:101
PetscScalar z
Definition variables.h:101
PetscScalar y
Definition variables.h:101
PetscInt64 PID
Definition variables.h:166
Cmpnts vertices[8]
Coordinates of the eight vertices of the cell.
Definition variables.h:161
Defines the vertices of a single hexahedral grid cell.
Definition variables.h:160
Defines a particle's core properties for Lagrangian tracking.
Definition variables.h:165
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
Header file for particle location functions using the walking search algorithm.
PetscErrorCode RetrieveCurrentCell(UserCtx *user, PetscInt idx, PetscInt idy, PetscInt idz, Cell *cell)
Retrieves the coordinates of the eight vertices of the current cell.
PetscErrorCode InitializeTraversalParameters(UserCtx *user, Particle *particle, PetscInt *idx, PetscInt *idy, PetscInt *idz, PetscInt *traversal_steps)
Initializes traversal parameters for locating a particle.
PetscErrorCode CheckCellWithinLocalGrid(UserCtx *user, PetscInt idx, PetscInt idy, PetscInt idz, PetscBool *is_within)
Checks if the current cell indices are within the local grid boundaries.