PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
Loading...
Searching...
No Matches
Macros | Functions
momentumsolvers.c File Reference
#include "momentumsolvers.h"
Include dependency graph for momentumsolvers.c:

Go to the source code of this file.

Macros

#define __FUNCT__   "ComputeTotalResidual"
 
#define __FUNCT__   "Momentum_Solver_Explicit_RungeKutta4"
 
#define __FUNCT__   "MomentumSolver_DualTime_Picard_RK4"
 

Functions

static PetscErrorCode ComputeTotalResidual (UserCtx *user)
 Internal helper implementation: ComputeTotalResidual().
 
PetscErrorCode MomentumSolver_Explicit_RungeKutta4 (UserCtx *user, IBMNodes *ibm, FSInfo *fsi)
 Internal helper implementation: MomentumSolver_Explicit_RungeKutta4().
 
PetscErrorCode MomentumSolver_DualTime_Picard_RK4 (UserCtx *user, IBMNodes *ibm, FSInfo *fsi)
 Internal helper implementation: MomentumSolver_DualTime_Picard_RK4().
 

Macro Definition Documentation

◆ __FUNCT__ [1/3]

#define __FUNCT__   "ComputeTotalResidual"

Definition at line 5 of file momentumsolvers.c.

◆ __FUNCT__ [2/3]

#define __FUNCT__   "Momentum_Solver_Explicit_RungeKutta4"

Definition at line 5 of file momentumsolvers.c.

◆ __FUNCT__ [3/3]

#define __FUNCT__   "MomentumSolver_DualTime_Picard_RK4"

Definition at line 5 of file momentumsolvers.c.

Function Documentation

◆ ComputeTotalResidual()

static PetscErrorCode ComputeTotalResidual ( UserCtx user)
static

Internal helper implementation: ComputeTotalResidual().

Local to this translation unit.

Definition at line 10 of file momentumsolvers.c.

11{
12 PetscErrorCode ierr;
13 SimCtx *simCtx = user->simCtx;
14
15 // Extract Time Parameters from Context
16 const PetscInt ti = simCtx->step;
17 const PetscInt tistart = simCtx->StartStep;
18 const PetscReal dt = simCtx->dt;
19
20 PetscFunctionBeginUser;
22 // 1. Calculate Spatial Terms (stored in user->Rhs)
23 // Rhs = -Div(Flux) + Viscous + Source
24 ierr = ComputeRHS(user, user->Rhs); CHKERRQ(ierr);
25
26 // 2. Add Physical Time Derivative Terms (BDF Discretization)
27 // The equation solved is: dU/dtau = RHS_Spatial + RHS_Temporal
28
29 if (COEF_TIME_ACCURACY > 1.1 && ti != tistart && ti != 1) {
30 // --- BDF2 (Second Order Backward Difference) ---
31 // Formula: (1.5*U^{n} - 2.0*U^{n-1} + 0.5*U^{n-2}) / dt
32 // Moved to RHS: -1.5/dt * U^{n} + ...
33 ierr = VecAXPY(user->Rhs, -COEF_TIME_ACCURACY/dt, user->Ucont);
34 ierr = VecAXPY(user->Rhs, +2.0/dt, user->Ucont_o);
35 ierr = VecAXPY(user->Rhs, -0.5/dt, user->Ucont_rm1);
36 } else {
37 // --- BDF1 (First Order / Euler Implicit) ---
38 // Formula: (U^{n} - U^{n-1}) / dt
39 // Moved to RHS: -1.0/dt * U^{n} + ...
40 ierr = VecAXPY(user->Rhs, -1.0/dt, user->Ucont);
41 ierr = VecAXPY(user->Rhs, +1.0/dt, user->Ucont_o);
42 }
43
44 // 3. Enforce Boundary Conditions on the Residual
45 ierr = EnforceRHSBoundaryConditions(user); CHKERRQ(ierr);
46
48 PetscFunctionReturn(0);
49}
PetscErrorCode EnforceRHSBoundaryConditions(UserCtx *user)
Enforces boundary conditions on the momentum equation's Right-Hand-Side (RHS) vector.
Definition Boundaries.c:591
#define PROFILE_FUNCTION_END
Marks the end of a profiled code block.
Definition logging.h:779
#define PROFILE_FUNCTION_BEGIN
Marks the beginning of a profiled code block (typically a function).
Definition logging.h:770
PetscErrorCode ComputeRHS(UserCtx *user, Vec Rhs)
Computes the Right-Hand Side (RHS) of the momentum equations.
Definition rhs.c:1185
Vec Rhs
Definition variables.h:845
SimCtx * simCtx
Back-pointer to the master simulation context.
Definition variables.h:814
PetscReal dt
Definition variables.h:658
Vec Ucont
Definition variables.h:837
PetscInt StartStep
Definition variables.h:653
Vec Ucont_o
Definition variables.h:844
Vec Ucont_rm1
Definition variables.h:845
PetscInt step
Definition variables.h:651
#define COEF_TIME_ACCURACY
Coefficient controlling the temporal accuracy scheme (e.g., 1.5 for 2nd Order Backward Difference).
Definition variables.h:57
The master context for the entire simulation.
Definition variables.h:643
Here is the call graph for this function:
Here is the caller graph for this function:

◆ MomentumSolver_Explicit_RungeKutta4()

PetscErrorCode MomentumSolver_Explicit_RungeKutta4 ( UserCtx user,
IBMNodes ibm,
FSInfo fsi 
)

Internal helper implementation: MomentumSolver_Explicit_RungeKutta4().

Advances the momentum equations using an explicit 4th-order Runge-Kutta scheme.

Local to this translation unit.

Definition at line 58 of file momentumsolvers.c.

59{
60 PetscErrorCode ierr;
61 (void)fsi;
62 // --- Context Acquisition ---
63 // Get the master simulation context from the first block's UserCtx.
64 // This is the bridge to access all former global variables.
65 SimCtx *simCtx = user[0].simCtx;
66 PetscReal dt = simCtx->dt;
67 //PetscReal st = simCtx->st;
68 PetscInt istage;
69 PetscReal alfa[] = {0.25, 1.0/3.0, 0.5, 1.0};
70
71 PetscFunctionBeginUser;
72
74
75 LOG_ALLOW(GLOBAL, LOG_INFO, "Executing explicit momentum solver (Runge-Kutta) for %d block(s).\n",simCtx->block_number);
76
77 // --- 1. Pre-Loop Initialization (Legacy Logic) ---
78 // This block prepares boundary conditions and allocates the RHS vector for all blocks
79 // before the main RK loop begins. This logic is preserved from the original code.
80 LOG_ALLOW(GLOBAL, LOG_DEBUG, "Preparing all blocks for Runge-Kutta solve...\n");
81 for (PetscInt bi = 0; bi < simCtx->block_number; bi++) {
82
83 ierr = ApplyBoundaryConditions(&user[bi]); CHKERRQ(ierr);
84
85 /*
86 // Immersed boundary interpolation (if enabled)
87 if (simCtx->immersed) {
88 LOG_ALLOW(LOCAL, LOG_DEBUG, " Performing pre-RK IBM interpolation for block %d.\n", bi);
89 for (PetscInt ibi = 0; ibi < simCtx->NumberOfBodies; ibi++) {
90 // The 'ibm' and 'fsi' pointers are passed directly from FlowSolver
91 ierr = ibm_interpolation_advanced(&user[bi], &ibm[ibi], ibi, 1); CHKERRQ(ierr);
92 }
93 }
94 */
95
96 // Allocate the persistent RHS vector for this block's context
97 ierr = VecDuplicate(user[bi].Ucont, &user[bi].Rhs); CHKERRQ(ierr);
98 }
99 LOG_ALLOW(GLOBAL, LOG_DEBUG, "Pre-loop initialization complete.\n");
100
101 // --- 2. Main Runge-Kutta Loop ---
102 // The legacy code had an outer `pseudot` loop that only ran once. We preserve it.
103 for (PetscInt pseudot = 0; pseudot < 1; pseudot++) {
104 // Loop over each block to perform the RK stages
105 for (PetscInt bi = 0; bi < simCtx->block_number; bi++) {
106 for (istage = 0; istage < 4; istage++) {
107 LOG_ALLOW(LOCAL, LOG_DEBUG, " Block %d, RK Stage %d (alpha=%.4f)...\n", bi, istage, alfa[istage]);
108
109 // a. Calculate the Right-Hand Side (RHS) of the momentum equation.
110 ierr = ComputeRHS(&user[bi], user[bi].Rhs); CHKERRQ(ierr);
111
112 // b. Advance Ucont to the next intermediate stage using the RK coefficient.
113 // Ucont_new = Ucont_old + alpha * dt * RHS
114 ierr = VecWAXPY(user[bi].Ucont, alfa[istage] * dt, user[bi].Rhs, user[bi].Ucont_o); CHKERRQ(ierr);
115
116 // c. Update local ghost values for the new intermediate Ucont.
117 ierr = DMGlobalToLocalBegin(user[bi].fda, user[bi].Ucont, INSERT_VALUES, user[bi].lUcont); CHKERRQ(ierr);
118 ierr = DMGlobalToLocalEnd(user[bi].fda, user[bi].Ucont, INSERT_VALUES, user[bi].lUcont); CHKERRQ(ierr);
119
120 // d. Re-apply boundary conditions for the new intermediate velocity.
121 // This is crucial for the stability and accuracy of the multi-stage scheme.
122 ierr = ApplyBoundaryConditions(&user[bi]); CHKERRQ(ierr);
123
124 } // End of RK stages for one block
125
126 /*
127 // Final IBM Interpolation for the block (if enabled)
128 if (simCtx->immersed) {
129 LOG_ALLOW(LOCAL, LOG_DEBUG, " Performing post-RK IBM interpolation for block %d.\n", bi);
130 for (PetscInt ibi = 0; ibi < simCtx->NumberOfBodies; ibi++) {
131 ierr = ibm_interpolation_advanced(&user[bi], &ibm[ibi], ibi, 1); CHKERRQ(ierr);
132 }
133 }
134 */
135
136 } // End loop over blocks
137
138 // --- 3. Inter-Block Communication (Legacy Logic) ---
139 // This is called after all blocks have completed their RK stages.
140 if (simCtx->block_number > 1) {
141 // LOG_ALLOW(GLOBAL, LOG_DEBUG, "Updating multi-block interfaces after RK stages.\n");
142 // ierr = Block_Interface_U(user); CHKERRQ(ierr);
143 }
144
145 } // End of pseudo-time loop
146
147 // --- 4. Cleanup ---
148 // Destroy the RHS vectors that were created at the start of this function.
149 for (PetscInt bi = 0; bi < simCtx->block_number; bi++) {
150 ierr = VecDestroy(&user[bi].Rhs); CHKERRQ(ierr);
151 }
152
153 LOG_ALLOW(GLOBAL, LOG_INFO, "Runge-Kutta solve completed for all blocks.\n");
154
156
157 PetscFunctionReturn(0);
158}
PetscErrorCode ApplyBoundaryConditions(UserCtx *user)
Main boundary-condition orchestrator executed during solver timestepping.
#define LOCAL
Logging scope definitions for controlling message output.
Definition logging.h:44
#define GLOBAL
Scope for global logging across all processes.
Definition logging.h:45
#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:199
@ LOG_INFO
Informational messages about program execution.
Definition logging.h:30
@ LOG_DEBUG
Detailed debugging information.
Definition logging.h:31
PetscInt block_number
Definition variables.h:712
Here is the call graph for this function:
Here is the caller graph for this function:

◆ MomentumSolver_DualTime_Picard_RK4()

PetscErrorCode MomentumSolver_DualTime_Picard_RK4 ( UserCtx user,
IBMNodes ibm,
FSInfo fsi 
)

Internal helper implementation: MomentumSolver_DualTime_Picard_RK4().

Solves the Momentum Equations using Dual-Time Stepping with a Fixed-Point RK4 Smoother.

Local to this translation unit.

Definition at line 166 of file momentumsolvers.c.

167{
168 (void)ibm;
169 (void)fsi;
170 // --- CONTEXT ACQUISITION BLOCK ---
171 SimCtx *simCtx = user[0].simCtx;
172 const PetscInt block_number = simCtx->block_number;
173
174 // Legacy Names (Physics/Time) - Kept for brevity in formulas
175 const PetscInt ti = simCtx->step;
176 const PetscReal dt = simCtx->dt; // Physical Time Step
177 const PetscReal cfl = simCtx->pseudo_cfl; // Initial CFL Scaling
178 const PetscReal alfa[] = {0.25, 1.0/3.0, 0.5, 1.0}; // RK4 Coefficients
179 Vec *pRhs; // Per-block backup of Rhs at last accepted pseudo-state.
180
181 // State flags
182 PetscBool force_restart = PETSC_FALSE;
183
184 // Renamed Solver Parameters
185 const PetscInt max_pseudo_steps = simCtx->mom_max_pseudo_steps; // Max Dual-Time Iterations
186 const PetscReal tol_abs_delta = simCtx->mom_atol; // Stop if |dU| < tol
187 const PetscReal tol_rtol_delta = simCtx->mom_rtol; // Stop if |dU|/|dU0| < tol
188 // --- END CONTEXT ACQUISITION BLOCK ---
189
190 PetscErrorCode ierr;
191 PetscMPIInt rank;
192 PetscInt istage, pseudo_iter;
193 PetscReal ts, te, cput;
194
195 // --- Global Convergence Metrics ---
196 PetscReal global_norm_delta = 10.0;
197 PetscReal global_rel_delta = 1.0;
198 PetscReal global_rel_resid = 1.0;
199
200 // --- Local Metric Arrays (Per Block) ---
201 PetscReal *delta_sol_norm_init, *delta_sol_norm_prev, *delta_sol_norm_curr, *delta_sol_rel_curr;
202 PetscReal *resid_norm_init, *resid_norm_prev, *resid_norm_curr, *resid_rel_curr;
203 PetscReal *pseudo_dt_scaling;
204
205 PetscFunctionBeginUser;
207 ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank); CHKERRQ(ierr);
208 LOG_ALLOW(GLOBAL, LOG_INFO, "Executing Dual-Time Momentum Solver for %d block(s).\n", block_number);
209
210 // --- Allocate metric arrays ---
211 ierr = PetscMalloc2(block_number, &delta_sol_norm_init, block_number, &delta_sol_norm_prev); CHKERRQ(ierr);
212 ierr = PetscMalloc2(block_number, &delta_sol_norm_curr, block_number, &delta_sol_rel_curr); CHKERRQ(ierr);
213 ierr = PetscMalloc2(block_number, &resid_norm_init, block_number, &resid_norm_prev); CHKERRQ(ierr);
214 ierr = PetscMalloc2(block_number, &resid_norm_curr, block_number, &pseudo_dt_scaling); CHKERRQ(ierr);
215 ierr = PetscMalloc1(block_number, &resid_rel_curr); CHKERRQ(ierr);
216 ierr = PetscMalloc1(block_number, &pRhs); CHKERRQ(ierr);
217
218 ierr = PetscTime(&ts); CHKERRQ(ierr);
219
220 // --- 1. Pre-Loop Initialization ---
221
222 if (block_number > 1) {
223 // ierr = Block_Interface_U(user); CHKERRQ(ierr);
224 }
225
226 for (PetscInt bi = 0; bi < block_number; bi++) {
227 ierr = ApplyBoundaryConditions(&user[bi]); CHKERRQ(ierr);
228
229 // Immersed boundary interpolation (if enabled)
230
231 // Allocate workspace
232 ierr = VecDuplicate(user[bi].Ucont, &user[bi].Rhs); CHKERRQ(ierr);
233 ierr = VecDuplicate(user[bi].Rhs, &pRhs[bi]); CHKERRQ(ierr);
234 ierr = VecDuplicate(user[bi].Ucont, &user[bi].dUcont); CHKERRQ(ierr);
235 ierr = VecDuplicate(user[bi].Ucont, &user[bi].pUcont); CHKERRQ(ierr);
236
237 // Initialize Backup (pUcont) with current state
238 ierr = VecCopy(user[bi].Ucont, user[bi].pUcont); CHKERRQ(ierr);
239
240 // --- Calculate Initial Total Residual (Spatial + Temporal) ---
241 ierr = ComputeTotalResidual(&user[bi]); CHKERRQ(ierr);
242
243 // Backup Rhs with current state
244 ierr = VecCopy(user[bi].Rhs, pRhs[bi]); CHKERRQ(ierr);
245
246 // Compute Initial Norms
247 ierr = VecNorm(user[bi].Rhs, NORM_INFINITY, &resid_norm_init[bi]); CHKERRQ(ierr);
248
249 // Initialize history for backtracking logic
250 resid_norm_prev[bi] = resid_norm_init[bi]; // Meaningful ratio on first iteration
251 delta_sol_norm_prev[bi] = 1000.0;
252 pseudo_dt_scaling[bi] = cfl; // Initialize adaptive scalar
253
254 LOG_ALLOW(GLOBAL,LOG_INFO," Block %d | Max RHS = %.6f | Pseudo-CFL = %.4f .\n", bi, resid_norm_init[bi], cfl);
255 }
256
257 // --- 2. Main Pseudo-Time Iteration Loop ---
258 pseudo_iter = 0;
259 while ((force_restart || (global_norm_delta > tol_abs_delta || global_rel_delta > tol_rtol_delta) || pseudo_iter < 1)
260 && pseudo_iter < max_pseudo_steps)
261 {
262 pseudo_iter++;
263 force_restart = PETSC_FALSE;
264
265 for (PetscInt bi = 0; bi < block_number; bi++) {
266
267 // === 4-Stage Explicit Runge-Kutta Loop ===
268 for (istage = 0; istage < 4; istage++) {
269
270 LOG_ALLOW(GLOBAL,LOG_TRACE," Pseudo-Iter: %d | RK-Stage: %d\n", pseudo_iter, istage);
271
272 // RK Update: U_new = U_old + (Scaler * Alpha * dt) * Total_Residual
273 ierr = VecWAXPY(user[bi].Ucont,
274 pseudo_dt_scaling[bi] * alfa[istage] * dt,
275 user[bi].Rhs,
276 user[bi].pUcont); CHKERRQ(ierr);
277
278 // Sync Ghosts & Re-apply BCs for intermediate stage
279 ierr = UpdateLocalGhosts(&user[bi],"Ucont"); CHKERRQ(ierr);
280
281 //ierr = DMGlobalToLocalBegin(user[bi].fda, user[bi].Ucont, INSERT_VALUES, user[bi].lUcont); CHKERRQ(ierr);
282 //ierr = DMGlobalToLocalEnd(user[bi].fda, user[bi].Ucont, INSERT_VALUES, user[bi].lUcont); CHKERRQ(ierr);
283
284 ierr = ApplyBoundaryConditions(&user[bi]); CHKERRQ(ierr);
285
286 // --- Re-calculate Total Residual for next stage ---
287 ierr = ComputeTotalResidual(&user[bi]); CHKERRQ(ierr);
288
289 } // End RK Stages
290
291 // Immersed boundary interpolation (if enabled)
292
293 // === Convergence Metrics Calculation ===
294
295 // Calculate dU = U_current - U_backup
296 ierr = VecWAXPY(user[bi].dUcont, -1.0, user[bi].pUcont, user[bi].Ucont); CHKERRQ(ierr);
297
298 // Compute Infinity Norms
299 ierr = VecNorm(user[bi].dUcont, NORM_INFINITY, &delta_sol_norm_curr[bi]); CHKERRQ(ierr);
300 ierr = VecNorm(user[bi].Rhs, NORM_INFINITY, &resid_norm_curr[bi]); CHKERRQ(ierr);
301
302 // Normalize relative metrics
303 if (pseudo_iter == 1) {
304 delta_sol_norm_init[bi] = delta_sol_norm_curr[bi];
305 delta_sol_rel_curr[bi] = 1.0;
306 resid_rel_curr[bi] = 1.0;
307 // resid_norm_init[bi] set correctly before the loop — do not overwrite
308 } else {
309 if (delta_sol_norm_init[bi] > 1.0e-10)
310 delta_sol_rel_curr[bi] = delta_sol_norm_curr[bi] / delta_sol_norm_init[bi];
311 else
312 delta_sol_rel_curr[bi] = 0.0;
313
314 if(resid_norm_init[bi] > 1.0e-10)
315 resid_rel_curr[bi] = resid_norm_curr[bi] / resid_norm_init[bi];
316 else
317 resid_rel_curr[bi] = 0.0;
318 }
319
320 // --- File Logging ---
321 if (!rank) {
322 FILE *f;
323 char filen[PETSC_MAX_PATH_LEN + 128];
324 ierr = PetscSNPrintf(filen, sizeof(filen), "%s/Momentum_Solver_Convergence_History_Block_%1d.log", simCtx->log_dir, bi); CHKERRQ(ierr);
325 if(simCtx->step == simCtx->StartStep + 1 && pseudo_iter == 1 && !simCtx->continueMode) f = fopen(filen, "w");
326 else f = fopen(filen, "a");
327
328 PetscFPrintf(PETSC_COMM_WORLD, f, "Step: %d | PseudoIter(k): %d| | Pseudo-cfl: %.4f |dUk|: %le | |dUk|/|dUprev|: %le | |Rk|: %le | |Rk|/|Rprev|: %le \n",
329 (int)ti, (int)pseudo_iter, pseudo_dt_scaling[bi], delta_sol_norm_curr[bi], delta_sol_rel_curr[bi], resid_norm_curr[bi],resid_rel_curr[bi]);
330 fclose(f);
331 }
332 } // End loop over blocks
333
334 // --- Update Global Convergence Criteria ---
335 global_norm_delta = -1.0e20;
336 global_rel_delta = -1.0e20;
337 global_rel_resid = -1.0e20;
338
339 for (PetscInt bi = 0; bi < block_number; bi++) {
340 global_norm_delta = PetscMax(delta_sol_norm_curr[bi], global_norm_delta);
341 global_rel_delta = PetscMax(delta_sol_rel_curr[bi], global_rel_delta);
342 global_rel_resid = PetscMax(resid_rel_curr[bi], global_rel_resid);
343 }
344 ierr = PetscTime(&te); CHKERRQ(ierr);
345 cput = te - ts;
346 LOG_ALLOW(GLOBAL, LOG_INFO, " Pseudo-Iter(k) %d: |dUk|=%e, |dUk|/|dU0| = %e, |Rk|/|R0| = %e, CPU=%.2fs\n",
347 pseudo_iter, global_norm_delta, global_rel_delta, global_rel_resid, cput);
348
349 // === Backtracking Line Search ===
350 for (PetscInt bi = 0; bi < block_number; bi++) {
351
352 PetscReal resid_ratio = resid_norm_curr[bi] / resid_norm_prev[bi];
353
354 PetscBool is_sol_nan, is_resid_nan, is_diverging;
355 // Check for NaNs in Solution and Residual
356 is_sol_nan = (PetscBool)PetscIsNanReal(delta_sol_norm_curr[bi]);
357 is_resid_nan = (PetscBool)PetscIsNanReal(resid_norm_curr[bi]);
358
359 // TRIGGER: Divergence Detected
360 is_diverging = (PetscBool)(pseudo_iter > 1 &&
361 resid_norm_curr[bi] > simCtx->mom_dt_rk4_residual_norm_noise_allowance_factor*resid_norm_prev[bi] &&
362 delta_sol_norm_curr[bi] > delta_sol_norm_prev[bi]
363 );
364
365 if(is_diverging || is_sol_nan || is_resid_nan) {
366
367 if(is_sol_nan){
368 LOG_ALLOW(LOCAL, LOG_WARNING, " Block %d: NaN detected in solution update norm.\n", bi);
369 }
370 if(is_resid_nan){
371 LOG_ALLOW(LOCAL, LOG_WARNING, " Block %d: NaN detected in residual norm.\n", bi);
372 }
373 if(is_diverging){
374 LOG_ALLOW(LOCAL, LOG_WARNING, " Block %d: Divergence detected (|R| and |dU| increasing).\n", bi);
375 }
376
377 if(pseudo_dt_scaling[bi] < simCtx->min_pseudo_cfl) {
378 SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_CONV_FAILED,
379 "Block %d: Solver diverged (NaN or explosion) and min time-step reached. Simulation aborted.", bi);
380 }
381
382 // For nan situations, we always backtrack.
383 force_restart = PETSC_TRUE;
384
385 // 1. Reduce Step Size
386 pseudo_dt_scaling[bi] *= 0.5*simCtx->pseudo_cfl_reduction_factor; // Aggressive deceleration (1/2 x reduction factor).
387 // NOTE: pseudo_iter-- is handled once after this loop to avoid multi-decrement for block_number > 1.
388
389 LOG_ALLOW(LOCAL, LOG_DEBUG, " Block %d: Backtracking triggered. Reducing scaling to %f.\n", bi, pseudo_dt_scaling[bi]);
390
391 // 3. RESTORE from Backup
392 ierr = VecCopy(user[bi].pUcont, user[bi].Ucont); CHKERRQ(ierr);
393
394 ierr = VecCopy(pRhs[bi], user[bi].Rhs); CHKERRQ(ierr);
395
396 // 4. SYNC Ghost Cells
397 ierr = UpdateLocalGhosts(&user[bi],"Ucont"); CHKERRQ(ierr);
398 //ierr = DMGlobalToLocalBegin(user[bi].fda, user[bi].Ucont, INSERT_VALUES, user[bi].lUcont); CHKERRQ(ierr);
399 //ierr = DMGlobalToLocalEnd(user[bi].fda, user[bi].Ucont, INSERT_VALUES, user[bi].lUcont); CHKERRQ(ierr);
400
401 // 5. Reset current norms to avoid polluting next loop check.
402 delta_sol_norm_curr[bi] = 1000.0;
403 resid_norm_curr[bi] = 1000.0;
404
405 }
406 else {
407 // SUCCESS: Step accepted. Update Backup and History.
408 ierr = VecCopy(user[bi].Ucont, user[bi].pUcont); CHKERRQ(ierr);
409 ierr = VecCopy(user[bi].Rhs, pRhs[bi]); CHKERRQ(ierr);
410 resid_norm_prev[bi] = resid_norm_curr[bi];
411 delta_sol_norm_prev[bi] = delta_sol_norm_curr[bi];
412
413 // Adaptive Ramping of solution speed (Pseudo-CFL)
414
415 PetscBool is_converging_fast;
416 is_converging_fast = (PetscBool)(resid_ratio < 0.9);
417 PetscBool is_stagnating;
418 is_stagnating = (PetscBool)(resid_ratio >= 0.98 && resid_ratio < 1.0);
419 PetscBool is_noisy;
420 is_noisy = (PetscBool)(resid_ratio >= 1.0 && resid_ratio < simCtx->mom_dt_rk4_residual_norm_noise_allowance_factor);
421
422 PetscReal current_cfl = pseudo_dt_scaling[bi];
423 PetscReal cfl_floor = 10*simCtx->min_pseudo_cfl; // Define a "low CFL" threshold for more aggressive acceleration.
424 PetscReal cfl_ceiling = 0.75*simCtx->max_pseudo_cfl; // Define a "high CFL" threshold for mild acceleration.
425
426 if(is_converging_fast){
427 // Residual is falling fast (good convergence)
428 pseudo_dt_scaling[bi] *= simCtx->pseudo_cfl_growth_factor;
429 } else if(is_stagnating && current_cfl < cfl_floor){
430 // Residual is falling very slow (weak convergence.) we accelerate at low CFL
431 pseudo_dt_scaling[bi] *= PetscMin(1.2*simCtx->pseudo_cfl_growth_factor, 1.1); // more aggressive acceleration at low CFL
432 } else if(is_stagnating && current_cfl >= cfl_floor && current_cfl < cfl_ceiling && !is_noisy){
433 // moderate CFL - normal acceleration
434 pseudo_dt_scaling[bi] *= simCtx->pseudo_cfl_growth_factor;
435 } else if(is_stagnating && current_cfl >= cfl_ceiling && !is_noisy){
436 // high CFL - mild acceleration (creep)
437 pseudo_dt_scaling[bi] *= PetscMax(0.95*simCtx->pseudo_cfl_growth_factor,1.005);
438 } else if(is_noisy){
439 // Residual is rising slightly (within noise allowance). we decelerate mildly.
440 if(current_cfl > cfl_ceiling){
441 // High CFL - more aggressive deceleration
442 pseudo_dt_scaling[bi] *= PetscMax(0.9*simCtx->pseudo_cfl_reduction_factor,0.9);
443 } else if(current_cfl > cfl_floor && current_cfl <= cfl_ceiling){
444 // Moderate CFL - mild deceleration
445 pseudo_dt_scaling[bi] *= simCtx->pseudo_cfl_reduction_factor;
446 } else if(current_cfl <= cfl_floor){
447 // Low CFL - gentle acceleration
448 // rise in residual likely due to noise, so we back off slightly.
449 pseudo_dt_scaling[bi] *= PetscMax(0.95*simCtx->pseudo_cfl_growth_factor,1.01);
450 }
451 }
452 // Clamp to max bound.
453 pseudo_dt_scaling[bi] = PetscMin(pseudo_dt_scaling[bi], simCtx->max_pseudo_cfl);
454 pseudo_dt_scaling[bi] = PetscMax(pseudo_dt_scaling[bi], simCtx->min_pseudo_cfl);
455 }
456 }
457
458 // Decrement once regardless of how many blocks backtracked.
459 if (force_restart) pseudo_iter--;
460
461 if (block_number > 1) {
462 // ierr = Block_Interface_U(user); CHKERRQ(ierr);
463 }
464 } // End while loop
465
466
467 // Smart CFL Carry Over
468
469 PetscReal local_min_cfl = 1.0e20;
470 PetscReal global_min_cfl;
471
472 PetscReal cfl_floor = 10*simCtx->min_pseudo_cfl;
473
474 // 1. Find the "Weakest Link" (lowest CFL) across local blocks
475 for (PetscInt bi = 0; bi < block_number; bi++) {
476 if (pseudo_dt_scaling[bi] < local_min_cfl) {
477 local_min_cfl = pseudo_dt_scaling[bi];
478 }
479 }
480
481 // 2. MPI Reduction (Crucial if blocks are distributed across ranks)
482 // We want the minimum CFL across the entire simulation, not just this processor.
483 ierr = MPI_Allreduce(&local_min_cfl, &global_min_cfl, 1, MPIU_REAL, MPI_MIN, PETSC_COMM_WORLD); CHKERRQ(ierr);
484
485 // 3. Apply Safety Factor & Update simCtx
486 // We use a factor (e.g., 0.85) to back off slightly for the start of the next timestep
487 // to accommodate the initial residual spike from the BDF term.
488
489 PetscReal cfl_carry_over_relaxation_factor = 1.0; // Could be parameterized in simCtx if desired.
490
491 PetscReal next_start_cfl;
492
493 if(global_min_cfl < cfl_floor){
494 // If the minimum CFL is very low, we back off more aggressively.
495 next_start_cfl = cfl_floor * 1.5;
496 } else{
497 next_start_cfl = global_min_cfl * cfl_carry_over_relaxation_factor;
498 }
499
500 next_start_cfl = PetscMin(next_start_cfl, simCtx->max_pseudo_cfl);
501 next_start_cfl = PetscMax(next_start_cfl, simCtx->min_pseudo_cfl);
502
503 // 5. Save for next Physical Time Step
504 if (!rank) {
505 LOG_ALLOW(GLOBAL, LOG_INFO, " CFL Adaptation: Step %d finished at CFL %.3f. Next step will start at CFL %.3f.\n",
506 simCtx->step, global_min_cfl, next_start_cfl);
507 }
508
509 // Update the context so the next call to MomentumSolver reads this new value
510 simCtx->pseudo_cfl = next_start_cfl;
511
512 // --- Final Cleanup ---
513 for (PetscInt bi = 0; bi < block_number; bi++) {
514 ierr = VecDestroy(&user[bi].Rhs); CHKERRQ(ierr);
515 ierr = VecDestroy(&user[bi].dUcont); CHKERRQ(ierr);
516 ierr = VecDestroy(&user[bi].pUcont); CHKERRQ(ierr);
517 ierr = VecDestroy(&pRhs[bi]); CHKERRQ(ierr);
518 }
519 ierr = PetscFree(pRhs); CHKERRQ(ierr);
520
521 ierr = PetscFree2(delta_sol_norm_init, delta_sol_norm_prev);CHKERRQ(ierr);
522 ierr = PetscFree2(delta_sol_norm_curr, delta_sol_rel_curr);CHKERRQ(ierr);
523 ierr = PetscFree2(resid_norm_init, resid_norm_prev);CHKERRQ(ierr);
524 ierr = PetscFree2(resid_norm_curr, pseudo_dt_scaling); CHKERRQ(ierr);
525 ierr = PetscFree(resid_rel_curr); CHKERRQ(ierr);
526
527 LOG_ALLOW(GLOBAL, LOG_INFO, "Momentum Solver converged/completed after %d pseudo-iterations.\n", pseudo_iter);
529 PetscFunctionReturn(0);
530}
@ LOG_TRACE
Very fine-grained tracing information for in-depth debugging.
Definition logging.h:32
@ LOG_WARNING
Non-critical issues that warrant attention.
Definition logging.h:29
static PetscErrorCode ComputeTotalResidual(UserCtx *user)
Internal helper implementation: ComputeTotalResidual().
PetscErrorCode UpdateLocalGhosts(UserCtx *user, const char *fieldName)
Updates the local vector (including ghost points) from its corresponding global vector.
Definition setup.c:1361
PetscBool continueMode
Definition variables.h:660
PetscReal mom_rtol
Definition variables.h:685
PetscReal pseudo_cfl_reduction_factor
Definition variables.h:692
PetscReal min_pseudo_cfl
Definition variables.h:693
PetscReal mom_atol
Definition variables.h:685
PetscReal max_pseudo_cfl
Definition variables.h:693
char log_dir[PETSC_MAX_PATH_LEN]
Definition variables.h:668
PetscReal pseudo_cfl_growth_factor
Definition variables.h:692
PetscReal mom_dt_rk4_residual_norm_noise_allowance_factor
Definition variables.h:694
PetscInt mom_max_pseudo_steps
Definition variables.h:684
PetscReal pseudo_cfl
Definition variables.h:691
Here is the call graph for this function:
Here is the caller graph for this function: