PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
Loading...
Searching...
No Matches
simulation.h File Reference
#include <petscpf.h>
#include <petscdmswarm.h>
#include <stdlib.h>
#include <time.h>
#include <math.h>
#include <petsctime.h>
#include <petscsys.h>
#include <petscdmcomposite.h>
#include <petscsystypes.h>
#include "variables.h"
#include "ParticleSwarm.h"
#include "walkingsearch.h"
#include "grid.h"
#include "logging.h"
#include "io.h"
#include "interpolation.h"
#include "initialcondition.h"
#include "ParticleMotion.h"
#include "Boundaries.h"
#include "setup.h"
#include "solvers.h"
Include dependency graph for simulation.h:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Functions

PetscErrorCode UpdateSolverHistoryVectors (UserCtx *user)
 Copies the current time step's solution fields into history vectors (e.g., U(t_n) -> U_o, U_o -> U_rm1) for the next time step's calculations.
 
PetscErrorCode AdvanceSimulation (SimCtx *simCtx)
 Executes the main time-marching loop for the particle simulation.
 
PetscErrorCode PerformInitialSetup (SimCtx *simCtx)
 Finalizes the simulation setup at t=0, ensuring a consistent state before time marching.
 

Function Documentation

◆ UpdateSolverHistoryVectors()

PetscErrorCode UpdateSolverHistoryVectors ( UserCtx user)

Copies the current time step's solution fields into history vectors (e.g., U(t_n) -> U_o, U_o -> U_rm1) for the next time step's calculations.

This function is critical for multi-step time integration schemes (like BDF2) used by the legacy solver. It must be called at the end of every time step, after the new solution has been fully computed.

The order of operations is important to avoid overwriting data prematurely.

Parameters
userThe UserCtx for a single block. The function modifies the history vectors (Ucont_o, Ucont_rm1, etc.) within this context.
Returns
PetscErrorCode 0 on success.

Definition at line 23 of file simulation.c.

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}
#define LOCAL
Logging scope definitions for controlling message output.
Definition logging.h:44
#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:207
@ LOG_DEBUG
Detailed debugging information.
Definition logging.h:33
PetscMPIInt rank
Definition variables.h:516
SimCtx * simCtx
Back-pointer to the master simulation context.
Definition variables.h:633
PetscInt rans
Definition variables.h:578
Vec K_Omega_o
Definition variables.h:680
Vec K_Omega
Definition variables.h:680
Vec lUcont_rm1
Definition variables.h:661
PetscInt _this
Definition variables.h:643
Vec Ucont
Definition variables.h:657
Vec lUcont_o
Definition variables.h:660
Vec Ucat_o
Definition variables.h:660
Vec lK_Omega_o
Definition variables.h:680
Vec Ucat
Definition variables.h:657
Vec Ucont_o
Definition variables.h:660
Vec Nvert_o
Definition variables.h:660
Vec Ucont_rm1
Definition variables.h:661
Vec Nvert
Definition variables.h:657
Vec lNvert_o
Definition variables.h:660
PetscInt immersed
Definition variables.h:535
Vec P_o
Definition variables.h:660
The master context for the entire simulation.
Definition variables.h:513
Here is the caller graph for this function:

◆ AdvanceSimulation()

PetscErrorCode AdvanceSimulation ( SimCtx simCtx)

Executes the main time-marching loop for the particle simulation.

[TEST VERSION]

This version uses the new, integrated LocateAllParticlesInGrid_TEST orchestrator and the ResetAllParticleStatuses helper for a clean, robust, and understandable workflow.

For each timestep, it performs:

  1. Sets the background fluid velocity field (Ucat) for the current step.
  2. Updates particle positions using velocity from the previous step's interpolation.
  3. Removes any particles that have left the global domain.
  4. A single call to LocateAllParticlesInGrid_TEST, which handles all particle location and migration until the swarm is fully settled.
  5. Interpolates the current fluid velocity to the newly settled particle locations.
  6. Scatters particle data back to Eulerian fields.
  7. Outputs data at specified intervals.
Parameters
userPointer to the UserCtx structure..
Returns
PetscErrorCode 0 on success, non-zero on failure.

Executes the main time-marching loop for the particle simulation.

This function orchestrates the advancement of the simulation from the configured StartStep to the final step. It does NOT perform the initial t=0 setup, as that is handled by InitializeEulerianState and PerformInitialSetup in main().

For each timestep, it performs the following sequence:

  1. Pre-Solver Actions: Updates time-dependent boundary conditions (e.g., fluxin) and resets particle states for the new step.
  2. Eulerian Solve: Calls the refactored legacy Flow_Solver to advance the entire fluid field by one time step.
  3. Lagrangian Update: Executes the full particle workflow: advection based on the previous step's velocity, followed by settling (location/migration) in the new grid, and finally interpolation of the new fluid velocity.
  4. Two-Way Coupling: Scatters particle data back to the grid to act as source terms for the subsequent time step.
  5. History & I/O: Updates the solver's history vectors and writes output files at the specified frequency.
Parameters
userArray of UserCtx structures for the finest grid level.
Returns
PetscErrorCode 0 on success, non-zero on failure.

Definition at line 181 of file simulation.c.

182{
183 PetscErrorCode ierr;
184 // Get the master context from the first block. All blocks share it.
185 UserCtx *user = simCtx->usermg.mgctx[simCtx->usermg.mglevels-1].user;
186
187 // Retrieve control parameters from SimCtx for clarity
188 const PetscInt StartStep = simCtx->StartStep;
189 const PetscInt StepsToRun = simCtx->StepsToRun;
190 const PetscInt OutputFreq = simCtx->OutputFreq;
191 const PetscReal dt = simCtx->dt;
192
193 // Variables for particle removal statistics
194 PetscInt removed_local_ob, removed_global_ob;
195 PetscInt removed_local_lost, removed_global_lost;
196
197 PetscFunctionBeginUser;
199 LOG_ALLOW(GLOBAL, LOG_INFO, "Starting main time-marching loop: %d steps from step %d (t=%.4f), dt=%.4f\n",
200 StepsToRun, StartStep, simCtx->StartTime, dt);
201
202 // --- Main Time-Marching Loop ---
203 for (PetscInt step = StartStep; step < StartStep + StepsToRun; ++step) {
204
205 // =================================================================
206 // 1. PRE-STEP SETUP
207 // =================================================================
208
209 // Update simulation time and step counters in the master context
210 simCtx->step = step;
211 simCtx->ti = simCtx->StartTime + step * simCtx->dt;
212
213 LOG_ALLOW(GLOBAL, LOG_INFO, "--- Advancing Step %d (t=%.4f) ---\n", step, simCtx->ti);
214
215 // Update any time-dependent boundary conditions (e.g., pulsating inlet)
216 // if (simCtx->inletprofile == 3) {
217 // LOG_ALLOW(GLOBAL, LOG_DEBUG, "Updating time-dependent inlet flux.\n");
218 // fluxin(&user[0]); // Assumes block 0 drives the global flux value
219 // }
220
221 // For particles, reset their status to prepare for the new advection/location cycle
222 if (simCtx->np > 0) {
223 LOG_ALLOW(GLOBAL, LOG_DEBUG, "Resetting all particle statuses to NEEDS_LOCATION.\n");
224 ierr = ResetAllParticleStatuses(user); CHKERRQ(ierr);
225 }
226
227 // =================================================================
228 // 2. EULERIAN SOLVER STEP
229 // =================================================================
230
231 // Call the refactored, high-level legacy solver. This single function
232 // advances the entire multi-block fluid field from t_n to t_{n+1}.
233 LOG_ALLOW(GLOBAL, LOG_INFO, "Updating Eulerian Field ...\n");
234
235 ierr = Flow_Solver(simCtx); CHKERRQ(ierr);
236
237 LOG_ALLOW(GLOBAL, LOG_INFO, "Eulerian Field solved ...\n");
238 // =================================================================
239 // 3. LAGRANGIAN PARTICLE STEP
240 // =================================================================
241
242 if (simCtx->np > 0) {
243 LOG_ALLOW(GLOBAL, LOG_INFO, "Updating Lagrangian particle system...\n");
244
245 // a. Advect particles using the velocity interpolated from the *previous* step.
246 // P(t_{n+1}) = P(t_n) + V_p(t_n) * dt
247 ierr = UpdateAllParticlePositions(user); CHKERRQ(ierr);
248
249 // b. Settle all particles: find their new host cells and migrate them across ranks.
250 ierr = LocateAllParticlesInGrid_TEST(user, simCtx->bboxlist); CHKERRQ(ierr);
251
252 // c. Remove any particles that are now lost or out of the global domain.
253 ierr = CheckAndRemoveLostParticles(user, &removed_local_lost, &removed_global_lost); CHKERRQ(ierr);
254 ierr = CheckAndRemoveOutOfBoundsParticles(user, &removed_local_ob, &removed_global_ob, simCtx->bboxlist); CHKERRQ(ierr);
255 if (removed_global_lost + removed_global_ob > 0) {
256 LOG_ALLOW(GLOBAL, LOG_INFO, "Removed %d particles globally this step.\n", removed_global_lost + removed_global_ob);
257 }
258
259 // d. Interpolate the NEW fluid velocity (just computed by Flow_Solver) onto the
260 // particles' new positions. This gives them V_p(t_{n+1}) for the *next* advection step.
261 ierr = InterpolateAllFieldsToSwarm(user); CHKERRQ(ierr);
262
263 // e. (For Two-Way Coupling) Scatter particle data back to the grid to act as a source term.
264 ierr = CalculateParticleCountPerCell(user); CHKERRQ(ierr);
265 ierr = ScatterAllParticleFieldsToEulerFields(user); CHKERRQ(ierr);
266 }
267
268 // =================================================================
269 // 4. UPDATE HISTORY & I/O
270 // =================================================================
271
272 // Copy the newly computed fields (Ucont, P, etc.) to the history vectors
273 // (_o, _rm1) to prepare for the next time step.
274 for (PetscInt bi = 0; bi < simCtx->block_number; bi++) {
275 ierr = UpdateSolverHistoryVectors(&user[bi]); CHKERRQ(ierr);
276 }
277
278 // Handle periodic file output
279 if (OutputFreq > 0 && (step + 1) % OutputFreq == 0) {
280 LOG_ALLOW(GLOBAL, LOG_INFO, "Writing output for step %d.\n", step + 1);
281 for (PetscInt bi = 0; bi < simCtx->block_number; bi++) {
282 ierr = WriteSimulationFields(&user[bi]); CHKERRQ(ierr);
283 }
284 if (simCtx->np > 0) {
285 ierr = WriteAllSwarmFields(user); CHKERRQ(ierr);
286 if(get_log_level() >= LOG_INFO){
288 }
289 }
290 }
291 } // --- End of Time-Marching Loop ---
292
293 LOG_ALLOW(GLOBAL, LOG_INFO, "Time marching completed. Final time t=%.4f.\n", simCtx->ti + dt);
295 PetscFunctionReturn(0);
296}
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 LocateAllParticlesInGrid_TEST(UserCtx *user, BoundingBox *bboxlist)
Orchestrates the complete particle location and migration process for one timestep.
PetscErrorCode UpdateAllParticlePositions(UserCtx *user)
Loops over all local particles in the DMSwarm, updating their positions based on velocity and the glo...
PetscErrorCode ResetAllParticleStatuses(UserCtx *user)
This function is designed to be called at the end of a full timestep, after all particle-based calcul...
PetscErrorCode CheckAndRemoveLostParticles(UserCtx *user, PetscInt *removedCountLocal, PetscInt *removedCountGlobal)
Removes particles that have been definitively flagged as LOST by the location algorithm.
PetscErrorCode CalculateParticleCountPerCell(UserCtx *user)
Counts particles in each cell of the DMDA 'da' and stores the result in user->ParticleCount.
PetscErrorCode ScatterAllParticleFieldsToEulerFields(UserCtx *user)
Scatters a predefined set of particle fields to their corresponding Eulerian fields.
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:1568
PetscErrorCode WriteSimulationFields(UserCtx *user)
Writes simulation fields to files.
Definition io.c:1280
#define GLOBAL
Scope for global logging across all processes.
Definition logging.h:45
#define PROFILE_FUNCTION_END
Marks the end of a profiled code block.
Definition logging.h:767
LogLevel get_log_level()
Retrieves the current logging level from the environment variable LOG_LEVEL.
Definition logging.c:49
PetscErrorCode LOG_PARTICLE_FIELDS(UserCtx *user, PetscInt printInterval)
Prints particle fields in a table that automatically adjusts its column widths.
Definition logging.c:402
@ LOG_INFO
Informational messages about program execution.
Definition logging.h:32
#define PROFILE_FUNCTION_BEGIN
Marks the beginning of a profiled code block (typically a function).
Definition logging.h:758
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 Flow_Solver(SimCtx *simCtx)
Orchestrates a single time step of the Eulerian fluid solver.
Definition solvers.c:26
UserCtx * user
Definition variables.h:418
PetscInt block_number
Definition variables.h:562
PetscReal StartTime
Definition variables.h:526
UserMG usermg
Definition variables.h:599
PetscReal dt
Definition variables.h:527
PetscInt StepsToRun
Definition variables.h:524
PetscInt np
Definition variables.h:585
PetscInt StartStep
Definition variables.h:523
BoundingBox * bboxlist
Definition variables.h:588
PetscInt OutputFreq
Definition variables.h:531
PetscInt mglevels
Definition variables.h:425
PetscInt step
Definition variables.h:521
MGCtx * mgctx
Definition variables.h:428
PetscReal ti
Definition variables.h:522
PetscInt LoggingFrequency
Definition variables.h:604
User-defined context containing data specific to a single computational grid level.
Definition variables.h:630
Here is the call graph for this function:
Here is the caller graph for this function:

◆ PerformInitialSetup()

PetscErrorCode PerformInitialSetup ( SimCtx simCtx)

Finalizes the simulation setup at t=0, ensuring a consistent state before time marching.

This function is called from main() after the initial Eulerian and Lagrangian states have been created but before the main time loop begins. Its responsibilities are:

  1. Settling the particle swarm: Migrates particles to their correct owner ranks and finds their initial host cells. This includes handling special surface initializations.
  2. Coupling the fields: Interpolates the initial Eulerian fields to the settled particle locations.
  3. Preparing for the first step: Scatters particle data back to the grid.
  4. Writing the initial output for step 0.
Parameters
simCtxPointer to the main simulation context structure.
Returns
PetscErrorCode 0 on success, non-zero on failure.

Definition at line 95 of file simulation.c.

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 // --- Set simulation time and step for this specific phase ---
105 simCtx->step = 0;
106 simCtx->ti = simCtx->StartTime;
107
108 LOG_ALLOW(GLOBAL, LOG_INFO, "[T=%.4f, Step=%d] Performing initial particle setup procedures.\n", simCtx->ti, simCtx->step);
109
110 // --- 1. Initial Particle Settlement (Location and Migration) ---
111 LOG_ALLOW(GLOBAL, LOG_INFO, "[T=%.4f, Step=%d] Initial Settlement: Locating and migrating all particles...\n", simCtx->ti, simCtx->step);
112 ierr = LocateAllParticlesInGrid_TEST(user, bboxlist); CHKERRQ(ierr);
113
114 // --- 2. Re-initialize Particles on Inlet Surface (if applicable) ---
115 // Note: Use simCtx->ParticleInitialization for consistency
116 if (simCtx->ParticleInitialization == 0 && 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_TEST(user, bboxlist); CHKERRQ(ierr);
125 }
126
127 // --- 3. Finalize State for t=0 ---
128 LOG_ALLOW(GLOBAL, LOG_INFO, "[T=%.4f, Step=%d] Interpolating initial fields to settled particles.\n", simCtx->ti, simCtx->step);
129 ierr = InterpolateAllFieldsToSwarm(user); CHKERRQ(ierr);
130 ierr = ScatterAllParticleFieldsToEulerFields(user); CHKERRQ(ierr);
131
132 // --- 4. Initial History and Output ---
133 // Update solver history vectors with the t=0 state before the first real step
134 for (PetscInt bi = 0; bi < simCtx->block_number; bi++) {
135 ierr = UpdateSolverHistoryVectors(&user[bi]); CHKERRQ(ierr);
136 }
137
138 if (simCtx->OutputFreq > 0 || (simCtx->StepsToRun == 0 && simCtx->StartStep == 0)) {
139 LOG_ALLOW(GLOBAL, LOG_INFO, "[T=%.4f, Step=%d] Writing initial simulation data.\n", simCtx->ti, simCtx->step);
140
141 // --- Particle Output (assumes functions operate on the master user context) ---
142 ierr = LOG_PARTICLE_FIELDS(user, simCtx->LoggingFrequency); CHKERRQ(ierr);
143 ierr = WriteAllSwarmFields(user); CHKERRQ(ierr);
144
145 // --- Eulerian Field Output (MUST loop over all blocks) --- // <<< CHANGED/FIXED
146 for (PetscInt bi = 0; bi < simCtx->block_number; bi++) {
147 ierr = WriteSimulationFields(&user[bi]); CHKERRQ(ierr);
148 }
149 }
150
151 LOG_ALLOW(GLOBAL, LOG_INFO, "--- Initial setup complete. Ready for time marching. ---\n\n");
152 PetscFunctionReturn(0);
153}
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...
PetscBool inletFaceDefined
Definition variables.h:649
PetscInt ParticleInitialization
Definition variables.h:589
Defines a 3D axis-aligned bounding box.
Definition variables.h:153
Here is the call graph for this function:
Here is the caller graph for this function: