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

C unit tests for particle walking-search helper routines. More...

#include "test_support.h"
#include "walkingsearch.h"
Include dependency graph for test_particle_kernels.c:

Go to the source code of this file.

Functions

static PetscErrorCode TestCheckCellWithinLocalGrid (void)
 Tests whether a cell index lies within the local grid ownership range.
 
static PetscErrorCode TestInitializeTraversalParameters (void)
 Tests initialization of traversal parameters for particle search.
 
static PetscErrorCode TestRetrieveCurrentCell (void)
 Tests retrieval of the current cell from logical particle coordinates.
 
int main (int argc, char **argv)
 Runs the unit-particles PETSc test binary.
 

Detailed Description

C unit tests for particle walking-search helper routines.

Definition in file test_particle_kernels.c.

Function Documentation

◆ TestCheckCellWithinLocalGrid()

static PetscErrorCode TestCheckCellWithinLocalGrid ( void  )
static

Tests whether a cell index lies within the local grid ownership range.

Definition at line 13 of file test_particle_kernels.c.

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, 4, 1, 1, &within));
26 PetscCall(PicurvAssertBool((PetscBool)!within, "cell (4,1,1) should fall outside the valid local cell range"));
27
28 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
29 PetscFunctionReturn(0);
30}
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.
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
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.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestInitializeTraversalParameters()

static PetscErrorCode TestInitializeTraversalParameters ( void  )
static

Tests initialization of traversal parameters for particle search.

Definition at line 35 of file test_particle_kernels.c.

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}
PetscErrorCode PicurvAssertIntEqual(PetscInt expected, PetscInt actual, const char *context)
Asserts that two integer values are equal.
PetscInt cell[3]
Definition variables.h:167
PetscInt64 PID
Definition variables.h:166
Defines a particle's core properties for Lagrangian tracking.
Definition variables.h:165
PetscErrorCode InitializeTraversalParameters(UserCtx *user, Particle *particle, PetscInt *idx, PetscInt *idy, PetscInt *idz, PetscInt *traversal_steps)
Initializes traversal parameters for locating a particle.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestRetrieveCurrentCell()

static PetscErrorCode TestRetrieveCurrentCell ( void  )
static

Tests retrieval of the current cell from logical particle coordinates.

Definition at line 64 of file test_particle_kernels.c.

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(0.25, cell.vertices[0].x, 1.0e-12, "vertex 0 x coordinate"));
76 PetscCall(PicurvAssertRealNear(0.25, cell.vertices[0].y, 1.0e-12, "vertex 0 y coordinate"));
77 PetscCall(PicurvAssertRealNear(0.25, cell.vertices[0].z, 1.0e-12, "vertex 0 z coordinate"));
78 PetscCall(PicurvAssertRealNear(0.50, cell.vertices[5].x, 1.0e-12, "vertex 5 x coordinate"));
79 PetscCall(PicurvAssertRealNear(0.50, cell.vertices[5].y, 1.0e-12, "vertex 5 y coordinate"));
80 PetscCall(PicurvAssertRealNear(0.50, cell.vertices[5].z, 1.0e-12, "vertex 5 z coordinate"));
81
82 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
83 PetscFunctionReturn(0);
84}
PetscErrorCode PicurvAssertRealNear(PetscReal expected, PetscReal actual, PetscReal tol, const char *context)
Asserts that two real values agree within tolerance.
PetscScalar x
Definition variables.h:101
PetscScalar z
Definition variables.h:101
PetscScalar y
Definition variables.h:101
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
PetscErrorCode RetrieveCurrentCell(UserCtx *user, PetscInt idx, PetscInt idy, PetscInt idz, Cell *cell)
Retrieves the coordinates of the eight vertices of the current cell.
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-particles PETSc test binary.

Definition at line 89 of file test_particle_kernels.c.

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)
Tests whether a cell index lies within the local grid ownership range.
static PetscErrorCode TestInitializeTraversalParameters(void)
Tests initialization of traversal parameters for particle search.
static PetscErrorCode TestRetrieveCurrentCell(void)
Tests retrieval of the current cell from logical particle coordinates.
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: