PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
Loading...
Searching...
No Matches
solvers.h File Reference
#include "variables.h"
#include "implicitsolvers.h"
#include "rhs.h"
#include "logging.h"
#include "poisson.h"
#include "setup.h"
#include "les.h"
Include dependency graph for solvers.h:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Functions

PetscErrorCode FlowSolver (SimCtx *simCtx)
 Orchestrates a single time step of the Eulerian fluid solver.
 

Function Documentation

◆ FlowSolver()

PetscErrorCode FlowSolver ( SimCtx simCtx)

Orchestrates a single time step of the Eulerian fluid solver.

This is the refactored, high-level entry point for advancing the fluid state from time t_n to t_{n+1}. It takes the master SimCtx as its primary argument.

Parameters
simCtxThe master simulation context, containing all solver settings, multigrid structures, and data vectors.
Returns
PetscErrorCode 0 on success.

This function is the refactored, high-level entry point for advancing the fluid state from time t_n to t_{n+1}. It is a direct adaptation of the legacy FlowSolver, but it now takes the master SimCtx as its primary argument to eliminate dependencies on global variables.

The sequence of operations is:

  1. (Optional) Update turbulence models (RANS/LES) to compute eddy viscosity.
  2. Call the core momentum solver (explicit Runge-Kutta or an implicit scheme) to get an intermediate velocity field.
  3. Call the pressure-Poisson solver to compute the pressure correction.
  4. Call the projection step to correct the velocity field, ensuring it is divergence-free.
  5. Perform final state updates, diagnostics, and I/O.
Parameters
simCtxThe master simulation context, containing all solver settings, multigrid structures, and data vectors.
Returns
PetscErrorCode 0 on success.

Definition at line 26 of file solvers.c.

27{
28 PetscErrorCode ierr;
29 UserMG *usermg = &simCtx->usermg;
30 PetscInt level = usermg->mglevels - 1;
31 UserCtx *user = usermg->mgctx[level].user;
32
33 PetscFunctionBeginUser;
35 LOG_ALLOW(GLOBAL, LOG_INFO, "[Step %d] Entering orchestrator...\n", simCtx->step);
36
37 /*
38 // ========================================================================
39 // SECTION: O-Grid Specific Force Calculations (Legacy Feature)
40 // ========================================================================
41 // This was a specialized calculation for non-immersed O-grid cases.
42 if (simCtx->Ogrid && !simCtx->immersed) {
43 LOG_ALLOW(GLOBAL, LOG_DEBUG, "Calculating O-grid forces...\n");
44 // Calc_forces_Ogrid(&user[0], simCtx->step, 0);
45 }
46 */
47
48
49
50 // ========================================================================
51 // SECTION: Turbulence Models (RANS/LES)
52 // ========================================================================
53 // These models compute the turbulent eddy viscosity (Nu_t) which is then
54 // used by the momentum solver in the diffusion term.
55
56 /*
57 if (simCtx->rans) {
58 LOG_ALLOW(GLOBAL, LOG_INFO, "Updating RANS (k-omega) model...\n");
59 for (PetscInt bi = 0; bi < simCtx->block_number; bi++) {
60 K_Omega_Set_Constant(&user[bi]);
61 if (simCtx->step == simCtx->StartStep) {
62 LOG_ALLOW(LOCAL, LOG_DEBUG, " Initializing K-Omega field for block %d.\n", bi);
63 K_Omega_IC(&user[bi]);
64 }
65 // In a full implementation, the K-Omega transport equations would be solved here.
66 // Solve_K_Omega(&user[bi]);
67 }
68 }
69 */
70
71 if (simCtx->les) {
72 LOG_ALLOW(GLOBAL, LOG_INFO, "Updating LES (Smagorinsky) model...\n");
73 for (PetscInt bi = 0; bi < simCtx->block_number; bi++) {
74 // LES models require Cartesian velocity to compute strain rates
75 UpdateLocalGhosts(&user[bi], "Ucont");
76 ierr = Contra2Cart(&user[bi]); CHKERRQ(ierr);
77 UpdateLocalGhosts(&user[bi], "Ucat");
78 if (simCtx->step % simCtx->dynamic_freq == 0) {
79 LOG_ALLOW(LOCAL, LOG_DEBUG, " Computing Smagorinsky constant for block %d.\n", bi);
81 }
82 // LOG_ALLOW(LOCAL, LOG_DEBUG, " Computing eddy viscosity for block %d.\n", bi);
83 ComputeEddyViscosityLES(&user[bi]);
84 }
85 }
86
87
88 // ========================================================================
89 // SECTION: Momentum Equation Solver
90 // ========================================================================
91 // This is the core of the time step. It computes an intermediate velocity
92 // field by solving the momentum equations.
93
94 LOG_ALLOW(GLOBAL, LOG_INFO, "Beginning momentum solve (Implicit Flag = %d)...\n", simCtx->implicit);
95
96 if (simCtx->implicit > 0) {
97 // --- Implicit Path ---
98 LOG_ALLOW(GLOBAL, LOG_INFO, "Executing implicit momentum solver...\n");
99 // Since IBM is disabled, we pass NULL for ibm and fsi arguments.
100 ierr = ImpRK(user, NULL, NULL); CHKERRQ(ierr);
101 } else {
102 // --- Explicit Path (Default) ---
103 LOG_ALLOW(GLOBAL, LOG_INFO, "Executing Iterative implicit momentum solver (Implicit Runge-Kutta)...\n");
104 // Since IBM is disabled, we pass NULL for ibm and fsi arguments.
105 ierr = RungeKutta(user, NULL, NULL); CHKERRQ(ierr);
106 }
107
108// ========================================================================
109// SECTION: Pressure-Poisson Solver
110// ========================================================================
111// This step enforces the continuity equation (incompressibility) by solving
112// for a pressure correction field.
113
114 LOG_ALLOW(GLOBAL, LOG_INFO, "Beginning pressure-Poisson solve (Poisson Flag = %d)...\n", simCtx->poisson);
115
116 if (simCtx->poisson == 0) {
117 ierr = PoissonSolver_MG(usermg); CHKERRQ(ierr);
118 } else {
119 // Logic for other Poisson solvers (e.g., Hypre) would go here.
120 LOG_ALLOW(GLOBAL, LOG_WARNING, "Poisson solver type %d is not currently enabled.\n", simCtx->poisson);
121 }
122
123 // ========================================================================
124 // SECTION: Velocity Correction (Projection)
125 // ========================================================================
126 // The pressure correction is used to update the pressure field and project
127 // the intermediate velocity onto a divergence-free space.
128
129 LOG_ALLOW(GLOBAL, LOG_INFO, "Applying velocity correction/projection step...\n");
130 for (PetscInt bi = 0; bi < simCtx->block_number; bi++) {
131 ierr = UpdatePressure(&user[bi]); CHKERRQ(ierr);
132 LOG_ALLOW(GLOBAL,LOG_INFO," Pressure Updated for Block %d.\n",bi);
133
134 ierr = Projection(&user[bi]); CHKERRQ(ierr);
135
136 LOG_ALLOW(GLOBAL,LOG_INFO," Velocity corrected for Block %d.\n",bi);
137
138 // Ensure local ghost cells for the final pressure field are correct
139 ierr = UpdateLocalGhosts(&user[bi],"P");
140 }
141
142 // ========================================================================
143
144 // ========================================================================
145 // SECTION: Final Diagnostics and I/O
146 // ========================================================================
147
148 for (PetscInt bi = 0; bi < simCtx->block_number; bi++) {
149 LOG_ALLOW(GLOBAL, LOG_INFO, "Finalizing state & Diagnostics for block %d...\n", bi);
150
151 // --- Perform Divergence Check ---
152 // This is a diagnostic to verify the quality of the velocity correction.
153 ierr = ComputeDivergence(&user[bi]); CHKERRQ(ierr);
154
155 // -- Log Continuity metrics ----
156 ierr = LOG_CONTINUITY_METRICS(&user[bi]);
157 /*
158 // --- Immersed Boundary Interpolation (Post-Correction) ---
159 // This step would update the velocity values AT the IB nodes to match the
160 // newly corrected fluid field. Important for the next time step.
161 if (simCtx->immersed) {
162 for (PetscInt ibi = 0; ibi < simCtx->NumberOfBodies; ibi++) {
163 ibm_interpolation_advanced(&user[bi], &simCtx->ibm[ibi], ibi, 1);
164 }
165 }
166 */
167
168 // --- Averaging and Statistics (if enabled) ---
169 /*
170 if (simCtx->averaging) {
171 LOG_ALLOW(LOCAL, LOG_DEBUG, "Performing statistical averaging for block %d.\n", bi);
172 Do_averaging(&user[bi]);
173 }
174 */
175
176 // }
177 }
178
179 LOG_ALLOW(GLOBAL, LOG_INFO, "orchestrator finished for step %d.\n", simCtx->step);
181 PetscFunctionReturn(0);
182}
PetscErrorCode ImpRK(UserCtx *user, IBMNodes *ibm, FSInfo *fsi)
Advances the momentum equations using an iterative explicit Runge-Kutta scheme.
PetscErrorCode RungeKutta(UserCtx *user, IBMNodes *ibm, FSInfo *fsi)
Advances the momentum equations using an explicit 4th-order Runge-Kutta scheme.
PetscErrorCode ComputeEddyViscosityLES(UserCtx *user)
Computes the turbulent eddy viscosity (Nu_t) for the LES model.
Definition les.c:301
PetscErrorCode ComputeSmagorinskyConstant(UserCtx *user)
Computes the dynamic Smagorinsky constant (Cs) for the LES model.
Definition les.c:23
#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
PetscErrorCode LOG_CONTINUITY_METRICS(UserCtx *user)
Logs continuity metrics for a single block to a file.
Definition logging.c:862
@ 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
#define PROFILE_FUNCTION_BEGIN
Marks the beginning of a profiled code block (typically a function).
Definition logging.h:731
PetscErrorCode Projection(UserCtx *user)
Corrects the contravariant velocity field Ucont to be divergence-free using the gradient of the press...
Definition poisson.c:948
PetscErrorCode UpdatePressure(UserCtx *user)
Updates the pressure field P with the pressure correction Phi computed by the Poisson solver.
Definition poisson.c:1502
PetscErrorCode PoissonSolver_MG(UserMG *usermg)
Solves the pressure-Poisson equation using a geometric multigrid method.
Definition poisson.c:3959
PetscErrorCode Contra2Cart(UserCtx *user)
Reconstructs Cartesian velocity (Ucat) at cell centers from contravariant velocity (Ucont) defined on...
Definition setup.c:2084
PetscErrorCode ComputeDivergence(UserCtx *user)
Definition setup.c:2377
PetscErrorCode UpdateLocalGhosts(UserCtx *user, const char *fieldName)
Updates the local vector (including ghost points) from its corresponding global vector.
Definition setup.c:1091
UserCtx * user
Definition variables.h:427
PetscInt block_number
Definition variables.h:593
PetscInt implicit
Definition variables.h:575
UserMG usermg
Definition variables.h:631
PetscInt poisson
Definition variables.h:578
PetscInt mglevels
Definition variables.h:434
PetscInt dynamic_freq
Definition variables.h:610
PetscInt step
Definition variables.h:546
PetscInt les
Definition variables.h:609
MGCtx * mgctx
Definition variables.h:437
User-defined context containing data specific to a single computational grid level.
Definition variables.h:661
User-level context for managing the entire multigrid hierarchy.
Definition variables.h:433
Here is the call graph for this function:
Here is the caller graph for this function: