PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
Loading...
Searching...
No Matches
simulation.c
Go to the documentation of this file.
1/**
2 * @file simulation.c // code for simulation loop
3 * @brief Test program for DMSwarm interpolation using the fdf-curvIB method.
4 * Provides the setup to start any simulation with DMSwarm and DMDAs.
5 **/
6
7#include "simulation.h"
8
9/**
10 * @brief Copies the current time step's solution fields into history vectors
11 * (e.g., U(t_n) -> U_o, U_o -> U_rm1) for the next time step's calculations.
12 *
13 * This function is critical for multi-step time integration schemes (like BDF2)
14 * used by the legacy solver. It must be called at the end of every time step,
15 * after the new solution has been fully computed.
16 *
17 * The order of operations is important to avoid overwriting data prematurely.
18 *
19 * @param user The UserCtx for a single block. The function modifies the history
20 * vectors (Ucont_o, Ucont_rm1, etc.) within this context.
21 * @return PetscErrorCode 0 on success.
22 */
24{
25 PetscErrorCode ierr;
26 SimCtx *simCtx = user->simCtx; // Access global settings if needed
27
28 PetscFunctionBeginUser;
29 LOG_ALLOW(LOCAL, LOG_DEBUG, "Rank %d, Block %d: Updating solver history vectors.\n",
30 simCtx->rank, user->_this);
31
32 // --- Primary Contravariant Velocity History ---
33 // The order is critical here.
34 // 1. First, move the n-1 state (Ucont_o) to the n-2 slot (Ucont_rm1).
35 ierr = VecCopy(user->Ucont_o, user->Ucont_rm1); CHKERRQ(ierr);
36 // 2. Then, move the new n state (Ucont) to the n-1 slot (Ucont_o).
37 ierr = VecCopy(user->Ucont, user->Ucont_o); CHKERRQ(ierr);
38
39 LOG_ALLOW(LOCAL,LOG_DEBUG, "Rank %d, Block %d, Ucont history updated.\n",simCtx->rank,user->_this);
40
41 // --- Update History for Other Fields ---
42 // These are typically only needed at the n-1 state.
43 ierr = VecCopy(user->Ucat, user->Ucat_o); CHKERRQ(ierr);
44 ierr = VecCopy(user->P, user->P_o); CHKERRQ(ierr);
45 LOG_ALLOW(LOCAL,LOG_DEBUG, "Rank %d, Block %d, Ucat & P history updated.\n",simCtx->rank,user->_this);
46
47 if (simCtx->immersed) {
48 ierr = VecCopy(user->Nvert, user->Nvert_o); CHKERRQ(ierr);
49 }
50
51 // --- Update History for Turbulence Models (if active) ---
52 if (simCtx->rans) {
53 ierr = VecCopy(user->K_Omega, user->K_Omega_o); CHKERRQ(ierr);
54 }
55
56 // --- Synchronize Local Ghost Regions for the new history vectors ---
57 // This is essential so that stencils in the next time step's calculations
58 // have correct values from neighboring processes.
59 ierr = DMGlobalToLocalBegin(user->fda, user->Ucont_o, INSERT_VALUES, user->lUcont_o); CHKERRQ(ierr);
60 ierr = DMGlobalToLocalEnd(user->fda, user->Ucont_o, INSERT_VALUES, user->lUcont_o); CHKERRQ(ierr);
61
62 ierr = DMGlobalToLocalBegin(user->fda, user->Ucont_rm1, INSERT_VALUES, user->lUcont_rm1); CHKERRQ(ierr);
63 ierr = DMGlobalToLocalEnd(user->fda, user->Ucont_rm1, INSERT_VALUES, user->lUcont_rm1); CHKERRQ(ierr);
64
65 if (simCtx->immersed) {
66 ierr = DMGlobalToLocalBegin(user->da, user->Nvert_o, INSERT_VALUES, user->lNvert_o); CHKERRQ(ierr);
67 ierr = DMGlobalToLocalEnd(user->da, user->Nvert_o, INSERT_VALUES, user->lNvert_o); CHKERRQ(ierr);
68 }
69
70 if (simCtx->rans) {
71 ierr = DMGlobalToLocalBegin(user->fda2, user->K_Omega_o, INSERT_VALUES, user->lK_Omega_o); CHKERRQ(ierr);
72 ierr = DMGlobalToLocalEnd(user->fda2, user->K_Omega_o, INSERT_VALUES, user->lK_Omega_o); CHKERRQ(ierr);
73 }
74
75 PetscFunctionReturn(0);
76}
77
78#undef __FUNCT__
79#define __FUNCT__ "PerformInitialSetup"
80/**
81 * @brief Finalizes the simulation setup at t=0, ensuring a consistent state before time marching.
82 *
83 * This function is called from main() after the initial Eulerian and Lagrangian states have been
84 * created but before the main time loop begins. Its responsibilities are:
85 *
86 * 1. Settling the particle swarm: Migrates particles to their correct owner ranks and finds their
87 * initial host cells. This includes handling special surface initializations.
88 * 2. Coupling the fields: Interpolates the initial Eulerian fields to the settled particle locations.
89 * 3. Preparing for the first step: Scatters particle data back to the grid.
90 * 4. Writing the initial output for step 0.
91 *
92 * @param simCtx Pointer to the main simulation context structure.
93 * @return PetscErrorCode 0 on success, non-zero on failure.
94 */
96{
97 PetscErrorCode ierr;
98 // --- Get pointers from SimCtx instead of passing them as arguments ---
99 UserCtx *user = simCtx->usermg.mgctx[simCtx->usermg.mglevels-1].user;
100 BoundingBox *bboxlist = simCtx->bboxlist;
101
102 PetscFunctionBeginUser;
103
104 LOG_ALLOW(GLOBAL, LOG_INFO, "[T=%.4f, Step=%d] Performing initial particle setup procedures.\n", simCtx->ti, simCtx->step);
105
106 // --- 1. Initial Particle Settlement (Location and Migration) ---
107 LOG_ALLOW(GLOBAL, LOG_INFO, "[T=%.4f, Step=%d] Initial Settlement: Locating and migrating all particles...\n", simCtx->ti, simCtx->step);
108 ierr = LocateAllParticlesInGrid(user, bboxlist); CHKERRQ(ierr);
109
111 LOG_ALLOW(GLOBAL, LOG_DEBUG, "[T=%.4f, Step=%d] Particle field states after Initial settlement...\n", simCtx->ti, simCtx->step);
112 ierr = LOG_PARTICLE_FIELDS(user,simCtx->LoggingFrequency); CHKERRQ(ierr);
113 }
114
115 // --- 2. Re-initialize Particles on Inlet Surface (if applicable) ---
116 if ((simCtx->ParticleInitialization == 0 || simCtx->ParticleInitialization == 3) && user->inletFaceDefined) {
117 LOG_ALLOW(GLOBAL, LOG_INFO, "[T=%.4f, Step=%d] Re-initializing particles on inlet surface...\n", simCtx->ti, simCtx->step);
118 ierr = ReinitializeParticlesOnInletSurface(user, simCtx->ti, simCtx->step); CHKERRQ(ierr);
119
120 LOG_ALLOW(GLOBAL, LOG_INFO, "[T=%.4f, Step=%d] Resetting statuses for post-reinitialization settlement.\n", simCtx->ti, simCtx->step);
121 ierr = ResetAllParticleStatuses(user); CHKERRQ(ierr);
122
123 LOG_ALLOW(GLOBAL, LOG_INFO, "[T=%.4f, Step=%d] Post-Reinitialization Settlement...\n", simCtx->ti, simCtx->step);
124 ierr = LocateAllParticlesInGrid(user, bboxlist); CHKERRQ(ierr);
125
126 }
127
128 // --- 3. Finalize State for t=0 ---
129 LOG_ALLOW(GLOBAL, LOG_INFO, "[T=%.4f, Step=%d] Interpolating initial fields to settled particles.\n", simCtx->ti, simCtx->step);
130 ierr = InterpolateAllFieldsToSwarm(user); CHKERRQ(ierr);
131 //ierr = ScatterAllParticleFieldsToEulerFields(user); CHKERRQ(ierr);
132
133 // --- 4. Initial History and Output ---
134 // Update solver history vectors with the t=0 state before the first real step
135 for (PetscInt bi = 0; bi < simCtx->block_number; bi++) {
136 ierr = UpdateSolverHistoryVectors(&user[bi]); CHKERRQ(ierr);
137 }
138
139 if (simCtx->OutputFreq > 0 || (simCtx->StepsToRun == 0 && simCtx->StartStep == 0)) {
140 LOG_ALLOW(GLOBAL, LOG_INFO, "[T=%.4f, Step=%d] Writing initial simulation data.\n", simCtx->ti, simCtx->step);
141
142 // --- Particle Output (assumes functions operate on the master user context) ---
143 if(get_log_level() == LOG_DEBUG){
144 LOG(GLOBAL, LOG_DEBUG, "[T=%.4f, Step=%d] Initial Particle field.\n", simCtx->ti, simCtx->step);
145 ierr = LOG_PARTICLE_FIELDS(user,simCtx->LoggingFrequency); CHKERRQ(ierr);
146 }
147 ierr = WriteAllSwarmFields(user); CHKERRQ(ierr);
148
149 // --- Eulerian Field Output (MUST loop over all blocks) --- // <<< CHANGED/FIXED
150 for (PetscInt bi = 0; bi < simCtx->block_number; bi++) {
151 ierr = WriteSimulationFields(&user[bi]); CHKERRQ(ierr);
152 }
153 }
154
155 LOG_ALLOW(GLOBAL, LOG_INFO, "--- Initial setup complete. Ready for time marching. ---\n");
156 PetscFunctionReturn(0);
157}
158
159#undef __FUNCT__
160#define __FUNCT__ "PerformLoadedParticleSetup"
161/**
162 * @brief Finalizes the simulation state after particle and fluid data have been loaded from a restart.
163 *
164 * This helper function performs the critical sequence of operations required to ensure
165 * the loaded Lagrangian and Eulerian states are fully consistent and the solver is
166 * ready to proceed. This includes:
167 * 1. Verifying particle locations in the grid and building runtime links.
168 * 2. Synchronizing particle velocity with the authoritative grid velocity via interpolation.
169 * 3. Scattering particle source terms (e.g., volume fraction) back to the grid.
170 * 4. Updating the solver's history vectors with the final, fully-coupled state.
171 * 5. Writing the complete, consistent state to output files for the restart step.
172 *
173 * @param simCtx The main simulation context.
174 * @return PetscErrorCode 0 on success.
175 */
176PetscErrorCode PerformLoadedParticleSetup(SimCtx *simCtx)
177{
178 PetscErrorCode ierr;
179 PetscFunctionBeginUser;
180 UserCtx *user = simCtx->usermg.mgctx[simCtx->usermg.mglevels-1].user;
181
182 // 0. Reset all particle statuses to ensure proper location checks.
183 ierr = ResetAllParticleStatuses(user); CHKERRQ(ierr);
184
185 // 1. Rebuild grid-to-particle links based on loaded coordinates.
186 ierr = LocateAllParticlesInGrid(user, simCtx->bboxlist); CHKERRQ(ierr);
187
188 if(get_log_level() == LOG_DEBUG){
189 LOG(GLOBAL, LOG_DEBUG, "[T=%.4f, Step=%d] Particle field states after locating loaded particles...\n", simCtx->ti, simCtx->step);
190 ierr = LOG_PARTICLE_FIELDS(user,simCtx->LoggingFrequency); CHKERRQ(ierr);
191 }
192
193 LOG_ALLOW(GLOBAL, LOG_INFO, "[T=%.4f, Step=%d] Interpolating initial fields to settled particles.\n", simCtx->ti, simCtx->step);
194
195 // 2. Ensure particles have velocity from the authoritative loaded grid for consistency.
196 ierr = InterpolateAllFieldsToSwarm(user); CHKERRQ(ierr);
197
198 // 3. Update Eulerian source terms from the loaded particle data.
199 //ierr = ScatterAllParticleFieldsToEulerFields(user); CHKERRQ(ierr);
200
201 // --- 4. Initial History and Output ---
202 // Update solver history vectors with the t=0 state before the first real step
203 for (PetscInt bi = 0; bi < simCtx->block_number; bi++) {
204 ierr = UpdateSolverHistoryVectors(&user[bi]); CHKERRQ(ierr);
205 }
206
207 if (simCtx->OutputFreq > 0 || (simCtx->StepsToRun == 0 && simCtx->StartStep == 0)) {
208 LOG_ALLOW(GLOBAL, LOG_INFO, "[T=%.4f, Step=%d] Writing initial simulation data.\n", simCtx->ti, simCtx->step);
209
210 // --- Particle Output (assumes functions operate on the master user context) ---
211 if(get_log_level() >=LOG_INFO) ierr = LOG_PARTICLE_FIELDS(user, simCtx->LoggingFrequency); CHKERRQ(ierr);
212 ierr = WriteAllSwarmFields(user); CHKERRQ(ierr);
213
214 // --- Eulerian Field Output (MUST loop over all blocks) --- // <<< CHANGED/FIXED
215 for (PetscInt bi = 0; bi < simCtx->block_number; bi++) {
216 ierr = WriteSimulationFields(&user[bi]); CHKERRQ(ierr);
217 }
218 }
219
220 LOG_ALLOW(GLOBAL, LOG_INFO, "--- Initial setup complete. Ready for time marching. ---\n");
221 PetscFunctionReturn(0);
222}
223
224#undef __FUNCT__
225#define __FUNCT__ "FinalizeRestartState"
226/**
227 * @brief Performs post-load/post-init consistency checks for a restarted simulation.
228 *
229 * This function is called from main() ONLY when a restart is being performed
230 * (i.e., StartStep > 0). It inspects the particle restart mode to determine the
231 * correct finalization procedure for the Lagrangian swarm.
232 *
233 * - If particles were loaded from a file (`mode == "load"`), it verifies their
234 * locations within the grid to establish necessary runtime links.
235 * - If new particles were initialized into the restarted flow (`mode == "init"`),
236 * it runs the full `PerformInitialSetup` sequence to migrate, locate, and
237 * couple the new particles with the existing fluid state.
238 *
239 * @param simCtx The main simulation context.
240 * @return PetscErrorCode 0 on success.
241 */
242PetscErrorCode FinalizeRestartState(SimCtx *simCtx)
243{
244 PetscErrorCode ierr;
245 UserCtx *user = simCtx->usermg.mgctx[simCtx->usermg.mglevels - 1].user;
246
247 PetscFunctionBeginUser;
248
249 LOG_ALLOW(GLOBAL, LOG_INFO, "--- Finalizing RESTART from state (step=%d, t=%.4f) ---\n", simCtx->StartStep, simCtx->ti);
250
251 // This function only needs to handle the particle finalization logic.
252 // The Eulerian state is assumed to be fully loaded and consistent at this point.
253 if (simCtx->np > 0) {
254
255 // Use the particle restart mode to decide the workflow.
256 if (strcmp(simCtx->particleRestartMode, "load") == 0) {
257 // PARTICLES WERE LOADED: The state is complete, but we must verify
258 // the loaded CellIDs and build the in-memory grid-to-particle links.
259 LOG_ALLOW(GLOBAL, LOG_INFO, "Particle Mode 'load': Verifying particle locations and building grid links...\n");
260 ierr = PerformLoadedParticleSetup(simCtx); CHKERRQ(ierr);
261
262 } else { // Mode must be "init"
263 // PARTICLES WERE RE-INITIALIZED: They need to be fully settled and coupled
264 // to the surrounding (restarted) fluid state.
265 LOG_ALLOW(GLOBAL, LOG_INFO, "Particle Mode 'init': Running full initial setup for new particles in restarted flow.\n");
266 ierr = PerformInitializedParticleSetup(simCtx); CHKERRQ(ierr);
267 }
268 } else {
269 LOG_ALLOW(GLOBAL, LOG_INFO, "No particles in simulation, restart finalization is complete.\n");
270
271 // Write the initial eulerian fields (this is done in PerformInitialSetup if particles exist.)
272 for(PetscInt bi = 0; bi < simCtx->block_number; bi ++){
273 ierr = WriteSimulationFields(&user[bi]); CHKERRQ(ierr);
274 }
275 }
276
277 LOG_ALLOW(GLOBAL, LOG_INFO, "--- Restart state successfully finalized. --\n");
278
279 PetscFunctionReturn(0);
280}
281
282#undef __FUNCT__
283#define __FUNCT__ "AdvanceSimulation"
284/**
285 * @brief Executes the main time-marching loop for the coupled Euler-Lagrange simulation.
286 *
287 * This function orchestrates the advancement of the simulation from the configured
288 * StartStep to the final step. It does NOT perform the initial t=0 setup, as that
289 * is handled by InitializeEulerianState and PerformInitialSetup in main().
290 *
291 * For each timestep, it performs the following sequence:
292 * 1. **Pre-Solver Actions:** Updates time-dependent boundary conditions (e.g., fluxin)
293 * and resets particle states for the new step.
294 * 2. **Eulerian Solve:** Calls the refactored legacy FlowSolver to advance the
295 * entire fluid field by one time step.
296 * 3. **Lagrangian Update:** Executes the full particle workflow: advection based
297 * on the previous step's velocity, followed by settling (location/migration)
298 * in the new grid, and finally interpolation of the new fluid velocity.
299 * 4. **Two-Way Coupling:** Scatters particle data back to the grid to act as source
300 * terms for the subsequent time step.
301 * 5. **History & I/O:** Updates the solver's history vectors and writes output files
302 * at the specified frequency.
303 *
304 * @param user Array of UserCtx structures for the finest grid level.
305 * @return PetscErrorCode 0 on success, non-zero on failure.
306 */
307PetscErrorCode AdvanceSimulation(SimCtx *simCtx)
308{
309 PetscErrorCode ierr;
310 // Get the master context from the first block. All blocks share it.
311 UserCtx *user = simCtx->usermg.mgctx[simCtx->usermg.mglevels-1].user;
312
313 // Retrieve control parameters from SimCtx for clarity
314 const PetscInt StartStep = simCtx->StartStep;
315 const PetscInt StepsToRun = simCtx->StepsToRun;
316 const PetscInt OutputFreq = simCtx->OutputFreq;
317 const PetscReal dt = simCtx->dt;
318
319 // Variables for particle removal statistics
320 PetscInt removed_local_ob, removed_global_ob;
321 PetscInt removed_local_lost, removed_global_lost;
322
323 PetscFunctionBeginUser;
325 LOG_ALLOW(GLOBAL, LOG_INFO, "Starting main time-marching loop: %d steps from step %d (t=%.4f), dt=%.4f\n",
326 StepsToRun, StartStep, simCtx->StartTime, dt);
327
328 // --- Main Time-Marching Loop ---
329 for (PetscInt step = StartStep; step < StartStep + StepsToRun; step++) {
330
331 // =================================================================
332 // 1. PRE-STEP SETUP
333 // =================================================================
334
335 // Update simulation time and step counters in the master context
336 simCtx->step = step + 1;
337 simCtx->ti += simCtx->dt; //simCtx->StartTime + step * simCtx->dt;
338
339 LOG_ALLOW(GLOBAL, LOG_INFO, "--- Advancing Step %d (To t=%.4f) ---\n", simCtx->step, simCtx->ti);
340
341
342 // For particles, reset their status to prepare for the new advection/location cycle
343 if (simCtx->np > 0) {
344 LOG_ALLOW(GLOBAL, LOG_DEBUG, "Resetting all particle statuses to NEEDS_LOCATION.\n");
345 ierr = ResetAllParticleStatuses(user); CHKERRQ(ierr);
346 }
347
348 // =================================================================
349 // 2. EULERIAN SOLVER STEP
350 // =================================================================
352 ierr = LOG_FIELD_ANATOMY(&user[0],"Coordinates","PreFlowSolver"); CHKERRQ(ierr);
353 }
354 LOG_ALLOW(GLOBAL, LOG_INFO, "Updating Eulerian Field ...\n");
355 if(strcmp(simCtx->eulerianSource,"load")==0){
356 //LOAD mode: Read pre-computed fields for the current step.
357 LOG_ALLOW(GLOBAL,LOG_INFO,"Eulerian Source 'load': Reading fields (t=%.4f,step=%d)...\n",simCtx->ti,simCtx->step);
358 for(PetscInt bi = 0; bi < simCtx->block_number;bi++){
359 ierr = ReadSimulationFields(&user[bi],simCtx->step); CHKERRQ(ierr);
360 }
361 }else if(strcmp(simCtx->eulerianSource,"solve")==0){
362 // SOLVE mode:Call the refactored, high-level legacy solver. This single function
363 // advances the entire multi-block fluid field from t_n to t_{n+1}.
364 LOG_ALLOW(GLOBAL,LOG_INFO,"Eulerian Source 'solve'. Updating Eulerian field via Solver...\n");
365 ierr = FlowSolver(simCtx); CHKERRQ(ierr);
366 }
367 LOG_ALLOW(GLOBAL, LOG_INFO, "Eulerian Field Updated ...\n");
369 LOG_ALLOW(GLOBAL, LOG_VERBOSE, "Post FlowSolver field states:\n");
370 ierr = LOG_FIELD_ANATOMY(&user[0],"Ucat","PostFlowSolver"); CHKERRQ(ierr);
371 ierr = LOG_FIELD_ANATOMY(&user[0],"P","PostFlowSolver"); CHKERRQ(ierr);
372 ierr = LOG_FIELD_ANATOMY(&user[0],"Ucont","PostFlowSolver"); CHKERRQ(ierr);
373 }
374
375
376 // =================================================================
377 // 3. LAGRANGIAN PARTICLE STEP
378 // =================================================================
379
380 if (simCtx->np > 0) {
381 LOG_ALLOW(GLOBAL, LOG_INFO, "Updating Lagrangian particle system...\n");
382
383 // a. Advect particles using the velocity interpolated from the *previous* step.
384 // P(t_{n+1}) = P(t_n) + V_p(t_n) * dt
385 ierr = UpdateAllParticlePositions(user); CHKERRQ(ierr);
386
387 // b. Settle all particles: find their new host cells and migrate them across ranks.
388 ierr = LocateAllParticlesInGrid(user, simCtx->bboxlist); CHKERRQ(ierr);
389
390 // c. Remove any particles that are now lost or out of the global domain.
391 ierr = CheckAndRemoveLostParticles(user, &removed_local_lost, &removed_global_lost); CHKERRQ(ierr);
392 ierr = CheckAndRemoveOutOfBoundsParticles(user, &removed_local_ob, &removed_global_ob, simCtx->bboxlist); CHKERRQ(ierr);
393 if (removed_global_lost + removed_global_ob > 0) {
394 LOG_ALLOW(GLOBAL, LOG_INFO, "Removed %d particles globally this step.\n", removed_global_lost + removed_global_ob);
395 }
396
397 // d. Interpolate the NEW fluid velocity (just computed by FlowSolver) onto the
398 // particles' new positions. This gives them V_p(t_{n+1}) for the *next* advection step.
399 ierr = InterpolateAllFieldsToSwarm(user); CHKERRQ(ierr);
400
401 // e. Update the Particle Fields (e.g., temperature, concentration) if applicable.
402 // This can be extended to include reactions, growth, etc.
403 ierr = UpdateAllParticleFields(user); CHKERRQ(ierr);
404
405 // f. (For Two-Way Coupling) Scatter particle data back to the grid to act as a source term.
406 // ierr = CalculateParticleCountPerCell(user); CHKERRQ(ierr);
407 // ierr = ScatterAllParticleFieldsToEulerFields(user); CHKERRQ(ierr);
408
410 LOG_ALLOW(GLOBAL, LOG_VERBOSE, "Post Lagrangian update field states:\n");
411 ierr = LOG_FIELD_MIN_MAX(&user[0],"Psi"); CHKERRQ(ierr);
412 }
413 }
414
415 // =================================================================
416 // 4. UPDATE HISTORY & I/O
417 // =================================================================
418
419 // Copy the newly computed fields (Ucont, P, etc.) to the history vectors
420 // (_o, _rm1) to prepare for the next time step.
421 for (PetscInt bi = 0; bi < simCtx->block_number; bi++) {
422 ierr = UpdateSolverHistoryVectors(&user[bi]); CHKERRQ(ierr);
423 }
424
425 //ierr = LOG_UCAT_ANATOMY(&user[0],"Final"); CHKERRQ(ierr);
426
427 // Handle periodic file output
428 if (OutputFreq > 0 && (step + 1) % OutputFreq == 0) {
429 LOG_ALLOW(GLOBAL, LOG_INFO, "Writing output for step %d.\n",simCtx->step);
430 for (PetscInt bi = 0; bi < simCtx->block_number; bi++) {
431 ierr = WriteSimulationFields(&user[bi]); CHKERRQ(ierr);
432 }
433 if (simCtx->np > 0) {
434 ierr = WriteAllSwarmFields(user); CHKERRQ(ierr);
435 if(get_log_level() >= LOG_INFO){
436 LOG(GLOBAL, LOG_INFO, "Particle states at step %d:\n", simCtx->step);
438 }
439 }
440 }
441
443
444 // Update Progress Bar
445 if(simCtx->rank == 0) {
446 PrintProgressBar(step,StartStep,StepsToRun,simCtx->ti);
447 if(get_log_level()>=LOG_WARNING) PetscPrintf(PETSC_COMM_SELF,"\n");
448 }
449 } // --- End of Time-Marching Loop ---
450
451 // After the loop, print the 100% complete bar on rank 0 and add a newline
452 // to ensure subsequent terminal output starts on a fresh line.
453 if (simCtx->rank == 0 && StepsToRun > 0) {
454 PrintProgressBar(StartStep + StepsToRun - 1, StartStep, StepsToRun, simCtx->ti);
455 PetscPrintf(PETSC_COMM_SELF, "\n");
456 fflush(stdout);
457 }
458
459 LOG_ALLOW(GLOBAL, LOG_INFO, "Time marching completed. Final time t=%.4f.\n", simCtx->ti);
461 PetscFunctionReturn(0);
462}
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 UpdateAllParticlePositions(UserCtx *user)
Loops over all local particles in the DMSwarm, updating their positions based on velocity and the glo...
PetscErrorCode LocateAllParticlesInGrid(UserCtx *user, BoundingBox *bboxlist)
Orchestrates the complete particle location and migration process for one timestep.
PetscErrorCode ResetAllParticleStatuses(UserCtx *user)
This function is designed to be called at the end of a full timestep, after all particle-based calcul...
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 UpdateAllParticleFields(UserCtx *user)
Orchestrates the update of all physical properties for particles.
PetscErrorCode InterpolateAllFieldsToSwarm(UserCtx *user)
Interpolates all relevant fields from the DMDA to the DMSwarm.
PetscErrorCode WriteAllSwarmFields(UserCtx *user)
Writes a predefined set of PETSc Swarm fields to files.
Definition io.c:1836
PetscErrorCode ReadSimulationFields(UserCtx *user, PetscInt ti)
Reads binary field data for velocity, pressure, and other required vectors.
Definition io.c:1024
PetscErrorCode WriteSimulationFields(UserCtx *user)
Writes simulation fields to files.
Definition io.c:1536
PetscBool is_function_allowed(const char *functionName)
Checks if a given function is in the allow-list.
Definition logging.c:157
#define LOCAL
Logging scope definitions for controlling message output.
Definition logging.h:46
#define GLOBAL
Scope for global logging across all processes.
Definition logging.h:47
#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:201
#define PROFILE_FUNCTION_END
Marks the end of a profiled code block.
Definition logging.h:740
#define LOG(scope, level, fmt,...)
Logging macro for PETSc-based applications with scope control.
Definition logging.h:85
PetscErrorCode LOG_FIELD_MIN_MAX(UserCtx *user, const char *fieldName)
Computes and logs the local and global min/max values of a 3-component vector field.
Definition logging.c:1273
void PrintProgressBar(PetscInt step, PetscInt startStep, PetscInt totalSteps, PetscReal currentTime)
Prints a progress bar to the console.
Definition logging.c:1207
PetscErrorCode LOG_FIELD_ANATOMY(UserCtx *user, const char *field_name, const char *stage_name)
Logs the anatomy of a specified field at key boundary locations, respecting the solver's specific gri...
Definition logging.c:1435
LogLevel get_log_level()
Retrieves the current logging level from the environment variable LOG_LEVEL.
Definition logging.c:39
PetscErrorCode LOG_PARTICLE_FIELDS(UserCtx *user, PetscInt printInterval)
Prints particle fields in a table that automatically adjusts its column widths.
Definition logging.c:400
@ LOG_INFO
Informational messages about program execution.
Definition logging.h:32
@ LOG_WARNING
Non-critical issues that warrant attention.
Definition logging.h:30
@ LOG_DEBUG
Detailed debugging information.
Definition logging.h:33
@ LOG_VERBOSE
Extremely detailed logs, typically for development use only.
Definition logging.h:35
#define PROFILE_FUNCTION_BEGIN
Marks the beginning of a profiled code block (typically a function).
Definition logging.h:731
PetscErrorCode ProfilingLogTimestepSummary(PetscInt step)
Logs the performance summary for the current timestep and resets timers.
Definition logging.c:1047
PetscErrorCode AdvanceSimulation(SimCtx *simCtx)
Executes the main time-marching loop for the coupled Euler-Lagrange simulation.
Definition simulation.c:307
PetscErrorCode UpdateSolverHistoryVectors(UserCtx *user)
Copies the current time step's solution fields into history vectors (e.g., U(t_n) -> U_o,...
Definition simulation.c:23
PetscErrorCode FinalizeRestartState(SimCtx *simCtx)
Performs post-load/post-init consistency checks for a restarted simulation.
Definition simulation.c:242
PetscErrorCode PerformInitializedParticleSetup(SimCtx *simCtx)
Finalizes the simulation setup at t=0, ensuring a consistent state before time marching.
Definition simulation.c:95
PetscErrorCode PerformLoadedParticleSetup(SimCtx *simCtx)
Finalizes the simulation state after particle and fluid data have been loaded from a restart.
Definition simulation.c:176
#define __FUNCT__
Definition simulation.c:79
PetscErrorCode FlowSolver(SimCtx *simCtx)
Orchestrates a single time step of the Eulerian fluid solver.
Definition solvers.c:26
UserCtx * user
Definition variables.h:427
PetscBool inletFaceDefined
Definition variables.h:680
PetscMPIInt rank
Definition variables.h:541
char eulerianSource[64]
Definition variables.h:557
PetscInt block_number
Definition variables.h:593
SimCtx * simCtx
Back-pointer to the master simulation context.
Definition variables.h:664
PetscInt rans
Definition variables.h:609
PetscReal StartTime
Definition variables.h:551
Vec K_Omega_o
Definition variables.h:712
UserMG usermg
Definition variables.h:631
Vec K_Omega
Definition variables.h:712
Vec lUcont_rm1
Definition variables.h:692
PetscInt _this
Definition variables.h:674
PetscReal dt
Definition variables.h:552
PetscInt StepsToRun
Definition variables.h:549
PetscInt np
Definition variables.h:616
Vec Ucont
Definition variables.h:688
PetscInt StartStep
Definition variables.h:548
Vec lUcont_o
Definition variables.h:691
Vec Ucat_o
Definition variables.h:691
char particleRestartMode[16]
Definition variables.h:621
BoundingBox * bboxlist
Definition variables.h:619
PetscInt OutputFreq
Definition variables.h:555
Vec lK_Omega_o
Definition variables.h:712
Vec Ucat
Definition variables.h:688
Vec Ucont_o
Definition variables.h:691
PetscInt mglevels
Definition variables.h:434
Vec Nvert_o
Definition variables.h:691
Vec Ucont_rm1
Definition variables.h:692
PetscInt step
Definition variables.h:546
Vec Nvert
Definition variables.h:688
MGCtx * mgctx
Definition variables.h:437
Vec lNvert_o
Definition variables.h:691
PetscInt ParticleInitialization
Definition variables.h:620
PetscReal ti
Definition variables.h:547
PetscInt immersed
Definition variables.h:567
PetscInt LoggingFrequency
Definition variables.h:636
Vec P_o
Definition variables.h:691
Defines a 3D axis-aligned bounding box.
Definition variables.h:154
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