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

C test module for PICurv. More...

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

Go to the source code of this file.

Functions

static PetscErrorCode TestDistributeParticlesCollectiveConsistency (void)
 Test-local routine.
 
static PetscErrorCode TestBoundingBoxCollectivesMultiRank (void)
 Test-local routine.
 
int main (int argc, char **argv)
 Entry point for this unit-test binary.
 

Detailed Description

C test module for PICurv.

Definition in file test_mpi_kernels.c.

Function Documentation

◆ TestDistributeParticlesCollectiveConsistency()

static PetscErrorCode TestDistributeParticlesCollectiveConsistency ( void  )
static

Test-local routine.

Definition at line 16 of file test_mpi_kernels.c.

17{
18 PetscMPIInt rank = 0, size = 1;
19 PetscInt local_particles = 0;
20 PetscInt remainder = 0;
21 PetscInt global_particles = 0;
22 PetscInt remainder_min = 0, remainder_max = 0;
23 const PetscInt total_particles = 137;
24
25 PetscFunctionBeginUser;
26 PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
27 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
28 PetscCall(PicurvAssertBool((PetscBool)(size >= 2), "unit-mpi requires at least two MPI ranks"));
29
30 PetscCall(DistributeParticles(total_particles, rank, size, &local_particles, &remainder));
31
32 PetscCall(PicurvAssertIntEqual(
33 total_particles / size + (((PetscInt)rank < remainder) ? 1 : 0),
34 local_particles,
35 "local particle share should match quotient+remainder policy"));
36
37 PetscCallMPI(MPI_Allreduce(&local_particles, &global_particles, 1, MPIU_INT, MPI_SUM, PETSC_COMM_WORLD));
38 PetscCall(PicurvAssertIntEqual(total_particles, global_particles, "distributed particle count must conserve total particles"));
39
40 PetscCallMPI(MPI_Allreduce(&remainder, &remainder_min, 1, MPIU_INT, MPI_MIN, PETSC_COMM_WORLD));
41 PetscCallMPI(MPI_Allreduce(&remainder, &remainder_max, 1, MPIU_INT, MPI_MAX, PETSC_COMM_WORLD));
42 PetscCall(PicurvAssertIntEqual(remainder_min, remainder_max, "all ranks should report the same remainder"));
43 PetscFunctionReturn(0);
44}
PetscErrorCode DistributeParticles(PetscInt numParticles, PetscMPIInt rank, PetscMPIInt size, PetscInt *particlesPerProcess, PetscInt *remainder)
Distributes particles evenly across MPI processes, handling any remainders.
PetscErrorCode PicurvAssertIntEqual(PetscInt expected, PetscInt actual, const char *context)
Shared test-support routine.
PetscErrorCode PicurvAssertBool(PetscBool value, const char *context)
Shared test-support routine.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TestBoundingBoxCollectivesMultiRank()

static PetscErrorCode TestBoundingBoxCollectivesMultiRank ( void  )
static

Test-local routine.

Definition at line 49 of file test_mpi_kernels.c.

50{
51 SimCtx *simCtx = NULL;
52 UserCtx *user = NULL;
53 BoundingBox local_bbox;
54 BoundingBox *boxes = NULL;
55 PetscMPIInt rank = 0, size = 1;
56 PetscReal global_min_x = 0.0;
57 PetscReal global_max_x = 0.0;
58
59 PetscFunctionBeginUser;
60 PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
61 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
62 PetscCall(PicurvAssertBool((PetscBool)(size >= 2), "unit-mpi requires at least two MPI ranks"));
63
64 PetscCall(PicurvCreateMinimalContexts(&simCtx, &user, 8, 6, 4));
65 PetscCall(ComputeLocalBoundingBox(user, &local_bbox));
66 PetscCall(GatherAllBoundingBoxes(user, &boxes));
67 PetscCall(BroadcastAllBoundingBoxes(user, &boxes));
68 PetscCall(PicurvAssertBool((PetscBool)(boxes != NULL), "all ranks should hold the gathered bounding-box table"));
69
70 for (PetscMPIInt r = 0; r < size; ++r) {
71 PetscCall(PicurvAssertBool((PetscBool)(boxes[r].min_coords.x <= boxes[r].max_coords.x), "bbox x-range should be valid"));
72 PetscCall(PicurvAssertBool((PetscBool)(boxes[r].min_coords.y <= boxes[r].max_coords.y), "bbox y-range should be valid"));
73 PetscCall(PicurvAssertBool((PetscBool)(boxes[r].min_coords.z <= boxes[r].max_coords.z), "bbox z-range should be valid"));
74 }
75
76 PetscCall(PicurvAssertRealNear(local_bbox.min_coords.x, boxes[rank].min_coords.x, 1.0e-10, "local bbox min x preserved"));
77 PetscCall(PicurvAssertRealNear(local_bbox.max_coords.x, boxes[rank].max_coords.x, 1.0e-10, "local bbox max x preserved"));
78
79 PetscCallMPI(MPI_Allreduce(&local_bbox.min_coords.x, &global_min_x, 1, MPIU_REAL, MPI_MIN, PETSC_COMM_WORLD));
80 PetscCallMPI(MPI_Allreduce(&local_bbox.max_coords.x, &global_max_x, 1, MPIU_REAL, MPI_MAX, PETSC_COMM_WORLD));
81 PetscCall(PicurvAssertBool((PetscBool)(global_min_x <= 0.0), "global min x should include domain start"));
82 PetscCall(PicurvAssertBool((PetscBool)(global_max_x >= 7.0), "global max x should include domain end"));
83
84 free(boxes);
85 PetscCall(PicurvDestroyMinimalContexts(&simCtx, &user));
86 PetscFunctionReturn(0);
87}
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
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.
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
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
Here is the call graph for this function:
Here is the caller graph for this function:

◆ main()

int main ( int  argc,
char **  argv 
)

Entry point for this unit-test binary.

Definition at line 92 of file test_mpi_kernels.c.

93{
94 PetscErrorCode ierr;
95 const PicurvTestCase cases[] = {
96 {"distribute-particles-collective-consistency", TestDistributeParticlesCollectiveConsistency},
97 {"bounding-box-collectives-multi-rank", TestBoundingBoxCollectivesMultiRank},
98 };
99
100 ierr = PetscInitialize(&argc, &argv, NULL, "PICurv MPI-focused runtime tests");
101 if (ierr) {
102 return (int)ierr;
103 }
104
105 ierr = PicurvRunTests("unit-mpi", cases, sizeof(cases) / sizeof(cases[0]));
106 if (ierr) {
107 PetscFinalize();
108 return (int)ierr;
109 }
110
111 ierr = PetscFinalize();
112 return (int)ierr;
113}
static PetscErrorCode TestDistributeParticlesCollectiveConsistency(void)
Test-local routine.
static PetscErrorCode TestBoundingBoxCollectivesMultiRank(void)
Test-local routine.
PetscErrorCode PicurvRunTests(const char *suite_name, const PicurvTestCase *cases, size_t case_count)
Shared test-support routine.
Named test case descriptor consumed by PicurvRunTests.
Here is the call graph for this function: