PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
Loading...
Searching...
No Matches
ParticleSwarm.h
Go to the documentation of this file.
1/**
2 * @file ParticleSwarm.h
3 * @brief Header file for Particle Swarm management functions.
4 *
5 * This file contains declarations of functions responsible for creating, managing,
6 * initializing, migrating, and printing particle swarms within a simulation using PETSc's DMSwarm.
7 */
8
9#ifndef PARTICLE_SWARM_H
10#define PARTICLE_SWARM_H
11
12// Include necessary headers
13#include <petsc.h> // PETSc library header
14#include <petscdmswarm.h> // PETSc DMSwarm header
15#include <stdbool.h>
16#include <math.h>
17#include "variables.h" // Common type definitions
18#include "logging.h" // Logging macros and definitions
19#include "walkingsearch.h"
20#include "Metric.h"
21#include "io.h"
22// --------------------- Function Declarations ---------------------
23
24/**
25 * @brief Creates and initializes a Particle Swarm.
26 *
27 * This function sets up a DMSwarm within the provided UserCtx structure, initializes
28 * particle fields, and distributes particles across MPI processes. It ensures that
29 * the number of particles is evenly divided among the available MPI ranks. If the total
30 * number of particles isn't divisible by the number of processes, the remainder is distributed
31 * to the first few ranks.
32 *
33 * Additionally, it now takes a 'bboxlist' array as an input parameter and passes it on to
34 * AssignInitialProperties(), enabling particle initialization at the midpoint of each rank's
35 * bounding box if ParticleInitialization is set to 0.
36 *
37 * @param[in,out] user Pointer to the UserCtx structure containing the simulation context.
38 * @param[in] numParticles Total number of particles to create across all MPI processes.
39 * @param[in] bboxlist Pointer to an array of BoundingBox structures, one per rank.
40 *
41 * @return PetscErrorCode Returns 0 on success, non-zero on failure.
42 *
43 * @note
44 * - Ensure that `numParticles` is a positive integer.
45 * - The `control.dat` file should contain necessary PETSc options.
46 * - The `bboxlist` array should be properly populated before calling this function.
47 */
48PetscErrorCode CreateParticleSwarm(UserCtx *user, PetscInt numParticles, PetscInt *particlesPerProcess, BoundingBox *bboxlist);
49
50/**
51 * @brief Initializes the DMSwarm object within the UserCtx structure.
52 *
53 * This function creates the DMSwarm, sets its type and dimension, and configures basic swarm properties.
54 *
55 * @param[in,out] user Pointer to the UserCtx structure containing simulation context.
56 *
57 * @return PetscErrorCode Returns 0 on success, non-zero on failure.
58 */
59PetscErrorCode InitializeSwarm(UserCtx* user);
60
61/**
62 * @brief Registers a swarm field without finalizing registration.
63 *
64 * This function calls DMSwarmRegisterPetscDatatypeField for the given field,
65 * but does not finalize the registration. The finalization is deferred until
66 * all fields have been registered.
67 *
68 * @param swarm [in] The DMSwarm object.
69 * @param fieldName [in] Name of the field to register.
70 * @param fieldDim [in] Dimension of the field (1 for scalar, 3 for vector, etc.).
71 * @param dtype [in] The datatype of the swarm field being registered.
72 * @return PetscErrorCode Returns 0 on success, non-zero on failure.
73 */
74PetscErrorCode RegisterSwarmField(DM swarm, const char *fieldName, PetscInt fieldDim, PetscDataType dtype);
75
76/**
77 * @brief Registers necessary particle fields within the DMSwarm.
78 *
79 * This function registers fields such as position, velocity, CellID, and weight for each particle.
80 *
81 * @param[in,out] swarm The DMSwarm object managing the particle swarm.
82 *
83 * @return PetscErrorCode Returns 0 on success, non-zero on failure.
84 */
85PetscErrorCode RegisterParticleFields(DM swarm);
86
87/**
88 * @brief Initializes all particle properties in the swarm.
89 *
90 * This function orchestrates the initialization of particle properties.
91 * It first determines the inlet face if surface initialization (Mode 0) is selected
92 * by parsing "bcs.dat".
93 * Then, it initializes basic particle properties (physical position, Particle ID,
94 * and placeholder Cell IDs) by calling `InitializeParticleBasicProperties`. This call
95 * uses the provided `rand_logic_i/j/k` RNGs, which must be pre-initialized for [0,1).
96 * The `rand_phys_x/y/z` RNGs (physically bounded) are passed but may not be used by
97 * `InitializeParticleBasicProperties` for position setting if all initialization paths
98 * use logical-to-physical mapping.
99 * Finally, it calls helper functions to initialize other registered swarm fields
100 * like "velocity", "weight", and "P" (pressure) to default values.
101 *
102 * @param[in,out] user Pointer to the `UserCtx` structure.
103 * @param[in] particlesPerProcess Number of particles assigned to this MPI process.
104 * @param[in] rand_phys_x RNG for physical x-coordinates (from `InitializeRandomGenerators`).
105 * @param[in] rand_phys_y RNG for physical y-coordinates (from `InitializeRandomGenerators`).
106 * @param[in] rand_phys_z RNG for physical z-coordinates (from `InitializeRandomGenerators`).
107 * @param[in] rand_logic_i RNG for i-logical dimension tasks [0,1) (from `InitializeLogicalSpaceRNGs`).
108 * @param[in] rand_logic_j RNG for j-logical dimension tasks [0,1) (from `InitializeLogicalSpaceRNGs`).
109 * @param[in] rand_logic_k RNG for k-logical dimension tasks [0,1) (from `InitializeLogicalSpaceRNGs`).
110 * @param[in] bboxlist Array of BoundingBox structures (potentially unused by IPBP).
111 *
112 * @return PetscErrorCode Returns 0 on success, non-zero on failure.
113 */
114PetscErrorCode AssignInitialPropertiesToSwarm(UserCtx* user,
115 PetscInt particlesPerProcess,
116 PetscRandom *rand_phys_x, // RNG from original InitializeRandomGenerators
117 PetscRandom *rand_phys_y, // RNG from original InitializeRandomGenerators
118 PetscRandom *rand_phys_z, // RNG from original InitializeRandomGenerators
119 PetscRandom *rand_logic_i, // RNG from InitializeLogicalSpaceRNGs
120 PetscRandom *rand_logic_j, // RNG from InitializeLogicalSpaceRNGs
121 PetscRandom *rand_logic_k, // RNG from InitializeLogicalSpaceRNGs
122 BoundingBox *bboxlist);
123
124/**
125 * @brief Distributes particles evenly across MPI processes, handling any remainders.
126 *
127 * This function calculates the number of particles each MPI process should handle,
128 * distributing the remainder particles to the first few ranks if necessary.
129 *
130 * @param[in] numParticles Total number of particles to create across all MPI processes.
131 * @param[in] rank MPI rank of the current process.
132 * @param[in] size Total number of MPI processes.
133 * @param[out] particlesPerProcess Number of particles assigned to the current MPI process.
134 * @param[out] remainder Remainder particles when dividing numParticles by size.
135 *
136 * @return PetscErrorCode Returns 0 on success, non-zero on failure.
137 */
138PetscErrorCode DistributeParticles(PetscInt numParticles, PetscMPIInt rank, PetscMPIInt size, PetscInt* particlesPerProcess, PetscInt* remainder);
139
140/**
141 * @brief Finalizes the swarm setup by destroying random generators and logging completion.
142 *
143 * This function cleans up resources by destroying random number generators and LOG_ALLOWs the completion of swarm setup.
144 *
145 * @param[in] randx Random number generator for the x-coordinate.
146 * @param[in] randy Random number generator for the y-coordinate.
147 * @param[in] randz Random number generator for the z-coordinate.
148 * @param[in] rand_logic_i Random number generator for the xi-coordinate.
149 * @param[in] rand_logic_j Random number generator for the eta-coordinate.
150 * @param[in] rand_logic_k Random number generator for the zeta-coordinate.
151 *
152 * @return PetscErrorCode Returns 0 on success, non-zero on failure.
153 */
154PetscErrorCode FinalizeSwarmSetup(PetscRandom *randx, PetscRandom *randy, PetscRandom *randz, PetscRandom *rand_logic_i, PetscRandom *rand_logic_j, PetscRandom *rand_logic_k);
155
156/**
157 * @brief Initializes a Particle struct with data from DMSwarm fields.
158 *
159 * This helper function populates a Particle structure using data retrieved from DMSwarm fields.
160 *
161 * @param[in] i Index of the particle in the DMSwarm.
162 * @param[in] PIDs Pointer to the array of particle IDs.
163 * @param[in] weights Pointer to the array of particle weights.
164 * @param[in] positions Pointer to the array of particle positions.
165 * @param[in] cellIndices Pointer to the array of particle cell indices.
166 * @param[in] velocities Pointer to the array of particle velocities.
167 * @param[in] LocStatus Pointer to the array of cell location status indicators.
168 * @param[out] particle Pointer to the Particle struct to initialize.
169 *
170 * @return PetscErrorCode Returns `0` on success, non-zero on failure.
171 */
172PetscErrorCode UnpackSwarmFields(PetscInt i, const PetscInt64 *PIDs, const PetscReal *weights,
173 const PetscReal *positions, const PetscInt *cellIndices,
174 PetscReal *velocities,PetscInt *LocStatus,Particle *particle);
175
176/**
177 * @brief Updates DMSwarm fields with data from a Particle struct.
178 *
179 * This helper function writes back the modified Particle data to the corresponding DMSwarm fields.
180 *
181 * @param[in] i Index of the particle in the DMSwarm.
182 * @param[in] particle Pointer to the Particle struct containing updated data.
183 * @param[in,out] weights Pointer to the array of particle weights.
184 * @param[in,out] cellIndices Pointer to the array of particle cell indices.
185 * @param[in,out] LocStatus Pointer to the array of cell location status indicators.
186 *
187 * @return PetscErrorCode Returns `0` on success, non-zero on failure.
188 */
189PetscErrorCode UpdateSwarmFields(PetscInt i, const Particle *particle,
190 PetscReal *weights, PetscInt *cellIndices, PetscInt *status_field);
191
192/**
193 * @brief Checks if a particle's location is within a specified bounding box.
194 *
195 * This function determines whether the given particle's location lies inside the provided bounding box.
196 * It performs an axis-aligned bounding box (AABB) check by comparing the particle's coordinates to the
197 * minimum and maximum coordinates of the bounding box in each dimension (x, y, z).
198 *
199 * Logging statements are included to provide detailed information about the function's execution.
200 *
201 * @param[in] bbox Pointer to the BoundingBox structure containing minimum and maximum coordinates.
202 * @param[in] particle Pointer to the Particle structure containing the particle's location and identifier.
203 *
204 * @return PetscBool Returns `PETSC_TRUE` if the particle is inside the bounding box, `PETSC_FALSE` otherwise.
205 *
206 * @note
207 * - The function assumes that the `bbox` and `particle` pointers are valid and non-NULL.
208 * - The function includes logging statements that start with the function name.
209 * - Be cautious when logging in performance-critical code sections, especially if the function is called frequently.
210 */
211PetscBool IsParticleInsideBoundingBox(const BoundingBox *bbox, const Particle *particle);
212
213/**
214 * @brief Updates a particle's interpolation weights based on distances to cell faces.
215 *
216 * This function computes interpolation weights using distances to the six
217 * cell faces (`d`) and updates the `weight` field of the provided particle.
218 *
219 * @param[in] d Pointer to an array of distances to the six cell faces.
220 * @param[out] particle Pointer to the Particle structure whose weights are to be updated.
221 *
222 * @return PetscErrorCode Returns 0 on success, or a non-zero error code on failure.
223 */
224PetscErrorCode UpdateParticleWeights(PetscReal *d, Particle *particle);
225
226/**
227 * @brief Perform particle swarm initialization, particle-grid interaction, and related operations.
228 *
229 * This function handles the following tasks:
230 * 1. Initializes the particle swarm using the provided bounding box list (bboxlist) to determine initial placement
231 * if ParticleInitialization is 0.
232 * 2. Locates particles within the computational grid.
233 * 3. Updates particle positions based on grid interactions (if such logic exists elsewhere in the code).
234 * 4. Interpolates particle velocities from grid points using trilinear interpolation.
235 *
236 * @param[in,out] user Pointer to the UserCtx structure containing grid and particle swarm information.
237 * @param[in] np Number of particles to initialize in the swarm.
238 * @param[in] bboxlist Pointer to an array of BoundingBox structures, one per MPI rank.
239 *
240 * @return PetscErrorCode Returns 0 on success, non-zero on failure.
241 *
242 * @note
243 * - Ensure that `np` (number of particles) is positive.
244 * - The `bboxlist` array must be correctly computed and passed in before calling this function.
245 * - If ParticleInitialization == 0, particles will be placed at the midpoint of the local bounding box.
246 */
247PetscErrorCode InitializeParticleSwarm(SimCtx *simCtx);
248
249#endif // PARTICLE_SWARM_H
PetscErrorCode UpdateParticleWeights(PetscReal *d, Particle *particle)
Updates a particle's interpolation weights based on distances to cell faces.
PetscErrorCode CreateParticleSwarm(UserCtx *user, PetscInt numParticles, PetscInt *particlesPerProcess, BoundingBox *bboxlist)
Creates and initializes a Particle Swarm.
PetscErrorCode FinalizeSwarmSetup(PetscRandom *randx, PetscRandom *randy, PetscRandom *randz, PetscRandom *rand_logic_i, PetscRandom *rand_logic_j, PetscRandom *rand_logic_k)
Finalizes the swarm setup by destroying random generators and logging completion.
PetscErrorCode DistributeParticles(PetscInt numParticles, PetscMPIInt rank, PetscMPIInt size, PetscInt *particlesPerProcess, PetscInt *remainder)
Distributes particles evenly across MPI processes, handling any remainders.
PetscBool IsParticleInsideBoundingBox(const BoundingBox *bbox, const Particle *particle)
Checks if a particle's location is within a specified bounding box.
PetscErrorCode UnpackSwarmFields(PetscInt i, const PetscInt64 *PIDs, const PetscReal *weights, const PetscReal *positions, const PetscInt *cellIndices, PetscReal *velocities, PetscInt *LocStatus, Particle *particle)
Initializes a Particle struct with data from DMSwarm fields.
PetscErrorCode InitializeSwarm(UserCtx *user)
Initializes the DMSwarm object within the UserCtx structure.
PetscErrorCode UpdateSwarmFields(PetscInt i, const Particle *particle, PetscReal *weights, PetscInt *cellIndices, PetscInt *status_field)
Updates DMSwarm fields with data from a Particle struct.
PetscErrorCode InitializeParticleSwarm(SimCtx *simCtx)
Perform particle swarm initialization, particle-grid interaction, and related operations.
PetscErrorCode RegisterSwarmField(DM swarm, const char *fieldName, PetscInt fieldDim, PetscDataType dtype)
Registers a swarm field without finalizing registration.
PetscErrorCode RegisterParticleFields(DM swarm)
Registers necessary particle fields within the DMSwarm.
PetscErrorCode AssignInitialPropertiesToSwarm(UserCtx *user, PetscInt particlesPerProcess, PetscRandom *rand_phys_x, PetscRandom *rand_phys_y, PetscRandom *rand_phys_z, PetscRandom *rand_logic_i, PetscRandom *rand_logic_j, PetscRandom *rand_logic_k, BoundingBox *bboxlist)
Initializes all particle properties in the swarm.
Public interface for data input/output routines.
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
Defines a particle's core properties for Lagrangian tracking.
Definition variables.h:165
The master context for the entire simulation.
Definition variables.h:538
User-defined context containing data specific to a single computational grid level.
Definition variables.h:661
Header file for particle location functions using the walking search algorithm.