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 ParticleMotion.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/**
22 * @brief Generates two independent standard normal random variables N(0,1)
23 * using the Box-Muller transform.
24 *
25 * @param[in] rnd The PETSc Random context (Uniform [0,1)).
26 * @param[out] n1 First Gaussian number.
27 * @param[out] n2 Second Gaussian number.
28 *
29 * @return PetscErrorCode
30 */
31PetscErrorCode GenerateGaussianNoise(PetscRandom rnd, PetscReal *n1, PetscReal *n2);
32
33/**
34 * @brief Calculates the stochastic displacement vector (Brownian motion) for a single particle.
35 * Equation: dX_stoch = sqrt(2 * Gamma_eff * dt) * N(0,1)
36 *
37 * @param[in] user Pointer to UserCtx (access to dt and BrownianMotionRNG).
38 * @param[in] diff_eff The effective diffusivity (Gamma + Gamma_t) at the particle's location.
39 * @param[out] displacement Pointer to a Cmpnts struct to store the resulting (dx, dy, dz).
40 *
41 * @return PetscErrorCode
42 */
43PetscErrorCode CalculateBrownianDisplacement(UserCtx *user, PetscReal diff_eff, Cmpnts *displacement);
44
45/**
46 * @brief Updates a particle's position based on its velocity and the timestep dt (stored in user->dt).
47 *
48 * @param[in] user Pointer to your UserCtx (must contain user->dt).
49 * @param[in,out] particle Pointer to the particle struct (contains, pos,vel,diffusivity etc).
50 *
51 * @return PetscErrorCode Returns 0 on success, or an error code on failure.
52 */
53PetscErrorCode UpdateParticlePosition(UserCtx *user, Particle *particle);
54
55/**
56 * @brief Loops over all local particles in the DMSwarm, updating their positions
57 * based on velocity and the global timestep user->dt.
58 * @param[in,out] user Pointer to UserCtx (must contain dt).
59 *
60 * @return PetscErrorCode Returns 0 on success, or an error code on failure.
61 */
62 PetscErrorCode UpdateAllParticlePositions(UserCtx *user);
63
64/**
65 * @brief Checks for particles outside the physical domain boundaries and removes them
66 * using DMSwarmRemovePointAtIndex.
67 *
68 * This function iterates through all particles local to the current MPI rank.
69 * It checks if a particle's position (x, y, or z) is outside the specified
70 * physical domain boundaries [xMin, xMax], [yMin, yMax], [zMin, zMax].
71 *
72 * If a particle is found out of bounds, it is removed using DMSwarmRemovePointAtIndex.
73 * NOTE: Removing points changes the indices of subsequent points in the iteration.
74 * Therefore, it's crucial to iterate BACKWARDS or carefully manage indices
75 * after a removal. Iterating backwards is generally safer.
76 *
77 * @param user Pointer to the UserCtx structure.
78 * @param[out] removedCountLocal Pointer to store the number of particles removed *on this rank*.
79 * @param[out] removedCountGlobal Pointer to store the total number of particles removed *across all ranks*.
80 * @param[in] bboxlist An array of BoundingBox structures for ALL MPI ranks, indexed 0 to (size-1).
81 * This array must be up-to-date and available on all ranks.
82 *
83 * @return PetscErrorCode 0 on success, non-zero on failure.
84 */
86 PetscInt *removedCountLocal,
87 PetscInt *removedCountGlobal,
88 const BoundingBox *bboxlist);
89
90/**
91 * @brief Removes particles that have been definitively flagged as LOST by the location algorithm.
92 *
93 * This function is the designated cleanup utility. It should be called after the
94 * `LocateAllParticlesInGrid` orchestrator has run and every particle's status
95 * has been definitively determined.
96 *
97 * It iterates through all locally owned particles and checks their `DMSwarm_location_status`
98 * field. If a particle's status is `LOST`, it is permanently removed from the simulation
99 * using `DMSwarmRemovePointAtIndex`.
100 *
101 * This approach centralizes the removal logic, making the `DMSwarm_location_status`
102 * the single source of truth for a particle's validity, which is more robust than
103 * relying on secondary geometric checks (like bounding boxes).
104 *
105 * @param[in,out] user Pointer to the UserCtx structure containing the swarm.
106 * @param[out] removedCountLocal Pointer to store the number of particles removed on this rank.
107 * @param[out] removedCountGlobal Pointer to store the total number of particles removed across all ranks.
108 *
109 * @return PetscErrorCode 0 on success, or a non-zero PETSc error code on failure.
110 */
111PetscErrorCode CheckAndRemoveLostParticles(UserCtx *user,
112 PetscInt *removedCountLocal,
113 PetscInt *removedCountGlobal);
114
115/**
116 * @brief Defines the basic migration pattern for particles within the swarm.
117 *
118 * This function establishes the migration pattern that dictates how particles
119 * move between different MPI ranks in the simulation. It initializes a migration
120 * list where each particle is assigned a target rank based on predefined conditions.
121 * The migration pattern can be customized to implement various migration behaviors.
122 *
123 * @param[in,out] user Pointer to the UserCtx structure containing simulation context.
124 *
125 * @return PetscErrorCode Returns 0 on success, non-zero on failure.
126 */
128
129/**
130 * @brief Performs the basic migration of particles based on the defined migration pattern.
131 *
132 * This function updates the positions of particles within the swarm by migrating them
133 * to target MPI ranks as specified in the migration list. It handles the migration process
134 * by setting the 'DMSwarm_rank' field for each particle and invokes the DMSwarm migration
135 * mechanism to relocate particles across MPI processes. After migration, it cleans up
136 * allocated resources and ensures synchronization across all MPI ranks.
137 *
138 * @param[in,out] user Pointer to the UserCtx structure containing simulation context.
139 *
140 * @return PetscErrorCode Returns 0 on success, non-zero on failure.
141 */
142PetscErrorCode PerformBasicMigration(UserCtx* user);
143
144/**
145 * @brief Identifies particles leaving the local bounding box and finds their target neighbor rank.
146 *
147 * Iterates local particles, checks against local bounding box. If outside, checks
148 * the pre-computed immediate neighbors (user->neighbors) using the global bboxlist
149 * to see if the particle landed in one of them. Populates the migrationList.
150 * Does NOT handle particles leaving the global domain (assumes CheckAndRemove was called).
151 *
152 * @param user Pointer to the UserCtx (contains local bbox and neighbors).
153 * @param bboxlist Array of BoundingBox structs for all ranks (for checking neighbor boxes).
154 * @param migrationList Pointer to an array of MigrationInfo structs (output, allocated/reallocated by this func).
155 * @param migrationCount Pointer to the number of particles marked for migration (output).
156 * @param listCapacity Pointer to the current allocated capacity of migrationList (in/out).
157 *
158 * @return PetscErrorCode 0 on success, non-zero on failure.
159 */
161 const BoundingBox *bboxlist,
162 MigrationInfo **migrationList,
163 PetscInt *migrationCount,
164 PetscInt *listCapacity);
165
166/**
167 * @brief Writes migration destinations into the DMSwarm rank field for marked particles.
168 *
169 * This helper consumes the migration list produced by `IdentifyMigratingParticles`
170 * and updates each selected particle's destination rank so that a subsequent
171 * `PerformMigration` call can transfer ownership correctly.
172 *
173 * @param[in,out] user Context containing the swarm and migration rank field.
174 * @param[in] migrationList Array of migration directives (local index + destination rank).
175 * @param[in] migrationCount Number of valid entries in `migrationList`.
176 * @return PetscErrorCode 0 on success.
177 */
178PetscErrorCode SetMigrationRanks(UserCtx* user, const MigrationInfo *migrationList, PetscInt migrationCount);
179
180/**
181 * @brief Performs particle migration based on the pre-populated DMSwarmPICField_rank field.
182 *
183 * Assumes SetMigrationRanks has already been called to mark particles with their target ranks.
184 * Calls DMSwarmMigrate to execute the communication and removal of un-migrated particles.
185 *
186 * @param user Pointer to the UserCtx structure containing the swarm.
187 *
188 * @return PetscErrorCode 0 on success, non-zero on failure.
189 */
190PetscErrorCode PerformMigration(UserCtx *user);
191
192/**
193 * @brief Counts particles in each cell of the DMDA 'da' and stores the result in user->ParticleCount.
194 *
195 * Assumes user->ParticleCount is a pre-allocated global vector associated with user->da
196 * and initialized to zero before calling this function (though it resets it internally).
197 * Assumes particle 'DMSwarm_CellID' field contains local cell indices.
198 *
199 * @param[in,out] user Pointer to the UserCtx structure containing da, swarm, and ParticleCount.
200 * @return PetscErrorCode Returns 0 on success, non-zero on failure.
201 */
202PetscErrorCode CalculateParticleCountPerCell(UserCtx *user);
203
204/**
205 * @brief Resizes a swarm collectively to a target global particle count.
206 *
207 * If particles are removed, the current implementation trims from the global tail
208 * ordering. If particles are added, new entries are created according to PETSc
209 * DMSwarm resize semantics.
210 *
211 * @param[in,out] swarm Swarm object to resize.
212 * @param[in] N_target Target global particle count.
213 * @return PetscErrorCode 0 on success.
214 */
215PetscErrorCode ResizeSwarmGlobally(DM swarm, PetscInt N_target);
216
217/**
218 * @brief Checks particle count in the reference file and resizes the swarm if needed.
219 *
220 * Reads the specified field file (e.g., position) into a temporary Vec to determine
221 * the number of particles (`N_file`) represented in that file for the given timestep.
222 * Compares `N_file` with the current swarm size (`N_current`). If they differ,
223 * resizes the swarm globally (adds or removes particles) to match `N_file`.
224 * Removal assumes excess particles are the globally last ones.
225 *
226 * @param[in,out] user Pointer to the UserCtx structure containing the DMSwarm.
227 * @param[in] ti Time index for constructing the file name.
228 * @param[in] ext File extension (e.g., "dat").
229 *
230 * @return PetscErrorCode 0 on success, non-zero on critical failure.
231 */
232PetscErrorCode PreCheckAndResizeSwarm(UserCtx *user, PetscInt ti, const char *ext);
233
234/**
235 * @brief Performs one full cycle of particle migration: identify, set ranks, and migrate.
236 *
237 * This function encapsulates the three main steps of migrating particles between MPI ranks:
238 * 1. Identify particles on the local rank that need to move based on their current
239 * positions and the domain decomposition (`bboxlist`).
240 * 2. Determine the destination rank for each migrating particle.
241 * 3. Perform the actual migration using PETSc's `DMSwarmMigrate`.
242 * It also calculates and logs the global number of particles migrated.
243 *
244 * @param user Pointer to the UserCtx structure.
245 * @param bboxlist Array of BoundingBox structures defining the spatial domain of each MPI rank.
246 * @param migrationList_p Pointer to a pointer for the MigrationInfo array. This array will be
247 * allocated/reallocated by `IdentifyMigratingParticles` if necessary.
248 * The caller is responsible for freeing this list eventually.
249 * @param migrationCount_p Pointer to store the number of particles identified for migration
250 * on the local rank. This is reset to 0 after migration for the current cycle.
251 * @param migrationListCapacity_p Pointer to store the current capacity of the `migrationList_p` array.
252 * @param currentTime Current simulation time (used for logging).
253 * @param step Current simulation step number (used for logging).
254 * @param migrationCycleName A descriptive name for this migration cycle (e.g., "Preliminary Sort", "Main Loop")
255 * for logging purposes.
256 * @param[out] globalMigrationCount_out Pointer to store the total number of particles migrated
257 * across all MPI ranks during this cycle.
258 * @return PetscErrorCode 0 on success, non-zero on failure.
259 */
260PetscErrorCode PerformSingleParticleMigrationCycle(UserCtx *user, const BoundingBox *bboxlist,
261 MigrationInfo **migrationList_p, PetscInt *migrationCount_p,
262 PetscInt *migrationListCapacity_p,
263 PetscReal currentTime, PetscInt step, const char *migrationCycleName,
264 PetscInt *globalMigrationCount_out);
265
266
267/**
268 * @brief Re-initializes the positions of particles currently on this rank if this rank owns
269 * part of the designated inlet surface.
270 *
271 * This function is intended for `user->ParticleInitialization == 0` (Surface Initialization mode)
272 * and is typically called after an initial migration step (e.g., in `PerformInitialSetup`).
273 * It ensures that all particles that should originate from the inlet surface and are now
274 * on the correct MPI rank are properly distributed across that rank's portion of the inlet.
275 *
276 * @param user Pointer to the UserCtx structure, containing simulation settings and grid information.
277 * @param currentTime Current simulation time (used for logging).
278 * @param step Current simulation step number (used for logging).
279 * @return PetscErrorCode 0 on success, non-zero on failure.
280 */
281PetscErrorCode ReinitializeParticlesOnInletSurface(UserCtx *user, PetscReal currentTime, PetscInt step);
282
283/**
284 * @brief Creates a sorted snapshot of all Particle IDs (PIDs) from a raw data array.
285 * @ingroup ParticleUtils
286 *
287 * This function is a crucial helper for the migration process. It captures the state of
288 * which particles are on the current MPI rank *before* migration occurs by taking a
289 * pointer to the swarm's raw PID data array. The resulting sorted array can then be used
290 * with an efficient binary search to quickly identify newcomer particles after migration.
291 *
292 * This function does NOT call DMSwarmGetField/RestoreField. It is the caller's
293 * responsibility to acquire the `pid_field` pointer before calling and restore it afterward.
294 *
295 * @param[in] pid_field A read-only pointer to the raw array of PIDs for the local swarm.
296 * @param[in] n_local The number of particles currently on the local rank.
297 * @param[out] pids_snapshot_out A pointer to a `PetscInt64*` array. This function will
298 * allocate memory for this array, and the caller is
299 * responsible for freeing it with `PetscFree()` when it
300 * is no longer needed.
301 *
302 * @return PetscErrorCode 0 on success, or a non-zero PETSc error code on failure.
303 */
304PetscErrorCode GetLocalPIDSnapshot(const PetscInt64 pid_field[],
305 PetscInt n_local,
306 PetscInt64 **pids_snapshot_out);
307
308/**
309 * @brief Safely adds a new migration task to a dynamically sized list.
310 *
311 * This utility function manages a dynamic array of MigrationInfo structs. It appends
312 * a new entry to the list and automatically doubles the array's capacity using
313 * `PetscRealloc` if the current capacity is exceeded. This prevents buffer overflows
314 * and avoids the need to know the number of migrating particles in advance.
315 *
316 * @param[in,out] migration_list_p A pointer to the MigrationInfo array pointer. The function
317 * will update this pointer if the array is reallocated.
318 * @param[in,out] capacity_p A pointer to an integer holding the current allocated
319 * capacity of the list (in number of elements). This will be
320 * updated upon reallocation.
321 * @param[in,out] count_p A pointer to an integer holding the current number of
322 * items in the list. This will be incremented by one.
323 * @param[in] particle_local_idx The local index (from 0 to nlocal-1) of the particle
324 * that needs to be migrated.
325 * @param[in] destination_rank The target MPI rank for the particle.
326 *
327 * @return PetscErrorCode 0 on success, or a non-zero PETSc error code on failure (e.g., from memory allocation).
328 */
329PetscErrorCode AddToMigrationList(MigrationInfo **migration_list_p,
330 PetscInt *capacity_p,
331 PetscInt *count_p,
332 PetscInt particle_local_idx,
333 PetscMPIInt destination_rank);
334
335/**
336 * @brief Identifies newly arrived particles after migration and flags them for a location search.
337 * @ingroup ParticleMotion
338 *
339 * This function is a critical component of the iterative migration process managed by
340 * the main particle settlement orchestrator (e.g., `SettleParticles`). After a
341 * `DMSwarmMigrate` call, each rank's local particle list is a new mix of resident
342 * particles and newly received ones. This function's job is to efficiently identify
343 * these "newcomers" and set their `DMSwarm_location_status` field to `NEEDS_LOCATION`.
344 *
345 * This ensures that in the subsequent pass of the migration `do-while` loop, only the
346 * newly arrived particles are processed by the expensive location algorithm, preventing
347 * redundant work on particles that are already settled on the current rank.
348 *
349 * The identification is done by comparing the PIDs of particles currently on the rank
350 * against a "snapshot" of PIDs taken *before* the migration occurred.
351 *
352 * @param[in] swarm The DMSwarm object, which has just completed a migration.
353 * @param[in] n_local_before The number of particles that were on this rank *before* the
354 * migration was performed.
355 * @param[in] pids_before A pre-sorted array of the PIDs that were on this rank before
356 * the migration. This is used for fast lookups.
357 *
358 * @return PetscErrorCode 0 on success, or a non-zero PETSc error code on failure.
359 *
360 * @note This function assumes the `pids_before` array is sorted in ascending order to
361 * enable the use of an efficient binary search.
362 */
363PetscErrorCode FlagNewcomersForLocation(DM swarm,
364 PetscInt n_local_before,
365 const PetscInt64 pids_before[]);
366
367
368/**
369 * @brief Fast-path migration for restart particles using preloaded Cell IDs.
370 *
371 * This function provides an optimized migration path specifically for particles
372 * loaded from restart files. Unlike the standard `LocateAllParticlesInGrid()`
373 * which performs expensive walking searches, this function leverages the fact that
374 * restart particles already have valid global Cell IDs loaded from disk.
375 *
376 * **How It Works:**
377 * 1. Iterates through all local particles.
378 * 2. For each particle with a valid Cell ID (ci, cj, ck):
379 * - Calls `FindOwnerOfCell(ci, cj, ck)` to determine the correct rank.
380 * - If owner differs from current rank, adds to migration list.
381 * - If owner matches current rank, the existing `ACTIVE_AND_LOCATED` status is preserved.
382 * 3. Uses existing `SetMigrationRanks()` and `PerformMigration()` infrastructure.
383 * 4. Achieves **single-pass direct migration** (no multi-hop, no walking searches).
384 *
385 * @param[in,out] user Pointer to UserCtx containing the swarm and RankCellInfoMap.
386 * The function updates particle status fields and performs migration.
387 *
388 * @return PetscErrorCode 0 on success, non-zero on failure.
389 *
390 * @note Testing status:
391 * Direct coverage currently focuses on restart fast-path ownership transfer.
392 * Non-restart multi-pass migration behavior remains part of the next
393 * simulation-core test backlog.
394 */
395PetscErrorCode MigrateRestartParticlesUsingCellID(UserCtx *user);
396
397/**
398 * @brief Orchestrates the complete particle location and migration process for one timestep.
399 * @ingroup ParticleLocation
400 *
401 * This function is the master orchestrator for ensuring every particle is on its correct
402 * MPI rank and has a valid host cell index. It is designed to be called once per
403 * timestep after particle positions have been updated.
404 *
405 * The function uses a robust, iterative "Guess and Verify" strategy within a
406 * do-while loop to handle complex particle motion across processor boundaries,
407 * especially on curvilinear grids.
408 *
409 * 1. **State Snapshot:** At the start of each pass, it captures a list of all Particle IDs (PIDs)
410 * on the current rank.
411 * 2. **"Guess" (Heuristic):** For particles that are "lost" (no valid host cell),
412 * it first attempts a fast, bounding-box-based guess to find a potential new owner rank.
413 * 3. **"Verify" (Robust Walk):** For all other particles, or if the guess fails,
414 * it uses a robust cell-walking algorithm (`LocateParticleOrFindMigrationTarget`)
415 * that determines the particle's status: located locally, needs migration, or is lost.
416 * 4. **Migration:** After identifying all migrating particles on a pass, it performs the
417 * MPI communication using the `SetMigrationRanks` and `PerformMigration` helpers.
418 * 5. **Newcomer Flagging:** After migration, it uses the PID snapshot from step 1 to
419 * efficiently identify newly arrived particles and flag them for location on the next pass.
420 * 6. **Iteration:** The process repeats in a `do-while` loop until a pass occurs where
421 * no particles migrate, ensuring the entire swarm is in a stable, consistent state.
422 *
423 * @param[in,out] user Pointer to the UserCtx, containing the swarm and all necessary
424 * domain topology information (bboxlist, RankCellInfoMap, etc.).
425 * @param[in] bboxlist An array of BoundingBox structures for ALL MPI ranks, indexed 0 to (size-1).
426 * This array must be up-to-date and available on all ranks.
427 * @return PetscErrorCode 0 on success, or a non-zero PETSc error code on failure.
428 *
429 * @note Testing status:
430 * Direct unit coverage currently pins the prior-cell fast path and the
431 * local guess-then-verify path. Multi-pass migration, newcomer flagging,
432 * and several lost/migration edge cases are still targeted for future
433 * bespoke tests.
434 */
435PetscErrorCode LocateAllParticlesInGrid(UserCtx *user,BoundingBox *bboxlist);
436
437/**
438 * @brief Marks all local particles as `NEEDS_LOCATION` for the next settlement pass.
439 *
440 * This function is designed to be called at the end of a full timestep, after all
441 * particle-based calculations are complete. It prepares the swarm for the next
442 * timestep by ensuring that after the next position update, every particle will be
443 * re-evaluated by the LocateAllParticlesInGrid orchestrator.
444 *
445 * It iterates through all locally owned particles and sets their
446 * `DMSwarm_location_status` field to `NEEDS_LOCATION`.
447 *
448 * @param[in,out] user Pointer to the UserCtx containing the swarm.
449 * @return PetscErrorCode 0 on success, or a non-zero PETSc error code on failure.
450 */
451PetscErrorCode ResetAllParticleStatuses(UserCtx *user);
452
453 #endif // PARTICLE_MOTION_H
PetscErrorCode GenerateGaussianNoise(PetscRandom rnd, PetscReal *n1, PetscReal *n2)
Generates two independent standard normal random variables N(0,1) using the Box-Muller transform.
PetscErrorCode ResizeSwarmGlobally(DM swarm, PetscInt N_target)
Resizes a swarm collectively to a target global particle count.
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)
Writes migration destinations into the DMSwarm rank field for marked particles.
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 MigrateRestartParticlesUsingCellID(UserCtx *user)
Fast-path migration for restart particles using preloaded Cell IDs.
PetscErrorCode UpdateAllParticlePositions(UserCtx *user)
Loops over all local particles in the DMSwarm, updating their positions based on velocity and the glo...
PetscErrorCode CalculateParticleCountPerCell(UserCtx *user)
Counts particles in each cell of the DMDA 'da' and stores the result in user->ParticleCount.
PetscErrorCode CalculateBrownianDisplacement(UserCtx *user, PetscReal diff_eff, Cmpnts *displacement)
Calculates the stochastic displacement vector (Brownian motion) for a single particle.
PetscErrorCode LocateAllParticlesInGrid(UserCtx *user, BoundingBox *bboxlist)
Orchestrates the complete particle location and migration process for one timestep.
PetscErrorCode UpdateParticlePosition(UserCtx *user, Particle *particle)
Updates a particle's position based on its velocity and the timestep dt (stored in user->dt).
PetscErrorCode ResetAllParticleStatuses(UserCtx *user)
Marks all local particles as NEEDS_LOCATION for the next settlement pass.
PetscErrorCode DefineBasicMigrationPattern(UserCtx *user)
Defines the basic migration pattern for particles within the swarm.
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.
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:154
A 3D point or vector with PetscScalar components.
Definition variables.h:100
Information needed to migrate a single particle between MPI ranks.
Definition variables.h:192
Defines a particle's core properties for Lagrangian tracking.
Definition variables.h:165
User-defined context containing data specific to a single computational grid level.
Definition variables.h:811
Header file for particle location functions using the walking search algorithm.