PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
Loading...
Searching...
No Matches
ParticleSwarm.c
Go to the documentation of this file.
1 // ParticleSwarm.c
2
3#include "ParticleSwarm.h"
4
5#define INTERPOLATION_DISTANCE_TOLERANCE 1.0e-14
6
7#undef __FUNCT__
8#define __FUNCT__ "InitializeSwarm"
9/**
10 * @brief Initializes the DMSwarm object within the UserCtx structure.
11 *
12 * This function creates the DMSwarm, sets its type and dimension, and configures basic swarm properties.
13 *
14 * @param[in,out] user Pointer to the UserCtx structure containing simulation context.
15 *
16 * @return PetscErrorCode Returns 0 on success, non-zero on failure.
17 */
18PetscErrorCode InitializeSwarm(UserCtx* user) {
19 PetscErrorCode ierr; // Error code for PETSc functions
20
21 PetscFunctionBeginUser;
23 // Create the DMSwarm object for particle management
24 ierr = DMCreate(PETSC_COMM_WORLD, &user->swarm); CHKERRQ(ierr);
25 ierr = DMSetType(user->swarm, DMSWARM); CHKERRQ(ierr);
26 ierr = DMSetDimension(user->swarm, 3); CHKERRQ(ierr);
27 ierr = DMSwarmSetType(user->swarm, DMSWARM_BASIC); CHKERRQ(ierr);
28 LOG_ALLOW(LOCAL,LOG_INFO, "DMSwarm created and configured.\n");
29
31 PetscFunctionReturn(0);
32}
33
34#undef __FUNCT__
35#define __FUNCT__ "RegisterSwarmField"
36
37/**
38 * @brief Registers a swarm field without finalizing registration.
39 *
40 * This function calls DMSwarmRegisterPetscDatatypeField for the given field,
41 * but does not finalize the registration. The finalization is deferred until
42 * all fields have been registered.
43 *
44 * @param swarm [in] The DMSwarm object.
45 * @param fieldName [in] Name of the field to register.
46 * @param fieldDim [in] Dimension of the field (1 for scalar, 3 for vector, etc.).
47 * @param dtype [in] The datatype of the swarm field being registered.
48 * @return PetscErrorCode Returns 0 on success, non-zero on failure.
49 */
50PetscErrorCode RegisterSwarmField(DM swarm, const char *fieldName, PetscInt fieldDim, PetscDataType dtype)
51{
52 PetscErrorCode ierr;
53 PetscFunctionBeginUser;
54
55 ierr = DMSwarmRegisterPetscDatatypeField(swarm, fieldName, fieldDim, dtype); CHKERRQ(ierr);
56 // PetscDataTypes is an extern char* [] defined in petscsystypes.h that gives string names for PetscDataType enums
57 LOG_ALLOW(LOCAL,LOG_DEBUG,"Registered field '%s' with dimension=%d, type=%s.\n",
58 fieldName, fieldDim, PetscDataTypes[dtype]);
59
60 PetscFunctionReturn(0);
61}
62
63#undef __FUNCT__
64#define __FUNCT__ "RegisterParticleFields"
65
66/**
67 * @brief Registers necessary particle fields within the DMSwarm.
68 *
69 * This function registers fields such as position, velocity, CellID, and weight for each particle.
70 *
71 * @param[in,out] swarm The DMSwarm object managing the particle swarm.
72 *
73 * @return PetscErrorCode Returns 0 on success, non-zero on failure.
74 */
75
76PetscErrorCode RegisterParticleFields(DM swarm)
77{
78 PetscErrorCode ierr;
79 PetscFunctionBeginUser;
80
81 // Register each field using the helper function
82 ierr = RegisterSwarmField(swarm, "position", 3 ,PETSC_REAL); CHKERRQ(ierr);
83 LOG_ALLOW(LOCAL,LOG_DEBUG, "Registered field 'position'.\n");
84
85 ierr = RegisterSwarmField(swarm, "velocity", 3, PETSC_REAL); CHKERRQ(ierr);
86 LOG_ALLOW(LOCAL,LOG_DEBUG,"Registered field 'velocity'.\n");
87
88 ierr = RegisterSwarmField(swarm, "DMSwarm_CellID", 3, PETSC_INT); CHKERRQ(ierr);
89 LOG_ALLOW(LOCAL,LOG_DEBUG,"Registered field 'DMSwarm_CellID'.\n");
90
91 ierr = RegisterSwarmField(swarm, "weight", 3,PETSC_REAL); CHKERRQ(ierr);
92 LOG_ALLOW(LOCAL,LOG_DEBUG,"Registered field 'weight'.\n");
93
94 ierr = RegisterSwarmField(swarm,"Diffusivity", 1,PETSC_REAL); CHKERRQ(ierr);
95 LOG_ALLOW(LOCAL,LOG_DEBUG,"Registered field 'Diffusivity' - Scalar.\n");
96
97 ierr = RegisterSwarmField(swarm,"DiffusivityGradient", 3,PETSC_REAL); CHKERRQ(ierr);
98 LOG_ALLOW(LOCAL,LOG_DEBUG,"Registered field 'DiffusivityGradient' - Vector.\n");
99
100 ierr = RegisterSwarmField(swarm,"Psi", 1,PETSC_REAL); CHKERRQ(ierr);
101 LOG_ALLOW(LOCAL,LOG_DEBUG,"Registered field 'Psi' - Scalar.\n");
102
103 ierr = RegisterSwarmField(swarm,"DMSwarm_location_status",1,PETSC_INT);CHKERRQ(ierr);
104 LOG_ALLOW(LOCAL,LOG_DEBUG,"Registered field 'DMSwarm_location_status' - Status of Location of Particle(located,lost etc).\n");
105
106 // Finalize the field registration after all fields have been added
107 ierr = DMSwarmFinalizeFieldRegister(swarm); CHKERRQ(ierr);
108 LOG_ALLOW(LOCAL,LOG_INFO,"RegisterParticleFields - Finalized field registration.\n");
109
110 PetscFunctionReturn(0);
111}
112
113#undef __FUNCT__
114#define __FUNCT__ "DetermineVolumetricInitializationParameters"
115/**
116 * @brief Determines cell selection and intra-cell logical coordinates for volumetric initialization (Mode 1).
117 *
118 * This function is called when `simCtx->ParticleInitialization == 1`. It randomly selects
119 * an *owned cell* on the current MPI rank and then generates random intra-cell logical
120 * coordinates `[0,1)^3` within that chosen cell.
121 *
122 * The process involves:
123 * 1. Checking if the current MPI rank owns any 3D cells.
124 * 2. If it does, it randomly selects an *owned cell index* in each logical direction (i, j, k)
125 * by scaling a `[0,1)` random number with the number of owned cells in that direction.
126 * 3. These local owned cell indices are then converted to the *local node indices*
127 * (`ci/cj/ck_metric_lnode_out`) corresponding to the origin of the selected cell,
128 * for use with `MetricLogicalToPhysical`. This conversion uses `xs/ys/zs_gnode`.
129 * 4. All three intra-cell logical coordinates (`xi/eta/zta_metric_logic_out`) for
130 * `MetricLogicalToPhysical` are chosen randomly within the `[0,1)` range.
131 * 5. A flag (`can_place_in_volume_out`) indicates if a valid placement could be determined.
132 *
133 * **Important Note on `DMDALocalInfo info` members (same as for surface init):**
134 * - `info->xs, info->ys, info->ks`: Global starting indices of *owned cells*.
135 * - `info->mx, info->my, info->mz`: Number of *grid points (nodes)* in each local dimension on this process.
136 * Therefore, the number of *owned cells* in a dimension is `info->mX - 1` (if `info->mX > 0`).
137 *
138 * @param[in] user Pointer to `UserCtx`. (Currently not used in this specific helper, but kept for API consistency).
139 * @param[in] info Pointer to `DMDALocalInfo` for the current rank's grid portion.
140 * @param[in] xs_gnode, ys_gnode, zs_gnode Local indices (in the ghosted array) of the first *owned node*.
141 * @param[in] rand_logic_i_ptr Pointer to the RNG for i-dimension tasks [0,1).
142 * @param[in] rand_logic_j_ptr Pointer to the RNG for j-dimension tasks [0,1).
143 * @param[in] rand_logic_k_ptr Pointer to the RNG for k-dimension tasks [0,1).
144 * @param[out] ci_metric_lnode_out Pointer to store the local i-node index of the selected cell's origin.
145 * @param[out] cj_metric_lnode_out Pointer to store the local j-node index of the selected cell's origin.
146 * @param[out] ck_metric_lnode_out Pointer to store the local k-node index of the selected cell's origin.
147 * @param[out] xi_metric_logic_out Pointer to store the intra-cell logical xi-coordinate [0,1).
148 * @param[out] eta_metric_logic_out Pointer to store the intra-cell logical eta-coordinate [0,1).
149 * @param[out] zta_metric_logic_out Pointer to store the intra-cell logical zeta-coordinate [0,1).
150 * @param[out] can_place_in_volume_out PETSC_TRUE if placement parameters were successfully determined, PETSC_FALSE otherwise.
151 * @return PetscErrorCode 0 on success, or a PETSc error code.
152 */
154 UserCtx *user, DMDALocalInfo *info,
155 PetscInt xs_gnode, PetscInt ys_gnode, PetscInt zs_gnode,
156 PetscRandom *rand_logic_i_ptr, PetscRandom *rand_logic_j_ptr, PetscRandom *rand_logic_k_ptr, /* Pointers to RNGs */
157 PetscInt *ci_metric_lnode_out, PetscInt *cj_metric_lnode_out, PetscInt *ck_metric_lnode_out,
158 PetscReal *xi_metric_logic_out, PetscReal *eta_metric_logic_out, PetscReal *zta_metric_logic_out,
159 PetscBool *can_place_in_volume_out)
160{
161 PetscErrorCode ierr = 0;
162 PetscReal r_val; // Temporary for random numbers from [0,1) RNGs
163 PetscInt local_owned_cell_idx_i, local_owned_cell_idx_j, local_owned_cell_idx_k;
164 PetscMPIInt rank_for_logging; // For logging if needed
165
166 PetscFunctionBeginUser;
167
169
170 ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank_for_logging); CHKERRQ(ierr);
171
172 *can_place_in_volume_out = PETSC_FALSE; // Default to: cannot place
173
174 // Default intra-cell logicals and cell node indices (e.g. if placement fails)
175 *xi_metric_logic_out = 0.5; *eta_metric_logic_out = 0.5; *zta_metric_logic_out = 0.5;
176 *ci_metric_lnode_out = xs_gnode; *cj_metric_lnode_out = ys_gnode; *ck_metric_lnode_out = zs_gnode;
177
178 // Calculate number of owned cells in each direction from node counts in info
179 // Get number of cells this rank owns in each dimension (tangential to the face mainly)
180 PetscInt owned_start_cell_i, num_owned_cells_on_rank_i;
181 PetscInt owned_start_cell_j, num_owned_cells_on_rank_j;
182 PetscInt owned_start_cell_k, num_owned_cells_on_rank_k;
183
184 ierr = GetOwnedCellRange(info, 0, &owned_start_cell_i, &num_owned_cells_on_rank_i); CHKERRQ(ierr);
185 ierr = GetOwnedCellRange(info, 1, &owned_start_cell_j, &num_owned_cells_on_rank_j); CHKERRQ(ierr);
186 ierr = GetOwnedCellRange(info, 2, &owned_start_cell_k, &num_owned_cells_on_rank_k); CHKERRQ(ierr);
187
188 if (num_owned_cells_on_rank_i > 0 && num_owned_cells_on_rank_j > 0 && num_owned_cells_on_rank_k > 0) { // If rank owns any 3D cells
189 *can_place_in_volume_out = PETSC_TRUE;
190
191 // --- 1. Select a Random Owned Cell ---
192 // The selected index will be a 0-based index relative to the start of this rank's owned cells.
193
194 // Select random local owned cell index in I-direction
195 ierr = PetscRandomGetValueReal(*rand_logic_i_ptr, &r_val); CHKERRQ(ierr); // Dereference RNG pointer
196 local_owned_cell_idx_i = (PetscInt)(r_val * num_owned_cells_on_rank_i);
197 // Clamp to be safe: local_owned_cell_idx_i should be in [0, num_owned_cells_on_rank_i - 1]
198 local_owned_cell_idx_i = PetscMin(PetscMax(0, local_owned_cell_idx_i), num_owned_cells_on_rank_i - 1);
199 *ci_metric_lnode_out = xs_gnode + local_owned_cell_idx_i; // Convert to local node index for cell origin
200
201 // Select random local owned cell index in J-direction
202 ierr = PetscRandomGetValueReal(*rand_logic_j_ptr, &r_val); CHKERRQ(ierr); // Dereference RNG pointer
203 local_owned_cell_idx_j = (PetscInt)(r_val * num_owned_cells_on_rank_j);
204 local_owned_cell_idx_j = PetscMin(PetscMax(0, local_owned_cell_idx_j), num_owned_cells_on_rank_j - 1);
205 *cj_metric_lnode_out = ys_gnode + local_owned_cell_idx_j;
206
207 // Select random local owned cell index in K-direction
208 ierr = PetscRandomGetValueReal(*rand_logic_k_ptr, &r_val); CHKERRQ(ierr); // Dereference RNG pointer
209 local_owned_cell_idx_k = (PetscInt)(r_val * num_owned_cells_on_rank_k);
210 local_owned_cell_idx_k = PetscMin(PetscMax(0, local_owned_cell_idx_k), num_owned_cells_on_rank_k - 1);
211 *ck_metric_lnode_out = zs_gnode + local_owned_cell_idx_k;
212
213 LOG_ALLOW(LOCAL, LOG_DEBUG, "Rank %d: Selected Cell (Owned Idx: %d,%d,%d -> LNodeStart: %d,%d,%d). OwnedCells(i,j,k): (%d,%d,%d). GhostNodeStarts(xs,ys,zs): (%d,%d,%d) \n",
214 rank_for_logging, local_owned_cell_idx_i, local_owned_cell_idx_j, local_owned_cell_idx_k,
215 *ci_metric_lnode_out, *cj_metric_lnode_out, *ck_metric_lnode_out,
216 num_owned_cells_on_rank_i, num_owned_cells_on_rank_j, num_owned_cells_on_rank_k,
217 xs_gnode, ys_gnode, zs_gnode);
218
219
220 // --- 2. Generate Random Intra-Cell Logical Coordinates [0,1) for MetricLogicalToPhysical ---
221 ierr = PetscRandomGetValueReal(*rand_logic_i_ptr, xi_metric_logic_out); CHKERRQ(ierr); // Re-use RNGs
222 ierr = PetscRandomGetValueReal(*rand_logic_j_ptr, eta_metric_logic_out); CHKERRQ(ierr);
223 ierr = PetscRandomGetValueReal(*rand_logic_k_ptr, zta_metric_logic_out); CHKERRQ(ierr);
224
225 // Ensure logical coordinates are strictly within [0,1) for robustness with MetricLogicalToPhysical
226 *xi_metric_logic_out = PetscMin(*xi_metric_logic_out, 1.0 - 1.0e-7);
227 *eta_metric_logic_out = PetscMin(*eta_metric_logic_out, 1.0 - 1.0e-7);
228 *zta_metric_logic_out = PetscMin(*zta_metric_logic_out, 1.0 - 1.0e-7);
229 // Ensure they are not negative either (though [0,1) RNGs shouldn't produce this)
230 *xi_metric_logic_out = PetscMax(*xi_metric_logic_out, 0.0);
231 *eta_metric_logic_out = PetscMax(*eta_metric_logic_out, 0.0);
232 *zta_metric_logic_out = PetscMax(*zta_metric_logic_out, 0.0);
233
234 } else {
235 // This rank does not own any 3D cells (e.g., in a 1D or 2D decomposition,
236 // or if the global domain itself is not 3D in terms of cells).
237 // *can_place_in_volume_out remains PETSC_FALSE.
238 LOG_ALLOW(LOCAL, LOG_WARNING, "Rank %d: Cannot place particle volumetrically. Rank has zero owned cells in at least one dimension (owned cells i,j,k: %d,%d,%d).\n",
239 rank_for_logging, num_owned_cells_on_rank_i, num_owned_cells_on_rank_j, num_owned_cells_on_rank_k);
240 }
241
243 PetscFunctionReturn(0);
244}
245
246#undef __FUNCT__
247#define __FUNCT__ "InitializeParticleBasicProperties"
248
249/**
250 * @brief Initializes basic properties for particles on the local process.
251 *
252 * This function assigns initial physical positions, Particle IDs (PIDs), and placeholder
253 * cell IDs to particles. The method of position initialization depends on
254 * `simCtx->ParticleInitialization`:
255 * - Mode 0 (Surface): Particles are placed on a designated inlet surface if the
256 * current rank services that surface. Otherwise, they are placed at (0,0,0)
257 * to be migrated to the correct rank later.
258 * - Mode 1 (Volumetric): Particles are placed randomly within a cell owned by
259 * the current rank.
260 *
261 * The logical coordinates for placement are generated using provided random number
262 * generators. These logical coordinates are then transformed to physical coordinates
263 * using `MetricLogicalToPhysical`.
264 *
265 * @param user Pointer to the UserCtx structure, containing simulation settings and grid information.
266 * @param particlesPerProcess The number of particles to initialize on this MPI rank.
267 * @param rand_logic_i Pointer to a PetscRandom generator for the xi logical coordinate.
268 * @param rand_logic_j Pointer to a PetscRandom generator for the eta logical coordinate.
269 * @param rand_logic_k Pointer to a PetscRandom generator for the zeta logical coordinate.
270 * @param bboxlist (Unused in this function for placement) Pointer to the bounding box list;
271 * provided for API consistency but not used for determining initial positions here.
272 * Particle positions are determined by logical-to-physical mapping based on rank's owned cells.
273 * @return PetscErrorCode 0 on success, non-zero on failure.
274 */
275static PetscErrorCode InitializeParticleBasicProperties(UserCtx *user,
276 PetscInt particlesPerProcess,
277 PetscRandom *rand_logic_i,
278 PetscRandom *rand_logic_j,
279 PetscRandom *rand_logic_k,
280 BoundingBox *bboxlist) // bboxlist unused for placement
281{
282 PetscErrorCode ierr;
283 DM swarm = user->swarm;
284 PetscReal *positions_field = NULL; // Pointer to swarm field for physical positions (x,y,z)
285 PetscInt64 *particleIDs = NULL; // Pointer to swarm field for Particle IDs
286 PetscInt *cellIDs_petsc = NULL; // Pointer to swarm field for DMSwarm_CellID (i,j,k of containing cell)
287 PetscInt *status_field = NULL; // Pointer to swarm field for DMSwarm_location_status(NEEDS_LOCATION etc)
288 PetscMPIInt rank,size; // MPI rank of the current process, and total number of ranks.
289 const Cmpnts ***coor_nodes_local_array; // Read-only access to local node coordinates (from user->da)
290 Vec Coor_local; // Local vector for node coordinates
291 DMDALocalInfo info; // Local grid information (node-based) from user->da
292 PetscInt xs_gnode_rank, ys_gnode_rank, zs_gnode_rank; // Local starting node indices (incl. ghosts) of rank's DA patch
293 PetscInt IM_nodes_global, JM_nodes_global, KM_nodes_global; // Global node counts in each direction
294
295 // Variables for surface initialization (Mode 0)
296 PetscBool can_this_rank_service_inlet = PETSC_FALSE;
297
298 PetscFunctionBeginUser;
299
301
302 SimCtx *simCtx = user->simCtx;
303
304 // --- 1. Input Validation and Basic Setup ---
305 if (!user || !rand_logic_i || !rand_logic_j || !rand_logic_k) {
306 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Null user or RNG pointer.");
307 }
308 ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank); CHKERRQ(ierr);
309 ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size); CHKERRQ(ierr);
310
311 // Get DMDA information for the node-centered coordinate grid (user->da)
312 ierr = DMGetCoordinatesLocal(user->da, &Coor_local); CHKERRQ(ierr);
313 if (!Coor_local) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "DMGetCoordinatesLocal for user->da returned NULL Coor_local.");
314 ierr = DMDAVecGetArrayRead(user->fda, Coor_local, (void*)&coor_nodes_local_array); CHKERRQ(ierr);
315 ierr = DMDAGetLocalInfo(user->da, &info); CHKERRQ(ierr);
316 ierr = DMDAGetCorners(user->da, &xs_gnode_rank, &ys_gnode_rank, &zs_gnode_rank, NULL, NULL, NULL); CHKERRQ(ierr);
317 ierr = DMDAGetInfo(user->da, NULL, &IM_nodes_global, &JM_nodes_global, &KM_nodes_global, NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL); CHKERRQ(ierr);
318
319 // Modification to IM_nodes_global etc. to account for 1-cell halo in each direction.
320 IM_nodes_global -= 1; JM_nodes_global -= 1; KM_nodes_global -= 1;
321
322 const PetscInt IM_cells_global = IM_nodes_global > 0 ? IM_nodes_global - 1 : 0;
323 const PetscInt JM_cells_global = JM_nodes_global > 0 ? JM_nodes_global - 1 : 0;
324 const PetscInt KM_cells_global = KM_nodes_global > 0 ? KM_nodes_global - 1 : 0;
325
326 LOG_ALLOW(LOCAL, LOG_INFO, "Rank %d: Initializing %d particles. Mode: %s.\n",
327 rank, particlesPerProcess, ParticleInitializationToString(simCtx->ParticleInitialization));
328
329 // --- 2. Pre-computation for Surface Initialization (PARTICLE_INIT_SURFACE_RANDOM and PARTICLE_INIT_SURFACE_EDGES) ---
331 simCtx->ParticleInitialization == PARTICLE_INIT_SURFACE_EDGES) { // Surface initialization
332 ierr = CanRankServiceInletFace(user, &info, IM_nodes_global, JM_nodes_global, KM_nodes_global, &can_this_rank_service_inlet); CHKERRQ(ierr);
333 if (can_this_rank_service_inlet) {
334 LOG_ALLOW(LOCAL, LOG_INFO, "Rank %d: Will attempt to place particles on inlet face %s.\n", rank, BCFaceToString((BCFace)user->identifiedInletBCFace));
335 } else {
336 LOG_ALLOW(LOCAL, LOG_INFO, "Rank %d: Cannot service inlet face %s. Particles will be at Inlet Center (%.6f,%.6f,%.6f) and rely on migration.\n", rank, BCFaceToString((BCFace)user->identifiedInletBCFace),user->simCtx->CMx_c,user->simCtx->CMy_c,user->simCtx->CMz_c);
337 }
338 }
339
340 // --- 3. Get Access to Swarm Fields ---
341 ierr = DMSwarmGetField(swarm, "position", NULL, NULL, (void**)&positions_field); CHKERRQ(ierr);
342 ierr = DMSwarmGetField(swarm, "DMSwarm_pid", NULL, NULL, (void**)&particleIDs); CHKERRQ(ierr);
343 ierr = DMSwarmGetField(swarm, "DMSwarm_CellID", NULL, NULL, (void**)&cellIDs_petsc); CHKERRQ(ierr);
344 ierr = DMSwarmGetField(swarm, "DMSwarm_location_status",NULL,NULL,(void**)&status_field); CHKERRQ(ierr);
345
346 // --- 4. Determine Starting Global PID for this Rank ---
347 PetscInt particles_per_rank_ideal = simCtx->np / size; // Assumes user->size is PETSC_COMM_WORLD size
348 PetscInt remainder_particles = simCtx->np % size;
349 PetscInt base_pid_for_rank = rank * particles_per_rank_ideal + PetscMin(rank, remainder_particles);
350 // This calculation must match how particlesPerProcess was determined (e.g., in DistributeParticles).
351
352 // --- 5. Loop Over Particles to Initialize ---
353 for (PetscInt p = 0; p < particlesPerProcess; p++) {
354 PetscInt idx = p;
355 PetscInt ci_metric_lnode, cj_metric_lnode, ck_metric_lnode;
356 PetscReal xi_metric_logic, eta_metric_logic, zta_metric_logic;
357 Cmpnts phys_coords = {0.0, 0.0, 0.0};
358 PetscBool particle_placed_by_this_rank = PETSC_FALSE;
359
360 if (simCtx->ParticleInitialization == PARTICLE_INIT_SURFACE_RANDOM) { // --- 5.a. Surface Random Initialization ---
361 if (can_this_rank_service_inlet) {
362 ierr = GetRandomCellAndLogicalCoordsOnInletFace(user, &info, xs_gnode_rank, ys_gnode_rank, zs_gnode_rank,
363 IM_nodes_global, JM_nodes_global, KM_nodes_global,
364 rand_logic_i, rand_logic_j, rand_logic_k,
365 &ci_metric_lnode, &cj_metric_lnode, &ck_metric_lnode,
366 &xi_metric_logic, &eta_metric_logic, &zta_metric_logic); CHKERRQ(ierr);
367 ierr = MetricLogicalToPhysical(user, coor_nodes_local_array,
368 ci_metric_lnode, cj_metric_lnode, ck_metric_lnode,
369 xi_metric_logic, eta_metric_logic, zta_metric_logic,
370 &phys_coords); CHKERRQ(ierr);
371 particle_placed_by_this_rank = PETSC_TRUE;
372 }else{
373 // Rank cannot service inlet - place at inlet center to be migrated later
374 phys_coords.x = user->simCtx->CMx_c;
375 phys_coords.y = user->simCtx->CMy_c;
376 phys_coords.z = user->simCtx->CMz_c;
377 particle_placed_by_this_rank = PETSC_FALSE; // Relies on migration
378 }
379 }else if(simCtx->ParticleInitialization == PARTICLE_INIT_SURFACE_EDGES) { // --- 5.a1. Surface Edges Initialization (deterministic) ---
380 if(can_this_rank_service_inlet) {
381 PetscInt64 particle_global_id = (PetscInt64)(base_pid_for_rank + p);
382 ierr = GetDeterministicFaceGridLocation(user,&info,xs_gnode_rank, ys_gnode_rank, zs_gnode_rank,
383 IM_cells_global, JM_cells_global, KM_cells_global,
384 particle_global_id,
385 &ci_metric_lnode, &cj_metric_lnode, &ck_metric_lnode,
386 &xi_metric_logic, &eta_metric_logic, &zta_metric_logic,
387 &particle_placed_by_this_rank); CHKERRQ(ierr);
388 if(particle_placed_by_this_rank){
389 ierr = MetricLogicalToPhysical(user, coor_nodes_local_array,
390 ci_metric_lnode, cj_metric_lnode, ck_metric_lnode,
391 xi_metric_logic, eta_metric_logic, zta_metric_logic,
392 &phys_coords); CHKERRQ(ierr);
393 }else{
394 // Even if rank can service face, it may not own the portion of the face that a particle is placed in.
395 phys_coords.x = user->simCtx->CMx_c;
396 phys_coords.y = user->simCtx->CMy_c;
397 phys_coords.z = user->simCtx->CMz_c;
398 }
399 }else{
400 // Rank cannot service inlet - place at inlet center to be migrated later
401 phys_coords.x = user->simCtx->CMx_c;
402 phys_coords.y = user->simCtx->CMy_c;
403 phys_coords.z = user->simCtx->CMz_c;
404 particle_placed_by_this_rank = PETSC_FALSE; // Relies on migration
405 }
406 }else if(simCtx->ParticleInitialization == PARTICLE_INIT_VOLUME){ // --- 5.b. Volumetric Initialization ---
407 PetscBool can_place_volumetrically;
408 ierr = DetermineVolumetricInitializationParameters(user, &info, xs_gnode_rank, ys_gnode_rank, zs_gnode_rank,
409 rand_logic_i, rand_logic_j, rand_logic_k,
410 &ci_metric_lnode, &cj_metric_lnode, &ck_metric_lnode,
411 &xi_metric_logic, &eta_metric_logic, &zta_metric_logic,
412 &can_place_volumetrically); CHKERRQ(ierr);
413 if(can_place_volumetrically){
414 ierr = MetricLogicalToPhysical(user, coor_nodes_local_array,
415 ci_metric_lnode, cj_metric_lnode, ck_metric_lnode,
416 xi_metric_logic, eta_metric_logic, zta_metric_logic,
417 &phys_coords); CHKERRQ(ierr);
418 particle_placed_by_this_rank = PETSC_TRUE;
419 } else {
421 "Rank %d: PID %lld (idx %ld) (Volumetric Mode) - DetermineVolumetric... returned false. Default Phys: (%.2f,%.2f,%.2f).\n",
422 rank, (long long)(base_pid_for_rank + p), (long)p, phys_coords.x, phys_coords.y, phys_coords.z);
423 }
424 }else if(simCtx->ParticleInitialization == PARTICLE_INIT_POINT_SOURCE){ // --- 5.c. Point Source Initialization ---
425 // All particles placed at the user-specified fixed point (psrc_x, psrc_y, psrc_z).
426 // No random number generation or logical-to-physical conversion needed.
427 phys_coords.x = simCtx->psrc_x;
428 phys_coords.y = simCtx->psrc_y;
429 phys_coords.z = simCtx->psrc_z;
430 particle_placed_by_this_rank = PETSC_TRUE;
431 }else {
432 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unknown ParticleInitialization mode %d.", simCtx->ParticleInitialization);
433 }
434
435 // --- 5.c. Store Particle Properties ---
436 positions_field[3*p+0] = phys_coords.x;
437 positions_field[3*p+1] = phys_coords.y;
438 positions_field[3*p+2] = phys_coords.z;
439
440 particleIDs[p] = (PetscInt64)base_pid_for_rank + p;
441 cellIDs_petsc[3*p+0] = -1; cellIDs_petsc[3*p+1] = -1; cellIDs_petsc[3*p+2] = -1;
442 status_field[p] = UNINITIALIZED;
443
444 // --- 5.d. Logging for this particle ---
445 if (particle_placed_by_this_rank) {
446 LOG_LOOP_ALLOW(LOCAL, LOG_VERBOSE, idx, user->simCtx->LoggingFrequency,//(particlesPerProcess > 20 ? particlesPerProcess/10 : 1),
447 "Rank %d: PID %lld (idx %ld) PLACED. Mode %s. Embedded Cell:(%d,%d,%d). Logical Coords: (%.2e,%.2f,%.2f).\n Final Coords: (%.6f,%.6f,%.6f).\n",
448 rank, (long long)particleIDs[p], (long)p, ParticleInitializationToString(simCtx->ParticleInitialization),
449 ci_metric_lnode, cj_metric_lnode, ck_metric_lnode,
450 xi_metric_logic, eta_metric_logic, zta_metric_logic,
451 phys_coords.x, phys_coords.y, phys_coords.z);
452
453 } else {
454 LOG_LOOP_ALLOW(LOCAL, LOG_WARNING, idx, user->simCtx->LoggingFrequency, //(particlesPerProcess > 20 ? particlesPerProcess/10 : 1),
455 "Rank %d: PID %lld (idx %ld) Mode %s NOT placed by this rank's logic. Default Coor: (%.2f,%.2f,%.2f). Relies on migration.\n",
456 rank, (long long)particleIDs[p], (long)p, ParticleInitializationToString(simCtx->ParticleInitialization),
457 phys_coords.x, phys_coords.y, phys_coords.z);
458 }
459 }
460
461 // --- 6. Restore Pointers and Cleanup ---
462 ierr = DMSwarmRestoreField(swarm, "position", NULL, NULL, (void**)&positions_field); CHKERRQ(ierr);
463 ierr = DMSwarmRestoreField(swarm, "DMSwarm_pid", NULL, NULL, (void**)&particleIDs); CHKERRQ(ierr);
464 ierr = DMSwarmRestoreField(swarm, "DMSwarm_CellID", NULL, NULL, (void**)&cellIDs_petsc); CHKERRQ(ierr);
465 ierr = DMSwarmRestoreField(swarm, "DMSwarm_location_status", NULL, NULL, (void**)&status_field); CHKERRQ(ierr);
466 ierr = DMDAVecRestoreArrayRead(user->fda, Coor_local, (void*)&coor_nodes_local_array); CHKERRQ(ierr);
467
468 LOG_ALLOW(LOCAL, LOG_INFO, "Rank %d: Completed processing for %d particles.\n",
469 rank, particlesPerProcess);
470
472
473 PetscFunctionReturn(0);
474}
475
476#undef __FUNCT__
477#define __FUNCT__ "InitializeSwarmFieldValue"
478/**
479 * @brief Helper function to Initialize a given particle’s field value.
480 *
481 * This function performs conditional, point-level Initialization for a swarm field based on its name.
482 * For example, you might want to initialize the "velocity" field to 0.0, but the "temperature"
483 * field to a nonzero default (e.g., 300.0). This function can be extended for other fields.
484 *
485 * @param[in] fieldName Name of the swarm field.
486 * @param[in] p Particle index.
487 * @param[in] fieldDim Dimension of the field.
488 * @param[out] fieldData Pointer to the field’s data array.
489 *
490 * @return PetscErrorCode Returns 0 on success.
491 */
492static PetscErrorCode InitializeSwarmFieldValue(const char *fieldName, PetscInt p, PetscInt fieldDim, PetscReal *fieldData)
493{
494 PetscFunctionBeginUser;
495
497
498 if (strcmp(fieldName, "velocity") == 0) {
499 // For velocity, initialize all components to zero
500 for (PetscInt d = 0; d < fieldDim; d++) {
501 fieldData[fieldDim * p + d] = 0.0;
502 }
503 } else if (strcmp(fieldName, "Diffusivity") == 0) {
504 // For Diffusivity,initialize to a default value (e.g., 1.0)
505 for (PetscInt d = 0; d < fieldDim; d++) {
506 fieldData[fieldDim * p + d] = 1.0;
507 }
508 } else if (strcmp(fieldName, "DiffusivityGradient") == 0) {
509 // For DiffusivityGradient,initialize to a default value (e.g., 1.0)
510 for (PetscInt d = 0; d < fieldDim; d++) {
511 fieldData[fieldDim * p + d] = 1.0;
512 }
513 } else if (strcmp(fieldName, "Psi") == 0) {
514 for (PetscInt d = 0; d < fieldDim; d++) {
515 fieldData[fieldDim * p + d] = 0.0;
516 }
517 } else {
518 // Default: initialize all components to zero
519 for (PetscInt d = 0; d < fieldDim; d++) {
520 fieldData[fieldDim * p + d] = 0.0;
521 }
522 }
523
525 PetscFunctionReturn(0);
526}
527
528
529#undef __FUNCT__
530#define __FUNCT__ "AssignInitialFieldToSwarm"
531/**
532 * @brief Initializes a generic swarm field with point-level updates.
533 *
534 * This field-agnostic function retrieves the specified swarm field (which may be
535 * scalar or multi-component) and initializes each particle's entry using a helper
536 * that performs conditional updates based on the field name.
537 *
538 * @param[in,out] user Pointer to the UserCtx structure containing the swarm.
539 * @param[in] fieldName Name of the swarm field to initialize.
540 * @param[in] fieldDim Dimension of the field (e.g., 1 for scalar, 3 for vector).
541 *
542 * @return PetscErrorCode Returns 0 on success, non-zero on failure.
543 */
544static PetscErrorCode AssignInitialFieldToSwarm(UserCtx *user, const char *fieldName, PetscInt fieldDim)
545{
546 PetscErrorCode ierr;
547 DM swarm = user->swarm;
548 PetscReal *fieldData = NULL;
549 PetscInt nLocal;
550
551 PetscFunctionBeginUser;
552
554
555 // Get the number of local particles
556 ierr = DMSwarmGetLocalSize(swarm, &nLocal); CHKERRQ(ierr);
557 LOG_ALLOW(LOCAL,LOG_INFO, "%d local particles found.\n", nLocal);
558
559 // Retrieve the swarm field pointer for the specified fieldName
560 ierr = DMSwarmGetField(swarm, fieldName, NULL, NULL, (void**)&fieldData); CHKERRQ(ierr);
561 LOG_ALLOW(LOCAL,LOG_DEBUG, "Retrieved field '%s'.\n", fieldName);
562
563 // Loop over all particles and update the field using the helper function
564 for (PetscInt p = 0; p < nLocal; p++) {
565 ierr = InitializeSwarmFieldValue(fieldName, p, fieldDim, fieldData); CHKERRQ(ierr);
566 PetscReal disp_data[fieldDim];
567
568 for (PetscInt d = 0; d < fieldDim; d++) {
569 disp_data[d] = fieldData[fieldDim* p + d];
570 }
571 LOG_LOOP_ALLOW(LOCAL,LOG_VERBOSE,p, 100," Particle %d: %s[%d] = [%.6f, ...,%.6f].\n", p,fieldName,fieldDim,disp_data[0],disp_data[fieldDim-1]);
572 }
573
574 // Restore the swarm field pointer
575 ierr = DMSwarmRestoreField(swarm, fieldName, NULL, NULL, (void**)&fieldData); CHKERRQ(ierr);
576 LOG_ALLOW(LOCAL,LOG_INFO, "Initialization of field '%s' complete.\n", fieldName);
577
578
580
581 PetscFunctionReturn(0);
582}
583
584#undef __FUNCT__
585#define __FUNCT__ "AssignInitialPropertiesToSwarm"
586
587/**
588 * @brief Initializes all particle properties in the swarm.
589 *
590 * This function orchestrates the initialization of particle properties.
591 * It first determines the inlet face if surface initialization (Mode 0) is selected
592 * by parsing "bcs.dat".
593 * Then, it initializes basic particle properties (physical position, Particle ID,
594 * and placeholder Cell IDs) by calling `InitializeParticleBasicProperties`. This call
595 * uses the provided `rand_logic_i/j/k` RNGs, which must be pre-initialized for [0,1).
596 * The `rand_phys_x/y/z` RNGs (physically bounded) are passed but may not be used by
597 * `InitializeParticleBasicProperties` for position setting if all initialization paths
598 * use logical-to-physical mapping.
599 * Finally, it calls helper functions to initialize other registered swarm fields
600 * like "velocity", "weight", and "Psi" (scalar) to default values.
601 *
602 * @param[in,out] user Pointer to the `UserCtx` structure.
603 * @param[in] particlesPerProcess Number of particles assigned to this MPI process.
604 * @param[in] rand_phys_x RNG for physical x-coordinates (from `InitializeRandomGenerators`).
605 * @param[in] rand_phys_y RNG for physical y-coordinates (from `InitializeRandomGenerators`).
606 * @param[in] rand_phys_z RNG for physical z-coordinates (from `InitializeRandomGenerators`).
607 * @param[in] rand_logic_i RNG for i-logical dimension tasks [0,1) (from `InitializeLogicalSpaceRNGs`).
608 * @param[in] rand_logic_j RNG for j-logical dimension tasks [0,1) (from `InitializeLogicalSpaceRNGs`).
609 * @param[in] rand_logic_k RNG for k-logical dimension tasks [0,1) (from `InitializeLogicalSpaceRNGs`).
610 * @param[in] bboxlist Array of BoundingBox structures (potentially unused by IPBP).
611 *
612 * @return PetscErrorCode Returns 0 on success, non-zero on failure.
613 */
615 PetscInt particlesPerProcess,
616 PetscRandom *rand_phys_x, // RNG from original InitializeRandomGenerators
617 PetscRandom *rand_phys_y, // RNG from original InitializeRandomGenerators
618 PetscRandom *rand_phys_z, // RNG from original InitializeRandomGenerators
619 PetscRandom *rand_logic_i, // RNG from InitializeLogicalSpaceRNGs
620 PetscRandom *rand_logic_j, // RNG from InitializeLogicalSpaceRNGs
621 PetscRandom *rand_logic_k, // RNG from InitializeLogicalSpaceRNGs
622 BoundingBox *bboxlist)
623{
624 PetscErrorCode ierr;
625 PetscFunctionBeginUser;
626
628
629 SimCtx *simCtx = user->simCtx;
630
631 // --- 0. Input Validation ---
632 if (!user || !bboxlist || !rand_logic_i || !rand_logic_j || !rand_logic_k || !rand_phys_x || !rand_phys_y || !rand_phys_z) {
633 // Check all RNGs now as they are passed in
634 LOG_ALLOW(GLOBAL, LOG_ERROR, "Null user, bboxlist, or RNG pointer.\n");
635 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Null input detected.");
636 }
637
638 LOG_ALLOW(GLOBAL, LOG_INFO, "Initializing swarm with %d particles per process. Mode: %s.\n",
639 particlesPerProcess, ParticleInitializationToString(simCtx->ParticleInitialization));
640
641 // --- 1. Parse BCS File for Inlet Information (if surface initialization) ---
643 simCtx->ParticleInitialization == PARTICLE_INIT_SURFACE_EDGES) { // Surface initialization
644 if(user->inletFaceDefined == PETSC_FALSE){
645 LOG_ALLOW(GLOBAL, LOG_ERROR, "Particle Initialization on inlet surface selected, but no INLET face was identified from bcs.dat. Cannot proceed.\n");
646 SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONGSTATE, "ParticleInitialization Mode 0 requires an INLET face to be defined in bcs.dat.");
647 }else{
648 LOG_ALLOW(GLOBAL, LOG_INFO, "After Parsing BCS file for Inlet, Inlet face = %s\n", BCFaceToString((BCFace)user->identifiedInletBCFace));
649 }
650 }
651
652 // --- 2. Initialize Basic Particle Properties (Position, PID, Cell IDs placeholder) ---
653 // The rand_logic_i/j/k are now passed directly.
654 // The rand_phys_x/y/z are passed but InitializeParticleBasicProperties (refactored version)
655 // will not use them for setting positions if all its paths use logical-to-physical mapping.
656 LOG_ALLOW(GLOBAL, LOG_DEBUG, "Calling InitializeParticleBasicProperties.\n");
657 ierr = InitializeParticleBasicProperties(user, particlesPerProcess,
658 rand_logic_i, rand_logic_j, rand_logic_k,
659 bboxlist); // bboxlist passed along
660 CHKERRQ(ierr);
661 LOG_ALLOW(GLOBAL, LOG_INFO, "Successfully initialized basic particle properties.\n");
662
663 // Note: The logical RNGs (rand_logic_i/j/k) are NOT destroyed here.
664 // They were created externally (e.g., by InitializeLogicalSpaceRNGs) and
665 // should be destroyed externally (e.g., in FinalizeSwarmSetup).
666 // Same for rand_phys_x/y/z.
667
668 // --- 3. Initialize Other Swarm Fields (Velocity, Weight, Pressure, etc.) ---
669 LOG_ALLOW(GLOBAL, LOG_DEBUG, "Initializing 'velocity' field.\n");
670 ierr = AssignInitialFieldToSwarm(user, "velocity", 3); CHKERRQ(ierr);
671 LOG_ALLOW(LOCAL, LOG_INFO, "'velocity' field initialization complete.\n");
672
673 LOG_ALLOW(GLOBAL, LOG_DEBUG, "Initializing 'weight' field.\n");
674 ierr = AssignInitialFieldToSwarm(user, "weight", 3); CHKERRQ(ierr); // Assuming weight is vec3
675 LOG_ALLOW(LOCAL, LOG_INFO, "'weight' field initialization complete.\n");
676
677 LOG_ALLOW(GLOBAL, LOG_DEBUG, "Initializing 'Diffusivity' field.\n");
678 ierr = AssignInitialFieldToSwarm(user, "Diffusivity", 1); CHKERRQ(ierr);
679 LOG_ALLOW(GLOBAL, LOG_INFO, "'Diffusivity' field initialization complete.\n");
680
681 LOG_ALLOW(GLOBAL, LOG_DEBUG, "Initializing 'DiffusivityGradient' field.\n");
682 ierr = AssignInitialFieldToSwarm(user, "DiffusivityGradient", 3); CHKERRQ(ierr);
683 LOG_ALLOW(GLOBAL, LOG_INFO, "'DiffusivityGradient' field initialization complete.\n");
684
685 LOG_ALLOW(GLOBAL, LOG_DEBUG, "Initializing 'Psi' (Scalar) field.\n");
686 ierr = AssignInitialFieldToSwarm(user, "Psi", 1); CHKERRQ(ierr);
687 LOG_ALLOW(GLOBAL, LOG_INFO, "'P' field initialization complete.\n");
688
689 LOG_ALLOW(GLOBAL, LOG_INFO, "Successfully completed all swarm property initialization.\n");
690
691
693
694 PetscFunctionReturn(0);
695}
696
697
698#undef __FUNCT__
699#define __FUNCT__ "DistributeParticles"
700/**
701 * @brief Distributes particles evenly across MPI processes, handling any remainders.
702 *
703 * This function calculates the number of particles each MPI process should handle,
704 * distributing the remainder particles to the first few ranks if necessary.
705 *
706 * @param[in] numParticles Total number of particles to create across all MPI processes.
707 * @param[in] rank MPI rank of the current process.
708 * @param[in] size Total number of MPI processes.
709 * @param[out] particlesPerProcess Number of particles assigned to the current MPI process.
710 * @param[out] remainder Remainder particles when dividing numParticles by size.
711 *
712 * @return PetscErrorCode Returns 0 on success, non-zero on failure.
713 */
714PetscErrorCode DistributeParticles(PetscInt numParticles, PetscMPIInt rank, PetscMPIInt size, PetscInt* particlesPerProcess, PetscInt* remainder) {
715
716 PetscFunctionBeginUser;
717
719 // Calculate the base number of particles per process
720 *particlesPerProcess = numParticles / size;
721 *remainder = numParticles % size;
722
723 // Distribute the remainder particles to the first 'remainder' ranks
724 if (rank < *remainder) {
725 *particlesPerProcess += 1;
726 LOG_ALLOW_SYNC(GLOBAL,LOG_INFO,"Rank %d receives an extra particle. Total: %d\n", rank, *particlesPerProcess);
727 } else {
728 LOG_ALLOW_SYNC(GLOBAL,LOG_INFO, "Rank %d receives %d particles.\n", rank, *particlesPerProcess);
729 }
730
732 PetscFunctionReturn(0);
733}
734
735
736#undef __FUNCT__
737#define __FUNCT__ "FinalizeSwarmSetup"
738/**
739 * @brief Finalizes the swarm setup by destroying random generators and logging completion.
740 *
741 * This function cleans up resources by destroying random number generators and LOG_ALLOWs the completion of swarm setup.
742 *
743 * @param[in] randx Random number generator for the x-coordinate.
744 * @param[in] randy Random number generator for the y-coordinate.
745 * @param[in] randz Random number generator for the z-coordinate.
746 * @param[in] rand_logic_i Random number generator for the xi-coordinate.
747 * @param[in] rand_logic_j Random number generator for the eta-coordinate.
748 * @param[in] rand_logic_k Random number generator for the zeta-coordinate.
749 *
750 * @return PetscErrorCode Returns 0 on success, non-zero on failure.
751 */
752PetscErrorCode FinalizeSwarmSetup(PetscRandom *randx, PetscRandom *randy, PetscRandom *randz, PetscRandom *rand_logic_i, PetscRandom *rand_logic_j, PetscRandom *rand_logic_k) {
753 PetscErrorCode ierr; // Error code for PETSc functions
754 PetscFunctionBeginUser;
756 // Destroy random number generators to free resources
757 // Physical space
758 ierr = PetscRandomDestroy(randx); CHKERRQ(ierr);
759 ierr = PetscRandomDestroy(randy); CHKERRQ(ierr);
760 ierr = PetscRandomDestroy(randz); CHKERRQ(ierr);
761 // Logical space
762 ierr = PetscRandomDestroy(rand_logic_i); CHKERRQ(ierr);
763 ierr = PetscRandomDestroy(rand_logic_j); CHKERRQ(ierr);
764 ierr = PetscRandomDestroy(rand_logic_k); CHKERRQ(ierr);
765
766 LOG_ALLOW(LOCAL,LOG_DEBUG,"Destroyed all random number generators.\n");
767
769 PetscFunctionReturn(0);
770}
771
772#undef __FUNCT__
773#define __FUNCT__ "CreateParticleSwarm"
774/**
775 * @brief Creates and initializes a Particle Swarm.
776 *
777 * This function sets up a DMSwarm within the provided UserCtx structure, initializes
778 * particle fields, and distributes particles across MPI processes. It ensures that
779 * the number of particles is evenly divided among the available MPI ranks. If the total
780 * number of particles isn't divisible by the number of processes, the remainder is distributed
781 * to the first few ranks.
782 *.
783*
784* @param[in,out] user Pointer to the UserCtx structure containing the simulation context.
785* @param[in] numParticles Total number of particles to create across all MPI processes.
786* @param[in] bboxlist Pointer to an array of BoundingBox structures, one per rank.
787*
788* @param[in] particlesPerProcess
789* @return PetscErrorCode Returns 0 on success, non-zero on failure.
790*
791* @note
792* - Ensure that `numParticles` is a positive integer.
793* - The `control.dat` file should contain necessary PETSc options.
794* - The `bboxlist` array should be properly populated before calling this function.
795*/
796PetscErrorCode CreateParticleSwarm(UserCtx *user, PetscInt numParticles, PetscInt *particlesPerProcess, BoundingBox *bboxlist) {
797 PetscErrorCode ierr; // PETSc error handling variable
798 PetscMPIInt rank, size; // Variables to store MPI rank and size
799 PetscInt remainder = 0; // Remainder of particles after division
800
801 PetscFunctionBeginUser;
803 // Validate input parameters
804 if (numParticles <= 0) {
805 LOG_ALLOW(GLOBAL,LOG_DEBUG, "Number of particles must be positive. Given: %d\n", numParticles);
806 return PETSC_ERR_ARG_OUTOFRANGE;
807 }
808
809 // Retrieve MPI rank and size
810 ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank); CHKERRQ(ierr);
811 ierr = MPI_Comm_size(PETSC_COMM_WORLD, &size); CHKERRQ(ierr);
812 LOG_ALLOW(GLOBAL,LOG_INFO," Domain dimensions: [%.2f,%.2f],[%.2f,%.2f],[%.2f,%.2f] \n",
813 user->Min_X,user->Max_X,user->Min_Y,user->Max_Y, user->Min_Z,user->Max_Z);
814 LOG_ALLOW_SYNC(GLOBAL,LOG_DEBUG, "[Rank %d] Local Bounding Box: [%.2f,%.2f],[%.2f,%.2f],[%.2f,%.2f] \n",
815 rank,user->bbox.min_coords.x,user->bbox.max_coords.x,
816 user->bbox.min_coords.y,user->bbox.max_coords.y,
817 user->bbox.min_coords.z,user->bbox.max_coords.z);
818 // Distribute particles among MPI processes
819 ierr = DistributeParticles(numParticles, rank, size, particlesPerProcess, &remainder); CHKERRQ(ierr);
820
821 // Initialize the DMSwarm - creates the swarm, sets the type and dimension
822 ierr = InitializeSwarm(user); CHKERRQ(ierr);
823
824 if (user->da) {
825 ierr = DMSwarmSetCellDM(user->swarm, user->da); CHKERRQ(ierr);
826 LOG_ALLOW(LOCAL,LOG_INFO,"Associated DMSwarm with Cell DM (user->da).\n");
827 } else {
828 // If user->da is essential for your simulation logic with particles, this should be a fatal error.
829 LOG_ALLOW(GLOBAL, LOG_WARNING, "user->da (Cell DM for Swarm) is NULL. Cell-based swarm operations might fail.\n");
830 // SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONGSTATE, "user->da (Cell DM) is NULL but required.");
831 }
832
833 // Register particle fields (position, velocity, CellID, weight, etc.)
834 ierr = RegisterParticleFields(user->swarm); CHKERRQ(ierr);
835
836 // Set the local number of particles for this rank and additional buffer for particle migration
837 ierr = DMSwarmSetLocalSizes(user->swarm, *particlesPerProcess, numParticles); CHKERRQ(ierr);
838 LOG_ALLOW(GLOBAL,LOG_INFO, "Set local swarm size: %d particles.\n", *particlesPerProcess);
839
840 // Optionally, LOG_ALLOW detailed DM info in debug mode
841 if (get_log_level() == LOG_DEBUG && is_function_allowed(__func__)) {
842 LOG_ALLOW(GLOBAL,LOG_DEBUG,"Viewing DMSwarm:\n");
843 ierr = DMView(user->swarm, PETSC_VIEWER_STDOUT_WORLD); CHKERRQ(ierr);
844 }
845
846 LOG_ALLOW(GLOBAL,LOG_INFO, "Particle swarm creation and initialization complete.\n");
847
849 PetscFunctionReturn(0);
850}
851
852// NOTE: The following two functions are helpers for unpacking and updating particle data
853// between DMSwarm fields and the Particle struct used in simulation logic.
854// While Swarm fields store data in arrays, the Particle struct provides a convenient
855// way to manipulate individual particle properties during simulation steps.
856
857#undef __FUNCT__
858#define __FUNCT__ "UnpackSwarmFields"
859
860/**
861 * @brief Initializes a Particle struct with data from DMSwarm fields.
862 *
863 * This helper function populates a Particle structure using data retrieved from DMSwarm fields.
864 *
865 * @param[in] i Index of the particle in the DMSwarm.
866 * @param[in] PIDs Pointer to the array of particle IDs.
867 * @param[in] weights Pointer to the array of particle weights.
868 * @param[in] positions Pointer to the array of particle positions.
869 * @param[in] cellIndices Pointer to the array of particle cell indices.
870 * @param[in] velocities Pointer to the array of particle velocities.
871 * @param[in] LocStatus Pointer to the array of cell location status indicators.
872 * @param[in] diffusivity Pointer to the array of particle diffusivities.
873 * @param[in] diffusivitygradient Pointer to the array of particle diffusivity gradients.
874 * @param[in] psi Pointer to the array of particle psi values.
875 * @param[out] particle Pointer to the Particle struct to initialize.
876 *
877 * @return PetscErrorCode Returns `0` on success, non-zero on failure.
878 */
879PetscErrorCode UnpackSwarmFields(PetscInt i, const PetscInt64 *PIDs, const PetscReal *weights,
880 const PetscReal *positions, const PetscInt *cellIndices,
881 PetscReal *velocities,PetscInt *LocStatus,PetscReal *diffusivity, Cmpnts *diffusivitygradient, PetscReal *psi, Particle *particle) {
882 PetscFunctionBeginUser;
883
885
886 PetscMPIInt rank;
887 PetscErrorCode ierr;
888
889 ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank); CHKERRQ(ierr);
890
891 if (particle == NULL) {
892 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Output Particle pointer is NULL. \n");
893 }
894
895 // logging the start of particle initialization
896 LOG_ALLOW(LOCAL,LOG_DEBUG, "[Rank %d]Unpacking Particle [%d] with PID: %ld.\n",rank, i, PIDs[i]);
897
898 // Initialize PID
899 if(PIDs == NULL){
900 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Input PIDs pointer is NULL.\n");
901 }
902 particle->PID = PIDs[i];
903 LOG_ALLOW(LOCAL,LOG_VERBOSE, "[Rank %d]Particle [%d] PID set to: %ld.\n", rank,i, particle->PID);
904
905 // Initialize weights
906 if(weights == NULL){
907 particle->weights.x = 1.0;
908 particle->weights.y = 1.0;
909 particle->weights.z = 1.0;
910 LOG_ALLOW(LOCAL,LOG_WARNING, "[Rank %d]Particle [%d] weights pointer is NULL. Defaulting weights to (1.0, 1.0, 1.0).\n", rank,i);
911 }else{
912 particle->weights.x = weights[3 * i];
913 particle->weights.y = weights[3 * i + 1];
914 particle->weights.z = weights[3 * i + 2];
915 LOG_ALLOW(LOCAL,LOG_VERBOSE, "[Rank %d]Particle [%d] weights set to: (%.6f, %.6f, %.6f).\n",
916 rank,i, particle->weights.x, particle->weights.y, particle->weights.z);
917 }
918 // Initialize locations
919 if(positions == NULL){
920 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Input positions pointer is NULL.\n");
921 }
922 particle->loc.x = positions[3 * i];
923 particle->loc.y = positions[3 * i + 1];
924 particle->loc.z = positions[3 * i + 2];
925 LOG_ALLOW(LOCAL,LOG_VERBOSE, "[Rank %d]Particle [%d] location set to: (%.6f, %.6f, %.6f).\n",
926 rank,i, particle->loc.x, particle->loc.y, particle->loc.z);
927
928 // Initialize velocities (assuming default zero; modify if necessary)
929 if(velocities == NULL){
930 particle->vel.x = 0.0;
931 particle->vel.y = 0.0;
932 particle->vel.z = 0.0;
933 LOG_ALLOW(LOCAL,LOG_WARNING, "[Rank %d]Particle [%d] velocities pointer is NULL. Defaulting velocities to (0.0, 0.0, 0.0).\n", rank,i);
934 }else{
935 particle->vel.x = velocities[3 * i];
936 particle->vel.y = velocities[3 * i + 1];
937 particle->vel.z = velocities[3 * i + 2];
938 LOG_ALLOW(LOCAL,LOG_VERBOSE,"[Rank %d]Particle [%d] velocities unpacked to: [%.6f,%.6f,%.6f].\n",rank, i,particle->vel.x,particle->vel.y,particle->vel.z);
939 }
940
941 // Initialize diffusivity
942 if(diffusivity == NULL){
943 particle->diffusivity = 1.0; // Default diffusivity
944 }else{
945 particle->diffusivity = diffusivity[i];
946 }
947 LOG_ALLOW(LOCAL,LOG_VERBOSE,"[Rank %d]Particle [%d] diffusivity set to: %.6f.\n",rank,i, particle->diffusivity);
948
949 // Initialize diffusivity gradient
950 if(diffusivitygradient == NULL){
951 particle->diffusivitygradient.x = 0.0;
952 particle->diffusivitygradient.y = 0.0;
953 particle->diffusivitygradient.z = 0.0;
954 LOG_ALLOW(LOCAL,LOG_WARNING, "[Rank %d]Particle [%d] diffusivity gradient pointer is NULL. Defaulting to (0.0, 0.0, 0.0).\n", rank,i);
955 }else{
956 particle->diffusivitygradient.x = diffusivitygradient[i].x;
957 particle->diffusivitygradient.y = diffusivitygradient[i].y;
958 particle->diffusivitygradient.z = diffusivitygradient[i].z;
959 }
960 LOG_ALLOW(LOCAL,LOG_VERBOSE,"[Rank %d]Particle [%d] diffusivity gradient set to: (%.6f, %.6f, %.6f).\n", rank,i,particle->diffusivitygradient.x,particle->diffusivitygradient.y,particle->diffusivitygradient.z);
961
962 // Initialize psi
963 if(psi == NULL){
964 particle->psi = 0.0; // Default psi
965 }else{
966 particle->psi = psi[i];
967 }
968 LOG_ALLOW(LOCAL,LOG_VERBOSE,"[Rank %d]Particle [%d] psi set to: %.6f.\n",rank,i, particle->psi);
969
970 // Initialize cell indices
971 if(cellIndices == NULL){
972 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Input cellIndices pointer is NULL.\n");
973 }
974 particle->cell[0] = cellIndices[3 * i];
975 particle->cell[1] = cellIndices[3 * i + 1];
976 particle->cell[2] = cellIndices[3 * i + 2];
977 LOG_ALLOW(LOCAL,LOG_VERBOSE,"[Rank %d]Particle [%d] cell indices set to: [%d, %d, %d].\n",rank,i, particle->cell[0], particle->cell[1], particle->cell[2]);
978
979 if(LocStatus == NULL){
980 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Input LocStatus pointer is NULL.\n");
981 }
982 // Initialize location status
983 particle->location_status = (ParticleLocationStatus)LocStatus[i];
984 LOG_ALLOW(LOCAL,LOG_VERBOSE, "[Rank %d]Particle [%d] Status set to: %d.\n",rank, i, particle->location_status);
985
986 // The destination_rank is only set by the location search, not read from the swarm,
987 // so we initialize it to a known invalid state.
988 particle->destination_rank = MPI_PROC_NULL;
989
990 // logging the completion of particle initialization
991 LOG_ALLOW(LOCAL,LOG_DEBUG,"[Rank %d]Completed initialization of Particle [%d]. \n", rank,i);
992
993
995
996 PetscFunctionReturn(0);
997}
998
999#undef __FUNCT__
1000#define __FUNCT__ "UpdateSwarmFields"
1001/**
1002 * @brief Updates DMSwarm data arrays from a Particle struct.
1003 *
1004 * This function writes data from the `Particle` struct back into the raw DMSwarm arrays.
1005 * It is robust: if any array pointer is NULL, that specific field is skipped.
1006 * This allows selective updating (e.g., update position but not velocity).
1007 *
1008 * @param[in] i Index of the particle in the local swarm arrays.
1009 * @param[in] particle Pointer to the Particle struct containing updated data.
1010 * @param[in,out] positions (Optional) Array of particle positions (size 3*n).
1011 * @param[in,out] velocities (Optional) Array of particle velocities (size 3*n).
1012 * @param[in,out] weights (Optional) Array of particle weights (size 3*n).
1013 * @param[in,out] cellIndices (Optional) Array of particle cell indices (size 3*n).
1014 * @param[in,out] status (Optional) Array of location status (size 1*n).
1015 * @param[in,out] diffusivity (Optional) Array of diffusivity values (size 1*n).
1016 * @param[in,out] diffusivitygradient (Optional) Array of diffusivity gradient values (size 3*n).
1017 * @param[in,out] psi (Optional) Array of scalar Psi values (size 1*n).
1018 *
1019 * @return PetscErrorCode Returns 0 on success.
1020 */
1021PetscErrorCode UpdateSwarmFields(PetscInt i, const Particle *particle,
1022 PetscReal *positions,
1023 PetscReal *velocities,
1024 PetscReal *weights,
1025 PetscInt *cellIndices,
1026 PetscInt *status,
1027 PetscReal *diffusivity,
1028 Cmpnts *diffusivitygradient,
1029 PetscReal *psi)
1030{
1031 PetscFunctionBeginUser;
1033
1034 if (!particle) {
1035 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Input Particle pointer is NULL.\n");
1036 }
1037
1038 // --- 1. Position (x, y, z) ---
1039 if (positions) {
1040 positions[3 * i + 0] = particle->loc.x;
1041 positions[3 * i + 1] = particle->loc.y;
1042 positions[3 * i + 2] = particle->loc.z;
1043 }
1044
1045 // --- 2. Velocity (u, v, w) ---
1046 if (velocities) {
1047 velocities[3 * i + 0] = particle->vel.x;
1048 velocities[3 * i + 1] = particle->vel.y;
1049 velocities[3 * i + 2] = particle->vel.z;
1050 }
1051
1052 // --- 3. Weights (i, j, k) ---
1053 if (weights) {
1054 weights[3 * i + 0] = particle->weights.x;
1055 weights[3 * i + 1] = particle->weights.y;
1056 weights[3 * i + 2] = particle->weights.z;
1057 }
1058
1059 // --- 4. Cell Indices (i, j, k) ---
1060 if (cellIndices) {
1061 cellIndices[3 * i + 0] = particle->cell[0];
1062 cellIndices[3 * i + 1] = particle->cell[1];
1063 cellIndices[3 * i + 2] = particle->cell[2];
1064 }
1065
1066 // --- 5. Status ---
1067 if (status) {
1068 status[i] = (PetscInt)particle->location_status;
1069 }
1070
1071 // --- 6. Diffusivity ---
1072 if (diffusivity) {
1073 diffusivity[i] = particle->diffusivity;
1074 }
1075
1076 if(diffusivitygradient){
1077 diffusivitygradient[i].x = particle->diffusivitygradient.x;
1078 diffusivitygradient[i].y = particle->diffusivitygradient.y;
1079 diffusivitygradient[i].z = particle->diffusivitygradient.z;
1080 }
1081 // --- 7. Psi ---
1082 if (psi) {
1083 psi[i] = particle->psi;
1084 }
1085
1086 // LOG_LOOP_ALLOW(LOCAL, LOG_VERBOSE, i, 1000, "Updated fields for Particle [%d].\n", i);
1087
1089 PetscFunctionReturn(0);
1090}
1091
1092#undef __FUNCT__
1093#define __FUNCT__ "IsParticleInsideBoundingBox"
1094/**
1095 * @brief Checks if a particle's location is within a specified bounding box.
1096 *
1097 * This function determines whether the given particle's location lies inside the provided bounding box.
1098 * It performs an axis-aligned bounding box (AABB) check by comparing the particle's coordinates to the
1099 * minimum and maximum coordinates of the bounding box in each dimension (x, y, z).
1100 *
1101 * logging statements are included to provide detailed information about the function's execution.
1102 *
1103 * @param[in] bbox Pointer to the BoundingBox structure containing minimum and maximum coordinates.
1104 * @param[in] particle Pointer to the Particle structure containing the particle's location and identifier.
1105 *
1106 * @return PetscBool Returns `PETSC_TRUE` if the particle is inside the bounding box, `PETSC_FALSE` otherwise.
1107 *
1108 * @note
1109 * - The function assumes that the `bbox` and `particle` pointers are valid and non-NULL.
1110 * - The function includes logging statements that start with the function name.
1111 * - The `LOG_ALLOW_SCOPE` variable is used to distinguish between `GLOBAL` and `LOCAL` LOG_ALLOW outputs.
1112 * - Be cautious when logging in performance-critical code sections, especially if the function is called frequently.
1113 */
1114PetscBool IsParticleInsideBoundingBox(const BoundingBox *bbox, const Particle *particle)
1115{
1116 PetscFunctionBeginUser;
1118
1119 // Validate input pointers
1120 if (!bbox) {
1121 // LOG_ALLOW error message and return PETSC_FALSE
1122 LOG_ALLOW(LOCAL,LOG_ERROR, "Error - 'bbox' pointer is NULL.");
1124 return PETSC_FALSE;
1125 }
1126 if (!particle) {
1127 LOG_ALLOW(LOCAL,LOG_ERROR,"Error - 'particle' pointer is NULL.");
1129 return PETSC_FALSE;
1130 }
1131
1132 // Extract particle location and bounding box coordinates
1133 const Cmpnts loc = particle->loc;
1134 const Cmpnts min_coords = bbox->min_coords;
1135 const Cmpnts max_coords = bbox->max_coords;
1136
1137 // LOG_ALLOW the particle location and bounding box coordinates for debugging
1138 LOG_ALLOW_SYNC(LOCAL, LOG_VERBOSE, "Particle PID %ld location: (%.6f, %.6f, %.6f).\n",particle->PID, loc.x, loc.y, loc.z);
1139 LOG_ALLOW_SYNC(LOCAL, LOG_VERBOSE, "BoundingBox min_coords: (%.6f, %.6f, %.6f), max_coords: (%.6f, %.6f, %.6f).\n",
1140 min_coords.x, min_coords.y, min_coords.z, max_coords.x, max_coords.y, max_coords.z);
1141
1142 // Check if the particle's location is within the bounding box
1143 if ((loc.x >= min_coords.x && loc.x <= max_coords.x) &&
1144 (loc.y >= min_coords.y && loc.y <= max_coords.y) &&
1145 (loc.z >= min_coords.z && loc.z <= max_coords.z)) {
1146 // Particle is inside the bounding box
1147 LOG_ALLOW_SYNC(LOCAL,LOG_VERBOSE, "Particle PID %ld is inside the bounding box.\n",particle->PID);
1149 return PETSC_TRUE;
1150 }
1151
1152 // Particle is outside the bounding box
1153 LOG_ALLOW_SYNC(LOCAL, LOG_VERBOSE,"Particle PID %ld is outside the bounding box.\n",particle->PID);
1155 return PETSC_FALSE;
1156}
1157
1158
1159#undef __FUNCT__
1160#define __FUNCT__ "UpdateParticleWeights"
1161/**
1162 * @brief Updates a particle's interpolation weights based on distances to cell faces.
1163 *
1164 * This function computes interpolation weights using distances to the six
1165 * cell faces (`d`) and updates the `weight` field of the provided particle.
1166 *
1167 * @param[in] d Pointer to an array of distances to the six cell faces.
1168 * @param[out] particle Pointer to the Particle structure whose weights are to be updated.
1169 *
1170 * @return PetscErrorCode Returns 0 on success, or a non-zero error code on failure.
1171 */
1172PetscErrorCode UpdateParticleWeights(PetscReal *d, Particle *particle) {
1173
1174 PetscFunctionBeginUser;
1176
1177 // Validate input pointers
1178 if (!d || !particle) {
1179 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL,
1180 "Null pointer argument (d or particle).");
1181 }
1182
1183
1184 // Validate distances
1185 for (PetscInt i = LEFT; i < NUM_FACES; i++) {
1188 "face distance d[%d] = %f <= %f; "
1189 "clamping to 1e-14 to avoid zero/negative.\n",
1190 i, (double)d[i], INTERPOLATION_DISTANCE_TOLERANCE);
1192 }
1193 }
1194
1195 // LOG_ALLOW the input distances
1197 "Calculating weights with distances: "
1198 "[LEFT=%f, RIGHT=%f, BOTTOM=%f, TOP=%f, FRONT=%f, BACK=%f].\n",
1199 d[LEFT], d[RIGHT], d[BOTTOM], d[TOP], d[FRONT], d[BACK]);
1200
1201 // Compute and update the particle's weights
1202 particle->weights.x = d[LEFT] / (d[LEFT] + d[RIGHT]);
1203 particle->weights.y = d[BOTTOM] / (d[BOTTOM] + d[TOP]);
1204 particle->weights.z = d[BACK] / (d[FRONT] + d[BACK]);
1205
1206 // LOG_ALLOW the updated weights
1208 "Updated particle weights: x=%f, y=%f, z=%f.\n",
1209 particle->weights.x, particle->weights.y, particle->weights.z);
1210
1211
1213 PetscFunctionReturn(0);
1214}
1215
1216/**
1217 * @brief Initializes or loads the particle swarm based on the simulation context.
1218 *
1219 * This function is the central point for setting up the DMSwarm. Its behavior
1220 * depends on the simulation context (simCtx):
1221 *
1222 * 1. **Fresh Start (simCtx->StartStep == 0):** A new particle population is
1223 * generated according to the specified initial conditions.
1224 *
1225 * 2. **Restart (simCtx->StartStep > 0):**
1226 * - If `simCtx->particleRestartMode` is "init", a new particle population
1227 * is generated, just like a fresh start. This allows injecting fresh
1228 * particles into a pre-computed flow field.
1229 * - If `simCtx->particleRestartMode` is "load", the particle state is loaded
1230 * from restart files corresponding to the StartStep.
1231 *
1232 * @param[in,out] simCtx Pointer to the main SimulationContext, which contains all
1233 * configuration and provides access to the UserCtx.
1234 * @return PetscErrorCode Returns 0 on success, non-zero on failure.
1235 */
1236#undef __FUNCT__
1237#define __FUNCT__ "InitializeParticleSwarm"
1238PetscErrorCode InitializeParticleSwarm(SimCtx *simCtx)
1239{
1240 PetscErrorCode ierr;
1241 PetscInt particlesPerProcess = 0;
1242 UserCtx *user = simCtx->usermg.mgctx[simCtx->usermg.mglevels - 1].user;
1243
1244 PetscFunctionBeginUser;
1246
1247 LOG_ALLOW(GLOBAL, LOG_INFO, "Starting particle swarm setup for %d particles.\n", simCtx->np);
1248
1249 // --- Phase 1: Create the DMSwarm Object (Always required) ---
1250 // This creates the container and registers the fields. It does not add particles yet.
1251 ierr = CreateParticleSwarm(user, simCtx->np, &particlesPerProcess, simCtx->bboxlist); CHKERRQ(ierr);
1252 LOG_ALLOW(GLOBAL, LOG_INFO, "DMSwarm object and fields created successfully.\n");
1253
1254
1255 // --- Phase 2: Decide whether to Initialize new particles or Load existing ones ---
1256 PetscBool should_initialize_new_particles = PETSC_FALSE;
1257 if(simCtx->exec_mode == EXEC_MODE_POSTPROCESSOR){
1258 should_initialize_new_particles = PETSC_TRUE;
1259 }else{
1260 if (simCtx->StartStep == 0) {
1261 should_initialize_new_particles = PETSC_TRUE; // Standard fresh start
1262 } else {
1263 // It's a restart, so check the user's requested particle mode.
1264 if (strcmp(simCtx->particleRestartMode, "init") == 0) {
1265 should_initialize_new_particles = PETSC_TRUE; // User wants to re-initialize particles in a restarted flow.
1266 }
1267 }
1268 }
1269
1270 // --- Phase 3: Execute the chosen particle setup path ---
1271 if (should_initialize_new_particles) {
1272 // --- PATH A: Generate a fresh population of particles ---
1273 LOG_ALLOW(GLOBAL, LOG_INFO, "Mode: INITIALIZE. Generating new particle population.\n");
1274 PetscRandom randx, randy, randz;
1275 PetscRandom rand_logic_i, rand_logic_j, rand_logic_k;
1276
1277 ierr = InitializeRandomGenerators(user, &randx, &randy, &randz); CHKERRQ(ierr);
1278 ierr = InitializeLogicalSpaceRNGs(&rand_logic_i, &rand_logic_j, &rand_logic_k); CHKERRQ(ierr);
1279 ierr = AssignInitialPropertiesToSwarm(user, particlesPerProcess, &randx, &randy, &randz, &rand_logic_i, &rand_logic_j, &rand_logic_k, simCtx->bboxlist); CHKERRQ(ierr);
1280 ierr = FinalizeSwarmSetup(&randx, &randy, &randz, &rand_logic_i, &rand_logic_j, &rand_logic_k); CHKERRQ(ierr);
1281
1282 } else {
1283 // --- PATH B: Load particle population from restart files ---
1284 // This path is only taken if simCtx->StartStep > 0 AND simCtx->particleRestartMode == "load"
1285 LOG_ALLOW(GLOBAL, LOG_INFO, "Mode: LOAD. Loading particle population from files for step %d.\n", simCtx->StartStep);
1286
1287 ierr = PreCheckAndResizeSwarm(user, simCtx->StartStep, "dat"); CHKERRQ(ierr);
1288
1289 ierr = ReadAllSwarmFields(user, simCtx->StartStep); CHKERRQ(ierr);
1290 // Note: We check for file-open errors inside ReadAllSwarmFields now.
1291
1292 LOG_ALLOW(GLOBAL, LOG_INFO, "Particle data loaded. CellID and status are preserved from file.\n");
1293 }
1294
1295 LOG_ALLOW_SYNC(GLOBAL, LOG_INFO, "Particle swarm setup complete.\n");
1296
1298 PetscFunctionReturn(0);
1299}
1300
PetscErrorCode GetRandomCellAndLogicalCoordsOnInletFace(UserCtx *user, const DMDALocalInfo *info, PetscInt xs_gnode_rank, PetscInt ys_gnode_rank, PetscInt zs_gnode_rank, PetscInt IM_nodes_global, PetscInt JM_nodes_global, PetscInt KM_nodes_global, PetscRandom *rand_logic_i_ptr, PetscRandom *rand_logic_j_ptr, PetscRandom *rand_logic_k_ptr, PetscInt *ci_metric_lnode_out, PetscInt *cj_metric_lnode_out, PetscInt *ck_metric_lnode_out, PetscReal *xi_metric_logic_out, PetscReal *eta_metric_logic_out, PetscReal *zta_metric_logic_out)
Assuming the current rank services the inlet face, this function selects a random cell (owned by this...
Definition Boundaries.c:467
PetscErrorCode CanRankServiceInletFace(UserCtx *user, const DMDALocalInfo *info, PetscInt IM_nodes_global, PetscInt JM_nodes_global, PetscInt KM_nodes_global, PetscBool *can_service_inlet_out)
Determines if the current MPI rank owns any part of the globally defined inlet face,...
Definition Boundaries.c:25
PetscErrorCode GetDeterministicFaceGridLocation(UserCtx *user, const DMDALocalInfo *info, PetscInt xs_gnode_rank, PetscInt ys_gnode_rank, PetscInt zs_gnode_rank, PetscInt IM_cells_global, PetscInt JM_cells_global, PetscInt KM_cells_global, PetscInt64 particle_global_id, PetscInt *ci_metric_lnode_out, PetscInt *cj_metric_lnode_out, PetscInt *ck_metric_lnode_out, PetscReal *xi_metric_logic_out, PetscReal *eta_metric_logic_out, PetscReal *zta_metric_logic_out, PetscBool *placement_successful_out)
Places particles in a deterministic grid/raster pattern on a specified domain face.
Definition Boundaries.c:268
PetscErrorCode MetricLogicalToPhysical(UserCtx *user, const Cmpnts ***X, PetscInt i, PetscInt j, PetscInt k, PetscReal xi, PetscReal eta, PetscReal zta, Cmpnts *Xp)
Public wrapper: map (cell index, ξ,η,ζ) to (x,y,z).
Definition Metric.c:77
PetscErrorCode PreCheckAndResizeSwarm(UserCtx *user, PetscInt ti, const char *ext)
Checks particle count in the reference file and resizes the swarm if needed.
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.
static PetscErrorCode AssignInitialFieldToSwarm(UserCtx *user, const char *fieldName, PetscInt fieldDim)
Initializes a generic swarm field with point-level updates.
#define INTERPOLATION_DISTANCE_TOLERANCE
static PetscErrorCode InitializeSwarmFieldValue(const char *fieldName, PetscInt p, PetscInt fieldDim, PetscReal *fieldData)
Helper function to Initialize a given particle’s field value.
PetscBool IsParticleInsideBoundingBox(const BoundingBox *bbox, const Particle *particle)
Checks if a particle's location is within a specified bounding box.
static PetscErrorCode InitializeParticleBasicProperties(UserCtx *user, PetscInt particlesPerProcess, PetscRandom *rand_logic_i, PetscRandom *rand_logic_j, PetscRandom *rand_logic_k, BoundingBox *bboxlist)
Initializes basic properties for particles on the local process.
PetscErrorCode UnpackSwarmFields(PetscInt i, const PetscInt64 *PIDs, const PetscReal *weights, const PetscReal *positions, const PetscInt *cellIndices, PetscReal *velocities, PetscInt *LocStatus, PetscReal *diffusivity, Cmpnts *diffusivitygradient, PetscReal *psi, 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 *positions, PetscReal *velocities, PetscReal *weights, PetscInt *cellIndices, PetscInt *status, PetscReal *diffusivity, Cmpnts *diffusivitygradient, PetscReal *psi)
Updates DMSwarm data arrays 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.
static PetscErrorCode DetermineVolumetricInitializationParameters(UserCtx *user, DMDALocalInfo *info, PetscInt xs_gnode, PetscInt ys_gnode, PetscInt zs_gnode, PetscRandom *rand_logic_i_ptr, PetscRandom *rand_logic_j_ptr, PetscRandom *rand_logic_k_ptr, PetscInt *ci_metric_lnode_out, PetscInt *cj_metric_lnode_out, PetscInt *ck_metric_lnode_out, PetscReal *xi_metric_logic_out, PetscReal *eta_metric_logic_out, PetscReal *zta_metric_logic_out, PetscBool *can_place_in_volume_out)
Determines cell selection and intra-cell logical coordinates for volumetric initialization (Mode 1).
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.
Header file for Particle Swarm management functions.
PetscErrorCode ReadAllSwarmFields(UserCtx *user, PetscInt ti)
Reads multiple fields (positions, velocity, CellID, and weight) into a DMSwarm.
Definition io.c:1534
#define LOG_LOOP_ALLOW(scope, level, iterVar, interval, fmt,...)
Logs a message inside a loop, but only every interval iterations.
Definition logging.h:312
PetscBool is_function_allowed(const char *functionName)
Checks if a given function is in the allow-list.
Definition logging.c:157
#define LOG_ALLOW_SYNC(scope, level, fmt,...)
----— DEBUG ---------------------------------------— #define LOG_ALLOW(scope, level,...
Definition logging.h:267
#define LOCAL
Logging scope definitions for controlling message output.
Definition logging.h:45
#define GLOBAL
Scope for global logging across all processes.
Definition logging.h:46
const char * BCFaceToString(BCFace face)
Helper function to convert BCFace enum to a string representation.
Definition logging.c:643
#define LOG_ALLOW(scope, level, fmt,...)
Logging macro that checks both the log level and whether the calling function is in the allowed-funct...
Definition logging.h:200
#define PROFILE_FUNCTION_END
Marks the end of a profiled code block.
Definition logging.h:746
LogLevel get_log_level()
Retrieves the current logging level from the environment variable LOG_LEVEL.
Definition logging.c:39
@ LOG_ERROR
Critical errors that may halt the program.
Definition logging.h:28
@ LOG_INFO
Informational messages about program execution.
Definition logging.h:31
@ LOG_WARNING
Non-critical issues that warrant attention.
Definition logging.h:29
@ LOG_DEBUG
Detailed debugging information.
Definition logging.h:32
@ LOG_VERBOSE
Extremely detailed logs, typically for development use only.
Definition logging.h:34
#define PROFILE_FUNCTION_BEGIN
Marks the beginning of a profiled code block (typically a function).
Definition logging.h:737
const char * ParticleInitializationToString(ParticleInitializationType ParticleInitialization)
Helper function to convert ParticleInitialization to a string representation.
Definition logging.c:675
PetscErrorCode InitializeRandomGenerators(UserCtx *user, PetscRandom *randx, PetscRandom *randy, PetscRandom *randz)
Initializes random number generators for assigning particle properties.
Definition setup.c:2642
PetscErrorCode GetOwnedCellRange(const DMDALocalInfo *info_nodes, PetscInt dim, PetscInt *xs_cell_global, PetscInt *xm_cell_local)
Determines the global starting index and number of CELLS owned by the current processor in a specifie...
Definition setup.c:1756
PetscErrorCode InitializeLogicalSpaceRNGs(PetscRandom *rand_logic_i, PetscRandom *rand_logic_j, PetscRandom *rand_logic_k)
Initializes random number generators for logical space operations [0.0, 1.0).
Definition setup.c:2693
Cmpnts vel
Definition variables.h:169
UserCtx * user
Definition variables.h:474
PetscBool inletFaceDefined
Definition variables.h:747
BCFace identifiedInletBCFace
Definition variables.h:748
PetscInt cell[3]
Definition variables.h:167
SimCtx * simCtx
Back-pointer to the master simulation context.
Definition variables.h:731
@ PARTICLE_INIT_SURFACE_RANDOM
Random placement on the inlet face.
Definition variables.h:466
@ PARTICLE_INIT_SURFACE_EDGES
Deterministic placement at inlet face edges.
Definition variables.h:469
@ PARTICLE_INIT_POINT_SOURCE
All particles at a fixed (psrc_x,psrc_y,psrc_z) — for validation.
Definition variables.h:468
@ PARTICLE_INIT_VOLUME
Random volumetric distribution across the domain.
Definition variables.h:467
ParticleLocationStatus
Defines the state of a particle with respect to its location and migration status during the iterativ...
Definition variables.h:135
@ UNINITIALIZED
Definition variables.h:140
PetscReal CMy_c
Definition variables.h:642
PetscReal Min_X
Definition variables.h:738
UserMG usermg
Definition variables.h:698
PetscReal psrc_x
Definition variables.h:643
Cmpnts max_coords
Maximum x, y, z coordinates of the bounding box.
Definition variables.h:156
Cmpnts diffusivitygradient
Definition variables.h:174
PetscInt np
Definition variables.h:676
PetscReal Max_Y
Definition variables.h:738
PetscInt StartStep
Definition variables.h:595
Cmpnts min_coords
Minimum x, y, z coordinates of the bounding box.
Definition variables.h:155
PetscScalar x
Definition variables.h:101
Cmpnts loc
Definition variables.h:168
PetscMPIInt destination_rank
Definition variables.h:172
PetscReal psrc_z
Point source location for PARTICLE_INIT_POINT_SOURCE.
Definition variables.h:643
ParticleLocationStatus location_status
Definition variables.h:171
char particleRestartMode[16]
Definition variables.h:681
BoundingBox * bboxlist
Definition variables.h:679
PetscReal CMz_c
Definition variables.h:642
ParticleInitializationType ParticleInitialization
Definition variables.h:680
PetscScalar z
Definition variables.h:101
PetscInt mglevels
Definition variables.h:481
PetscReal Min_Z
Definition variables.h:738
PetscReal psrc_y
Definition variables.h:643
PetscReal Max_X
Definition variables.h:738
PetscReal Min_Y
Definition variables.h:738
PetscReal diffusivity
Definition variables.h:173
PetscReal psi
Definition variables.h:175
PetscScalar y
Definition variables.h:101
@ EXEC_MODE_POSTPROCESSOR
Definition variables.h:559
Cmpnts weights
Definition variables.h:170
@ TOP
Definition variables.h:145
@ FRONT
Definition variables.h:145
@ BOTTOM
Definition variables.h:145
@ BACK
Definition variables.h:145
@ LEFT
Definition variables.h:145
@ NUM_FACES
Definition variables.h:145
@ RIGHT
Definition variables.h:145
MGCtx * mgctx
Definition variables.h:484
ExecutionMode exec_mode
Definition variables.h:603
BoundingBox bbox
Definition variables.h:739
PetscReal Max_Z
Definition variables.h:738
PetscInt64 PID
Definition variables.h:166
PetscInt LoggingFrequency
Definition variables.h:703
PetscReal CMx_c
Definition variables.h:642
BCFace
Identifies the six logical faces of a structured computational block.
Definition variables.h:203
Defines a 3D axis-aligned bounding box.
Definition variables.h:154
A 3D point or vector with PetscScalar components.
Definition variables.h:100
Defines a particle's core properties for Lagrangian tracking.
Definition variables.h:165
The master context for the entire simulation.
Definition variables.h:585
User-defined context containing data specific to a single computational grid level.
Definition variables.h:728