PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
Loading...
Searching...
No Matches
momentumsolvers.h
Go to the documentation of this file.
1#ifndef MOMENTUNSOLVERS_H
2#define MOMENTUMSOLVERS_H
3
4#include "variables.h" // Provides definitions for UserCtx, SimCtx, IBMNodes, etc.
5#include "logging.h"
6#include "rhs.h"
7#include "Boundaries.h"
8
9/*================================================================================*
10 * MOMENTUM EQUATION SOLVERS *
11 *================================================================================*/
12
13/**
14 * @brief Advances the momentum equations using an explicit 4th-order Runge-Kutta scheme.
15 * @param user Array of UserCtx structs for all blocks.
16 * @param ibm (Optional) Pointer to IBM data. Pass NULL if disabled.
17 * @param fsi (Optional) Pointer to FSI data. Pass NULL if disabled.
18 * @return PetscErrorCode 0 on success.
19 *
20 * @note Testing status:
21 * The explicit RK path remains on the near-term backlog for direct
22 * positive-path bespoke coverage; today it is weaker than the dual-time
23 * path in the test surface.
24 */
25extern PetscErrorCode MomentumSolver_Explicit_RungeKutta4(UserCtx *user, IBMNodes *ibm, FSInfo *fsi);
26
27/**
28 * @brief Solves the Momentum Equations using Dual-Time Stepping with a Fixed-Point RK4 Smoother.
29 *
30 * =================================================================================================
31 * GLOSSARY & THEORETICAL BASIS
32 * =================================================================================================
33 * 1. METHODOLOGY: Dual-Time Stepping (Pseudo-Time Integration)
34 * We aim to solve the implicit BDF equation: R_spatial(U) + dU/dt_physical = 0.
35 * We do this by introducing a fictitious "Pseudo-Time" (tau) and iterating to steady state:
36 * dU/d(tau) = - [ R_spatial(U) + BDF_Terms(U) ]
37 * When dU/d(tau) -> 0, the physical time step is satisfied.
38 * 2. ALGORITHM: Fixed-Point Iteration with Explicit Runge-Kutta
39 * This is technically a Fixed-Point iteration on the operator:
40 * U_new = U_old + pseudo_dt_scaling * dt_pseudo * Total_Residual(U_old)
41 * We use a 4-Stage Explicit RK scheme (Jameson-Schmidt-Turkel coeffs) to smooth errors.
42 * 3. STABILITY: Backtracking Line Search
43 * If a pseudo-time step causes the Residual or Solution Error to GROW (Divergence),
44 * the solver "Backtracks": it restores the previous solution, cuts the pseudo-time step
45 * scaling factor (pseudo_dt_scaling) in half, and retries the iteration.
46 * =================================================================================================
47 * VARIABLE MAPPING
48 * =================================================================================================
49 * -- Physics Variables (Legacy Names Kept) --
50 * ti : Physical Time Step Index.
51 * dt : Physical Time Step size (Delta t).
52 * st : Pseudo-Time Step size (Delta tau).
53 * alfa : Runge-Kutta stage coefficients {1/4, 1/3, 1/2, 1}.
54 * -- Convergence & Solver Control (Renamed) --
55 * pseudo_iter : Counter for the inner dual-time loop.
56 * pseudo_dt_scaling : Adaptive scalar for the pseudo-time step (formerly lambda).
57 * delta_sol_norm : The L_inf norm of the change in solution (dU).
58 * resid_norm : The L_inf norm of the Total Residual (RHS).
59 * =================================================================================================
60 *
61 * @param user Primary `UserCtx` input for the operation.
62 * @param ibm Parameter `ibm` passed to `MomentumSolver_DualTime_Picard_RK4()`.
63 * @param fsi Parameter `fsi` passed to `MomentumSolver_DualTime_Picard_RK4()`.
64 * @return PetscErrorCode 0 on success.
65 *
66 * @note Testing status:
67 * This solver is covered primarily through runtime smoke and orchestration
68 * tests. A smaller direct invariant-style positive-path harness remains
69 * part of the next-gap backlog.
70 */
71PetscErrorCode MomentumSolver_DualTime_Picard_RK4(UserCtx *user, IBMNodes *ibm, FSInfo *fsi);
72
73#endif // MOMENTUMSOLVERS_H
Logging utilities and macros for PETSc-based applications.
PetscErrorCode MomentumSolver_DualTime_Picard_RK4(UserCtx *user, IBMNodes *ibm, FSInfo *fsi)
Solves the Momentum Equations using Dual-Time Stepping with a Fixed-Point RK4 Smoother.
PetscErrorCode MomentumSolver_Explicit_RungeKutta4(UserCtx *user, IBMNodes *ibm, FSInfo *fsi)
Advances the momentum equations using an explicit 4th-order Runge-Kutta scheme.
Main header file for a complex fluid dynamics solver.
Holds all data related to the state and motion of a body in FSI.
Definition variables.h:445
Represents a collection of nodes forming a surface for the IBM.
Definition variables.h:372
User-defined context containing data specific to a single computational grid level.
Definition variables.h:811