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

Go to the source code of this file.

Functions

PetscErrorCode RungeKutta (UserCtx *user, IBMNodes *ibm, FSInfo *fsi)
 Advances the momentum equations using an explicit 4th-order Runge-Kutta scheme.
 
PetscErrorCode ImplicitMomentumSolver (UserCtx *user, IBMNodes *ibm, FSInfo *fsi)
 Advances the momentum equations using an implicit scheme.
 

Function Documentation

◆ RungeKutta()

PetscErrorCode RungeKutta ( UserCtx user,
IBMNodes ibm,
FSInfo fsi 
)
extern

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

Parameters
userArray of UserCtx structs for all blocks.
ibm(Optional) Pointer to IBM data. Pass NULL if disabled.
fsi(Optional) Pointer to FSI data. Pass NULL if disabled.
Returns
PetscErrorCode 0 on success.

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

This function computes an intermediate, non-divergence-free contravariant velocity field (Ucont) at time t_{n+1} for all computational blocks.

This is a minimally-edited version of the legacy solver. It retains its internal loop over all blocks and is intended to be called once per time step from the main Flow_Solver orchestrator. All former global variables are now accessed via the SimCtx passed in through the first block's UserCtx.

Parameters
userThe array of UserCtx structs for all blocks.
ibm(Optional) Pointer to the full array of IBM data structures. Pass NULL if disabled.
fsi(Optional) Pointer to the full array of FSI data structures. Pass NULL if disabled.
Returns
PetscErrorCode 0 on success.

Definition at line 22 of file implicitsolvers.c.

23{
24 PetscErrorCode ierr;
25 // --- Context Acquisition ---
26 // Get the master simulation context from the first block's UserCtx.
27 // This is the bridge to access all former global variables.
28 SimCtx *simCtx = user[0].simCtx;
29 PetscReal dt = simCtx->dt;
30 PetscReal st = simCtx->st;
31 PetscInt istage;
32 PetscReal alfa[] = {0.25, 1.0/3.0, 0.5, 1.0};
33
34 PetscFunctionBeginUser;
35 LOG_ALLOW(GLOBAL, LOG_INFO, "Executing explicit momentum solver (Runge-Kutta) for %d block(s).\n",simCtx->block_number);
36
37 // --- 1. Pre-Loop Initialization (Legacy Logic) ---
38 // This block prepares boundary conditions and allocates the RHS vector for all blocks
39 // before the main RK loop begins. This logic is preserved from the original code.
40 LOG_ALLOW(GLOBAL, LOG_DEBUG, "Preparing all blocks for Runge-Kutta solve...\n");
41 for (PetscInt bi = 0; bi < simCtx->block_number; bi++) {
42 ierr = InflowFlux(&user[bi]); CHKERRQ(ierr);
43 ierr = OutflowFlux(&user[bi]); CHKERRQ(ierr);
44 ierr = FormBCS(&user[bi]); CHKERRQ(ierr);
45
46 /*
47 // Immersed boundary interpolation (if enabled)
48 if (simCtx->immersed) {
49 LOG_ALLOW(LOCAL, LOG_DEBUG, " Performing pre-RK IBM interpolation for block %d.\n", bi);
50 for (PetscInt ibi = 0; ibi < simCtx->NumberOfBodies; ibi++) {
51 // The 'ibm' and 'fsi' pointers are passed directly from Flow_Solver
52 ierr = ibm_interpolation_advanced(&user[bi], &ibm[ibi], ibi, 1); CHKERRQ(ierr);
53 }
54 }
55 */
56
57 // Allocate the persistent RHS vector for this block's context
58 ierr = VecDuplicate(user[bi].Ucont, &user[bi].Rhs); CHKERRQ(ierr);
59 }
60 LOG_ALLOW(GLOBAL, LOG_DEBUG, "Pre-loop initialization complete.\n");
61
62 // --- 2. Main Runge-Kutta Loop ---
63 // The legacy code had an outer `pseudot` loop that only ran once. We preserve it.
64 for (PetscInt pseudot = 0; pseudot < 1; pseudot++) {
65 // Loop over each block to perform the RK stages
66 for (PetscInt bi = 0; bi < simCtx->block_number; bi++) {
67 for (istage = 0; istage < 4; istage++) {
68 LOG_ALLOW(LOCAL, LOG_DEBUG, " Block %d, RK Stage %d (alpha=%.4f)...\n", bi, istage, alfa[istage]);
69
70 // a. Calculate the Right-Hand Side (RHS) of the momentum equation.
71 ierr = FormFunction1(&user[bi], user[bi].Rhs); CHKERRQ(ierr);
72
73 // b. Advance Ucont to the next intermediate stage using the RK coefficient.
74 // Ucont_new = Ucont_old + alpha * dt * RHS
75 ierr = VecWAXPY(user[bi].Ucont, alfa[istage] * dt * st, user[bi].Rhs, user[bi].Ucont_o); CHKERRQ(ierr);
76
77 // c. Update local ghost values for the new intermediate Ucont.
78 ierr = DMGlobalToLocalBegin(user[bi].fda, user[bi].Ucont, INSERT_VALUES, user[bi].lUcont); CHKERRQ(ierr);
79 ierr = DMGlobalToLocalEnd(user[bi].fda, user[bi].Ucont, INSERT_VALUES, user[bi].lUcont); CHKERRQ(ierr);
80
81 // d. Re-apply boundary conditions for the new intermediate velocity.
82 // This is crucial for the stability and accuracy of the multi-stage scheme.
83 ierr = InflowFlux(&user[bi]); CHKERRQ(ierr);
84 ierr = OutflowFlux(&user[bi]); CHKERRQ(ierr);
85 ierr = FormBCS(&user[bi]); CHKERRQ(ierr);
86
87 } // End of RK stages for one block
88
89 /*
90 // Final IBM Interpolation for the block (if enabled)
91 if (simCtx->immersed) {
92 LOG_ALLOW(LOCAL, LOG_DEBUG, " Performing post-RK IBM interpolation for block %d.\n", bi);
93 for (PetscInt ibi = 0; ibi < simCtx->NumberOfBodies; ibi++) {
94 ierr = ibm_interpolation_advanced(&user[bi], &ibm[ibi], ibi, 1); CHKERRQ(ierr);
95 }
96 }
97 */
98
99 } // End loop over blocks
100
101 // --- 3. Inter-Block Communication (Legacy Logic) ---
102 // This is called after all blocks have completed their RK stages.
103 if (simCtx->block_number > 1) {
104 // LOG_ALLOW(GLOBAL, LOG_DEBUG, "Updating multi-block interfaces after RK stages.\n");
105 // ierr = Block_Interface_U(user); CHKERRQ(ierr);
106 }
107
108 } // End of pseudo-time loop
109
110 // --- 4. Cleanup ---
111 // Destroy the RHS vectors that were created at the start of this function.
112 for (PetscInt bi = 0; bi < simCtx->block_number; bi++) {
113 ierr = VecDestroy(&user[bi].Rhs); CHKERRQ(ierr);
114 }
115
116 LOG_ALLOW(GLOBAL, LOG_INFO, "Runge-Kutta solve completed for all blocks.\n");
117 PetscFunctionReturn(0);
118}
PetscErrorCode OutflowFlux(UserCtx *user)
PetscErrorCode InflowFlux(UserCtx *user)
Definition Boundaries.c:819
PetscErrorCode FormBCS(UserCtx *user)
#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:207
@ LOG_INFO
Informational messages about program execution.
Definition logging.h:32
@ LOG_DEBUG
Detailed debugging information.
Definition logging.h:33
PetscErrorCode FormFunction1(UserCtx *user, Vec Rhs)
Computes the Right-Hand Side (RHS) of the momentum equations.
Definition rhs.c:1231
PetscInt block_number
Definition variables.h:562
SimCtx * simCtx
Back-pointer to the master simulation context.
Definition variables.h:633
PetscReal st
Definition variables.h:549
PetscReal dt
Definition variables.h:527
The master context for the entire simulation.
Definition variables.h:513
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ImplicitMomentumSolver()

PetscErrorCode ImplicitMomentumSolver ( UserCtx user,
IBMNodes ibm,
FSInfo fsi 
)
extern

Advances the momentum equations using an implicit scheme.

Parameters
userArray of UserCtx structs for all blocks.
ibm(Optional) Pointer to IBM data. Pass NULL if disabled.
fsi(Optional) Pointer to FSI data. Pass NULL if disabled.
Returns
PetscErrorCode 0 on success.