PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
Loading...
Searching...
No Matches
ParticleMotion.h
Go to the documentation of this file.
1/**
2 * @file ParticleSwarm.h
3 * @brief Header file for Particle Motion and migration related functions.
4 *
5 * This file contains declarations of functions responsible for moving and migrating particle swarms within a simulation using PETSc's DMSwarm.
6 */
7
8 #ifndef PARTICLE_MOTION_H
9 #define PARTICLE_MOTION_H
10
11// Include necessary headers
12#include <petsc.h> // PETSc library header
13#include <petscdmswarm.h> // PETSc DMSwarm header
14#include <stdbool.h>
15#include <petscsys.h> // For PetscRealloc
16#include <math.h>
17#include "variables.h" // Common type definitions
18#include "logging.h" // Logging macros and definitions
19#include "walkingsearch.h" // Walking search function for particle migration
20/**
21 * @brief Updates a particle's position based on its velocity and the timestep dt (stored in user->dt).
22 *
23 * @param[in] user Pointer to your UserCtx (must contain user->dt).
24 * @param[in,out] position Pointer to the particle's current position (Cmpnts).
25 * @param[in] velocity Pointer to the particle's velocity (Cmpnts).
26 *
27 * @return PetscErrorCode Returns 0 on success, or an error code on failure.
28 */
29 PetscErrorCode UpdateParticlePosition(UserCtx *user, Cmpnts *position, const Cmpnts *velocity);
30
31/**
32 * @brief Loops over all local particles in the DMSwarm, updating their positions
33 * based on velocity and the global timestep user->dt.
34 * @param[in,out] user Pointer to UserCtx (must contain dt).
35 *
36 * @return PetscErrorCode Returns 0 on success, or an error code on failure.
37 */
38 PetscErrorCode UpdateAllParticlePositions(UserCtx *user);
39
40/**
41 * @brief Checks for particles outside the physical domain boundaries and removes them
42 * using DMSwarmRemovePointAtIndex.
43 *
44 * This function iterates through all particles local to the current MPI rank.
45 * It checks if a particle's position (x, y, or z) is outside the specified
46 * physical domain boundaries [xMin, xMax], [yMin, yMax], [zMin, zMax].
47 *
48 * If a particle is found out of bounds, it is removed using DMSwarmRemovePointAtIndex.
49 * NOTE: Removing points changes the indices of subsequent points in the iteration.
50 * Therefore, it's crucial to iterate BACKWARDS or carefully manage indices
51 * after a removal. Iterating backwards is generally safer.
52 *
53 * @param user Pointer to the UserCtx structure.
54 * @param[out] removedCountLocal Pointer to store the number of particles removed *on this rank*.
55 * @param[out] removedCountGlobal Pointer to store the total number of particles removed *across all ranks*.
56 * @param[in] bboxlist An array of BoundingBox structures for ALL MPI ranks, indexed 0 to (size-1).
57 * This array must be up-to-date and available on all ranks.
58 *
59 * @return PetscErrorCode 0 on success, non-zero on failure.
60 */
62 PetscInt *removedCountLocal,
63 PetscInt *removedCountGlobal,
64 const BoundingBox *bboxlist);
65
66/**
67 * @brief Removes particles that have been definitively flagged as LOST by the location algorithm.
68 *
69 * This function is the designated cleanup utility. It should be called after the
70 * `LocateAllParticlesInGrid` orchestrator has run and every particle's status
71 * has been definitively determined.
72 *
73 * It iterates through all locally owned particles and checks their `DMSwarm_location_status`
74 * field. If a particle's status is `LOST`, it is permanently removed from the simulation
75 * using `DMSwarmRemovePointAtIndex`.
76 *
77 * This approach centralizes the removal logic, making the `DMSwarm_location_status`
78 * the single source of truth for a particle's validity, which is more robust than
79 * relying on secondary geometric checks (like bounding boxes).
80 *
81 * @param[in,out] user Pointer to the UserCtx structure containing the swarm.
82 * @param[out] removedCountLocal Pointer to store the number of particles removed on this rank.
83 * @param[out] removedCountGlobal Pointer to store the total number of particles removed across all ranks.
84 *
85 * @return PetscErrorCode 0 on success, or a non-zero PETSc error code on failure.
86 */
87PetscErrorCode CheckAndRemoveLostParticles(UserCtx *user,
88 PetscInt *removedCountLocal,
89 PetscInt *removedCountGlobal);
90
91/**
92 * @brief Defines the basic migration pattern for particles within the swarm.
93 *
94 * This function establishes the migration pattern that dictates how particles
95 * move between different MPI ranks in the simulation. It initializes a migration
96 * list where each particle is assigned a target rank based on predefined conditions.
97 * The migration pattern can be customized to implement various migration behaviors.
98 *
99 * @param[in,out] user Pointer to the UserCtx structure containing simulation context.
100 *
101 * @return PetscErrorCode Returns 0 on success, non-zero on failure.
102 */
104
105/**
106 * @brief Performs the basic migration of particles based on the defined migration pattern.
107 *
108 * This function updates the positions of particles within the swarm by migrating them
109 * to target MPI ranks as specified in the migration list. It handles the migration process
110 * by setting the 'DMSwarm_rank' field for each particle and invokes the DMSwarm migration
111 * mechanism to relocate particles across MPI processes. After migration, it cleans up
112 * allocated resources and ensures synchronization across all MPI ranks.
113 *
114 * @param[in,out] user Pointer to the UserCtx structure containing simulation context.
115 *
116 * @return PetscErrorCode Returns 0 on success, non-zero on failure.
117 */
118PetscErrorCode PerformBasicMigration(UserCtx* user);
119
120/**
121 * @brief Identifies particles leaving the local bounding box and finds their target neighbor rank.
122 *
123 * Iterates local particles, checks against local bounding box. If outside, checks
124 * the pre-computed immediate neighbors (user->neighbors) using the global bboxlist
125 * to see if the particle landed in one of them. Populates the migrationList.
126 * Does NOT handle particles leaving the global domain (assumes CheckAndRemove was called).
127 *
128 * @param user Pointer to the UserCtx (contains local bbox and neighbors).
129 * @param bboxlist Array of BoundingBox structs for all ranks (for checking neighbor boxes).
130 * @param migrationList Pointer to an array of MigrationInfo structs (output, allocated/reallocated by this func).
131 * @param migrationCount Pointer to the number of particles marked for migration (output).
132 * @param listCapacity Pointer to the current allocated capacity of migrationList (in/out).
133 *
134 * @return PetscErrorCode 0 on success, non-zero on failure.
135 */
137 const BoundingBox *bboxlist,
138 MigrationInfo **migrationList,
139 PetscInt *migrationCount,
140 PetscInt *listCapacity);
141
142// --- Helper function to set migration rank field ---
143// This needs to be called AFTER identifying migrating particles
144PetscErrorCode SetMigrationRanks(UserCtx* user, const MigrationInfo *migrationList, PetscInt migrationCount);
145
146/**
147 * @brief Performs particle migration based on the pre-populated DMSwarmPICField_rank field.
148 *
149 * Assumes SetMigrationRanks has already been called to mark particles with their target ranks.
150 * Calls DMSwarmMigrate to execute the communication and removal of un-migrated particles.
151 *
152 * @param user Pointer to the UserCtx structure containing the swarm.
153 *
154 * @return PetscErrorCode 0 on success, non-zero on failure.
155 */
156PetscErrorCode PerformMigration(UserCtx *user);
157
158/**
159 * @brief Counts particles in each cell of the DMDA 'da' and stores the result in user->ParticleCount.
160 *
161 * Assumes user->ParticleCount is a pre-allocated global vector associated with user->da
162 * and initialized to zero before calling this function (though it resets it internally).
163 * Assumes particle 'DMSwarm_CellID' field contains local cell indices.
164 *
165 * @param[in,out] user Pointer to the UserCtx structure containing da, swarm, and ParticleCount.
166 * @return PetscErrorCode Returns 0 on success, non-zero on failure.
167 */
168PetscErrorCode CalculateParticleCountPerCell(UserCtx *user);
169
170// --- Helper function to resize swarm globally (add or remove) ---
171// This assumes removing excess particles means removing the globally last ones.
172PetscErrorCode ResizeSwarmGlobally(DM swarm, PetscInt N_target);
173
174/**
175 * @brief Checks particle count in the reference file and resizes the swarm if needed.
176 *
177 * Reads the specified field file (e.g., position) into a temporary Vec to determine
178 * the number of particles (`N_file`) represented in that file for the given timestep.
179 * Compares `N_file` with the current swarm size (`N_current`). If they differ,
180 * resizes the swarm globally (adds or removes particles) to match `N_file`.
181 * Removal assumes excess particles are the globally last ones.
182 *
183 * @param[in,out] user Pointer to the UserCtx structure containing the DMSwarm.
184 * @param[in] fieldName Name of the reference field (e.g., "position").
185 * @param[in] ti Time index for constructing the file name.
186 * @param[in] ext File extension (e.g., "dat").
187 * @param[out] skipStep Pointer to boolean flag, set to PETSC_TRUE if the step
188 * should be skipped (e.g., file not found), PETSC_FALSE otherwise.
189 *
190 * @return PetscErrorCode 0 on success, non-zero on critical failure.
191 * If the reference file is not found, returns 0 and sets skipStep = PETSC_TRUE.
192 */
193PetscErrorCode PreCheckAndResizeSwarm(UserCtx *user, PetscInt ti, const char *ext);
194
195/**
196 * @brief Performs one full cycle of particle migration: identify, set ranks, and migrate.
197 *
198 * This function encapsulates the three main steps of migrating particles between MPI ranks:
199 * 1. Identify particles on the local rank that need to move based on their current
200 * positions and the domain decomposition (`bboxlist`).
201 * 2. Determine the destination rank for each migrating particle.
202 * 3. Perform the actual migration using PETSc's `DMSwarmMigrate`.
203 * It also calculates and logs the global number of particles migrated.
204 *
205 * @param user Pointer to the UserCtx structure.
206 * @param bboxlist Array of BoundingBox structures defining the spatial domain of each MPI rank.
207 * @param migrationList_p Pointer to a pointer for the MigrationInfo array. This array will be
208 * allocated/reallocated by `IdentifyMigratingParticles` if necessary.
209 * The caller is responsible for freeing this list eventually.
210 * @param migrationCount_p Pointer to store the number of particles identified for migration
211 * on the local rank. This is reset to 0 after migration for the current cycle.
212 * @param migrationListCapacity_p Pointer to store the current capacity of the `migrationList_p` array.
213 * @param currentTime Current simulation time (used for logging).
214 * @param step Current simulation step number (used for logging).
215 * @param migrationCycleName A descriptive name for this migration cycle (e.g., "Preliminary Sort", "Main Loop")
216 * for logging purposes.
217 * @param[out] globalMigrationCount_out Pointer to store the total number of particles migrated
218 * across all MPI ranks during this cycle.
219 * @return PetscErrorCode 0 on success, non-zero on failure.
220 */
221PetscErrorCode PerformSingleParticleMigrationCycle(UserCtx *user, const BoundingBox *bboxlist,
222 MigrationInfo **migrationList_p, PetscInt *migrationCount_p,
223 PetscInt *migrationListCapacity_p,
224 PetscReal currentTime, PetscInt step, const char *migrationCycleName,
225 PetscInt *globalMigrationCount_out);
226
227
228/**
229 * @brief Re-initializes the positions of particles currently on this rank if this rank owns
230 * part of the designated inlet surface.
231 *
232 * This function is intended for `user->ParticleInitialization == 0` (Surface Initialization mode)
233 * and is typically called after an initial migration step (e.g., in `PerformInitialSetup`).
234 * It ensures that all particles that should originate from the inlet surface and are now
235 * on the correct MPI rank are properly distributed across that rank's portion of the inlet.
236 *
237 * @param user Pointer to the UserCtx structure, containing simulation settings and grid information.
238 * @param currentTime Current simulation time (used for logging).
239 * @param step Current simulation step number (used for logging).
240 * @return PetscErrorCode 0 on success, non-zero on failure.
241 */
242PetscErrorCode ReinitializeParticlesOnInletSurface(UserCtx *user, PetscReal currentTime, PetscInt step);
243
244/**
245 * @brief Creates a sorted snapshot of all Particle IDs (PIDs) from a raw data array.
246 * @ingroup ParticleUtils
247 *
248 * This function is a crucial helper for the migration process. It captures the state of
249 * which particles are on the current MPI rank *before* migration occurs by taking a
250 * pointer to the swarm's raw PID data array. The resulting sorted array can then be used
251 * with an efficient binary search to quickly identify newcomer particles after migration.
252 *
253 * This function does NOT call DMSwarmGetField/RestoreField. It is the caller's
254 * responsibility to acquire the `pid_field` pointer before calling and restore it afterward.
255 *
256 * @param[in] pid_field A read-only pointer to the raw array of PIDs for the local swarm.
257 * @param[in] n_local The number of particles currently on the local rank.
258 * @param[out] pids_snapshot_out A pointer to a `PetscInt64*` array. This function will
259 * allocate memory for this array, and the caller is
260 * responsible for freeing it with `PetscFree()` when it
261 * is no longer needed.
262 *
263 * @return PetscErrorCode 0 on success, or a non-zero PETSc error code on failure.
264 */
265PetscErrorCode GetLocalPIDSnapshot(const PetscInt64 pid_field[],
266 PetscInt n_local,
267 PetscInt64 **pids_snapshot_out);
268
269/**
270 * @brief Safely adds a new migration task to a dynamically sized list.
271 *
272 * This utility function manages a dynamic array of MigrationInfo structs. It appends
273 * a new entry to the list and automatically doubles the array's capacity using
274 * `PetscRealloc` if the current capacity is exceeded. This prevents buffer overflows
275 * and avoids the need to know the number of migrating particles in advance.
276 *
277 * @param[in,out] migration_list_p A pointer to the MigrationInfo array pointer. The function
278 * will update this pointer if the array is reallocated.
279 * @param[in,out] capacity_p A pointer to an integer holding the current allocated
280 * capacity of the list (in number of elements). This will be
281 * updated upon reallocation.
282 * @param[in,out] count_p A pointer to an integer holding the current number of
283 * items in the list. This will be incremented by one.
284 * @param[in] particle_local_idx The local index (from 0 to nlocal-1) of the particle
285 * that needs to be migrated.
286 * @param[in] destination_rank The target MPI rank for the particle.
287 *
288 * @return PetscErrorCode 0 on success, or a non-zero PETSc error code on failure (e.g., from memory allocation).
289 */
290PetscErrorCode AddToMigrationList(MigrationInfo **migration_list_p,
291 PetscInt *capacity_p,
292 PetscInt *count_p,
293 PetscInt particle_local_idx,
294 PetscMPIInt destination_rank);
295
296/**
297 * @brief Identifies newly arrived particles after migration and flags them for a location search.
298 * @ingroup ParticleMotion
299 *
300 * This function is a critical component of the iterative migration process managed by
301 * the main particle settlement orchestrator (e.g., `SettleParticles`). After a
302 * `DMSwarmMigrate` call, each rank's local particle list is a new mix of resident
303 * particles and newly received ones. This function's job is to efficiently identify
304 * these "newcomers" and set their `DMSwarm_location_status` field to `NEEDS_LOCATION`.
305 *
306 * This ensures that in the subsequent pass of the migration `do-while` loop, only the
307 * newly arrived particles are processed by the expensive location algorithm, preventing
308 * redundant work on particles that are already settled on the current rank.
309 *
310 * The identification is done by comparing the PIDs of particles currently on the rank
311 * against a "snapshot" of PIDs taken *before* the migration occurred.
312 *
313 * @param[in] swarm The DMSwarm object, which has just completed a migration.
314 * @param[in] n_local_before The number of particles that were on this rank *before* the
315 * migration was performed.
316 * @param[in] pids_before A pre-sorted array of the PIDs that were on this rank before
317 * the migration. This is used for fast lookups.
318 *
319 * @return PetscErrorCode 0 on success, or a non-zero PETSc error code on failure.
320 *
321 * @note This function assumes the `pids_before` array is sorted in ascending order to
322 * enable the use of an efficient binary search.
323 */
324PetscErrorCode FlagNewcomersForLocation(DM swarm,
325 PetscInt n_local_before,
326 const PetscInt64 pids_before[]);
327
328
329/**
330 * @brief Locates all particles within the grid and calculates their interpolation weights.
331 * @ingroup ParticleLocation
332 *
333 * This function iterates through all particles currently local to this MPI rank.
334 * For each particle, it first checks if the particle is within the rank's
335 * pre-calculated bounding box (`user->bbox`). If it is, it calls the
336 * `LocateParticleInGrid` function to perform the walking search.
337 *
338 * `LocateParticleInGrid` is responsible for finding the containing cell `(i,j,k)`
339 * and calculating the corresponding interpolation weights `(w1,w2,w3)`. It updates
340 * the `particle->cell` and `particle->weights` fields directly upon success.
341 * If the search fails (particle not found within MAX_TRAVERSAL, goes out of bounds,
342 * or gets stuck without resolution), `LocateParticleInGrid` sets the particle's
343 * `cell` to `{-1,-1,-1}` and `weights` to `{0.0, 0.0, 0.0}`.
344 *
345 * After attempting location, this function updates the corresponding entries in the
346 * DMSwarm's "DMSwarm_CellID" and "weight" fields using the potentially modified
347 * data from the `particle` struct.
348 *
349 * @param[in] user Pointer to the UserCtx structure containing grid, swarm, and bounding box info.
350 *
351 * @return PetscErrorCode Returns `0` on success, non-zero on failure (e.g., errors accessing DMSwarm fields).
352 *
353 * @note Assumes `user->bbox` is correctly initialized for the local rank.
354 * @note Assumes `InitializeParticle` correctly populates the temporary `particle` struct.
355 * @note Assumes `UpdateSwarmFields` correctly writes data back to the DMSwarm.
356 */
357PetscErrorCode LocateAllParticlesInGrid(UserCtx *user);
358
359/**
360 * @brief Orchestrates the complete particle location and migration process for one timestep.
361 * @ingroup ParticleLocation
362 *
363 * This function is the master orchestrator for ensuring every particle is on its correct
364 * MPI rank and has a valid host cell index. It is designed to be called once per
365 * timestep after particle positions have been updated.
366 *
367 * The function uses a robust, iterative "Guess and Verify" strategy within a
368 * do-while loop to handle complex particle motion across processor boundaries,
369 * especially on curvilinear grids.
370 *
371 * 1. **State Snapshot:** At the start of each pass, it captures a list of all Particle IDs (PIDs)
372 * on the current rank.
373 * 2. **"Guess" (Heuristic):** For particles that are "lost" (no valid host cell),
374 * it first attempts a fast, bounding-box-based guess to find a potential new owner rank.
375 * 3. **"Verify" (Robust Walk):** For all other particles, or if the guess fails,
376 * it uses a robust cell-walking algorithm (`LocateParticleOrFindMigrationTarget`)
377 * that determines the particle's status: located locally, needs migration, or is lost.
378 * 4. **Migration:** After identifying all migrating particles on a pass, it performs the
379 * MPI communication using the `SetMigrationRanks` and `PerformMigration` helpers.
380 * 5. **Newcomer Flagging:** After migration, it uses the PID snapshot from step 1 to
381 * efficiently identify newly arrived particles and flag them for location on the next pass.
382 * 6. **Iteration:** The process repeats in a `do-while` loop until a pass occurs where
383 * no particles migrate, ensuring the entire swarm is in a stable, consistent state.
384 *
385 * @param[in,out] user Pointer to the UserCtx, containing the swarm and all necessary
386 * domain topology information (bboxlist, RankCellInfoMap, etc.).
387 * @param[in] bboxlist An array of BoundingBox structures for ALL MPI ranks, indexed 0 to (size-1).
388 * This array must be up-to-date and available on all ranks.
389 * @return PetscErrorCode 0 on success, or a non-zero PETSc error code on failure.
390 */
391PetscErrorCode LocateAllParticlesInGrid_TEST(UserCtx *user,BoundingBox *bboxlist);
392
393/**
394 * This function is designed to be called at the end of a full timestep, after all
395 * particle-based calculations are complete. It prepares the swarm for the next
396 * timestep by ensuring that after the next position update, every particle will be
397 * re-evaluated by the LocateAllParticlesInGrid orchestrator.
398 *
399 * It iterates through all locally owned particles and sets their
400 * `DMSwarm_location_status` field to `NEEDS_LOCATION`.
401 *
402 * @param[in,out] user Pointer to the UserCtx containing the swarm.
403 * @return PetscErrorCode 0 on success, or a non-zero PETSc error code on failure.
404 */
405PetscErrorCode ResetAllParticleStatuses(UserCtx *user);
406
407 #endif // PARTICLE_MOTION_H
PetscErrorCode UpdateParticlePosition(UserCtx *user, Cmpnts *position, const Cmpnts *velocity)
Updates a particle's position based on its velocity and the timestep dt (stored in user->dt).
PetscErrorCode ResizeSwarmGlobally(DM swarm, PetscInt N_target)
PetscErrorCode AddToMigrationList(MigrationInfo **migration_list_p, PetscInt *capacity_p, PetscInt *count_p, PetscInt particle_local_idx, PetscMPIInt destination_rank)
Safely adds a new migration task to a dynamically sized list.
PetscErrorCode SetMigrationRanks(UserCtx *user, const MigrationInfo *migrationList, PetscInt migrationCount)
Sets the target rank field (DMSwarmPICField_rank) for particles scheduled for migration.
PetscErrorCode CheckAndRemoveOutOfBoundsParticles(UserCtx *user, PetscInt *removedCountLocal, PetscInt *removedCountGlobal, const BoundingBox *bboxlist)
Checks for particles outside the physical domain boundaries and removes them using DMSwarmRemovePoint...
PetscErrorCode GetLocalPIDSnapshot(const PetscInt64 pid_field[], PetscInt n_local, PetscInt64 **pids_snapshot_out)
Creates a sorted snapshot of all Particle IDs (PIDs) from a raw data array.
PetscErrorCode LocateAllParticlesInGrid_TEST(UserCtx *user, BoundingBox *bboxlist)
Orchestrates the complete particle location and migration process for one timestep.
PetscErrorCode UpdateAllParticlePositions(UserCtx *user)
Loops over all local particles in the DMSwarm, updating their positions based on velocity and the glo...
PetscErrorCode ResetAllParticleStatuses(UserCtx *user)
This function is designed to be called at the end of a full timestep, after all particle-based calcul...
PetscErrorCode DefineBasicMigrationPattern(UserCtx *user)
Defines the basic migration pattern for particles within the swarm.
PetscErrorCode LocateAllParticlesInGrid(UserCtx *user)
Locates all particles within the grid and calculates their interpolation weights.
PetscErrorCode PerformSingleParticleMigrationCycle(UserCtx *user, const BoundingBox *bboxlist, MigrationInfo **migrationList_p, PetscInt *migrationCount_p, PetscInt *migrationListCapacity_p, PetscReal currentTime, PetscInt step, const char *migrationCycleName, PetscInt *globalMigrationCount_out)
Performs one full cycle of particle migration: identify, set ranks, and migrate.
PetscErrorCode ReinitializeParticlesOnInletSurface(UserCtx *user, PetscReal currentTime, PetscInt step)
Re-initializes the positions of particles currently on this rank if this rank owns part of the design...
PetscErrorCode CheckAndRemoveLostParticles(UserCtx *user, PetscInt *removedCountLocal, PetscInt *removedCountGlobal)
Removes particles that have been definitively flagged as LOST by the location algorithm.
PetscErrorCode IdentifyMigratingParticles(UserCtx *user, const BoundingBox *bboxlist, MigrationInfo **migrationList, PetscInt *migrationCount, PetscInt *listCapacity)
Identifies particles leaving the local bounding box and finds their target neighbor rank.
PetscErrorCode FlagNewcomersForLocation(DM swarm, PetscInt n_local_before, const PetscInt64 pids_before[])
Identifies newly arrived particles after migration and flags them for a location search.
PetscErrorCode PreCheckAndResizeSwarm(UserCtx *user, PetscInt ti, const char *ext)
Checks particle count in the reference file and resizes the swarm if needed.
PetscErrorCode PerformBasicMigration(UserCtx *user)
Performs the basic migration of particles based on the defined migration pattern.
PetscErrorCode PerformMigration(UserCtx *user)
Performs particle migration based on the pre-populated DMSwarmPICField_rank field.
PetscErrorCode CalculateParticleCountPerCell(UserCtx *user)
Counts particles in each cell of the DMDA 'da' and stores the result in user->ParticleCount.
Logging utilities and macros for PETSc-based applications.
Main header file for a complex fluid dynamics solver.
Defines a 3D axis-aligned bounding box.
Definition variables.h:153
A 3D point or vector with PetscScalar components.
Definition variables.h:99
Information needed to migrate a single particle between MPI ranks.
Definition variables.h:188
User-defined context containing data specific to a single computational grid level.
Definition variables.h:630
Header file for particle location functions using the walking search algorithm.