PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
Loading...
Searching...
No Matches
implicitsolvers.c
Go to the documentation of this file.
1#include "implicitsolvers.h"
2
3#undef __FUNCT__
4#define __FUNCT__ "RungeKutta"
5/**
6 * @brief Advances the momentum equations for ALL blocks by one time step
7 * using an explicit 4-stage, 4th-order Runge-Kutta scheme.
8 *
9 * This function computes an intermediate, non-divergence-free contravariant
10 * velocity field (Ucont) at time t_{n+1} for all computational blocks.
11 *
12 * This is a minimally-edited version of the legacy solver. It retains its
13 * internal loop over all blocks and is intended to be called once per time step
14 * from the main FlowSolver orchestrator. All former global variables are now
15 * accessed via the SimCtx passed in through the first block's UserCtx.
16 *
17 * @param user The array of UserCtx structs for all blocks.
18 * @param ibm (Optional) Pointer to the full array of IBM data structures. Pass NULL if disabled.
19 * @param fsi (Optional) Pointer to the full array of FSI data structures. Pass NULL if disabled.
20 * @return PetscErrorCode 0 on success.
21 */
22PetscErrorCode RungeKutta(UserCtx *user, IBMNodes *ibm, FSInfo *fsi)
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
37
38 LOG_ALLOW(GLOBAL, LOG_INFO, "Executing explicit momentum solver (Runge-Kutta) for %d block(s).\n",simCtx->block_number);
39
40 // --- 1. Pre-Loop Initialization (Legacy Logic) ---
41 // This block prepares boundary conditions and allocates the RHS vector for all blocks
42 // before the main RK loop begins. This logic is preserved from the original code.
43 LOG_ALLOW(GLOBAL, LOG_DEBUG, "Preparing all blocks for Runge-Kutta solve...\n");
44 for (PetscInt bi = 0; bi < simCtx->block_number; bi++) {
45 ierr = InflowFlux(&user[bi]); CHKERRQ(ierr);
46 ierr = OutflowFlux(&user[bi]); CHKERRQ(ierr);
47 ierr = FormBCS(&user[bi]); CHKERRQ(ierr);
48
49 /*
50 // Immersed boundary interpolation (if enabled)
51 if (simCtx->immersed) {
52 LOG_ALLOW(LOCAL, LOG_DEBUG, " Performing pre-RK IBM interpolation for block %d.\n", bi);
53 for (PetscInt ibi = 0; ibi < simCtx->NumberOfBodies; ibi++) {
54 // The 'ibm' and 'fsi' pointers are passed directly from FlowSolver
55 ierr = ibm_interpolation_advanced(&user[bi], &ibm[ibi], ibi, 1); CHKERRQ(ierr);
56 }
57 }
58 */
59
60 // Allocate the persistent RHS vector for this block's context
61 ierr = VecDuplicate(user[bi].Ucont, &user[bi].Rhs); CHKERRQ(ierr);
62 }
63 LOG_ALLOW(GLOBAL, LOG_DEBUG, "Pre-loop initialization complete.\n");
64
65 // --- 2. Main Runge-Kutta Loop ---
66 // The legacy code had an outer `pseudot` loop that only ran once. We preserve it.
67 for (PetscInt pseudot = 0; pseudot < 1; pseudot++) {
68 // Loop over each block to perform the RK stages
69 for (PetscInt bi = 0; bi < simCtx->block_number; bi++) {
70 for (istage = 0; istage < 4; istage++) {
71 LOG_ALLOW(LOCAL, LOG_DEBUG, " Block %d, RK Stage %d (alpha=%.4f)...\n", bi, istage, alfa[istage]);
72
73 // a. Calculate the Right-Hand Side (RHS) of the momentum equation.
74 ierr = ComputeRHS(&user[bi], user[bi].Rhs); CHKERRQ(ierr);
75
76 // b. Advance Ucont to the next intermediate stage using the RK coefficient.
77 // Ucont_new = Ucont_old + alpha * dt * RHS
78 ierr = VecWAXPY(user[bi].Ucont, alfa[istage] * dt * st, user[bi].Rhs, user[bi].Ucont_o); CHKERRQ(ierr);
79
80 // c. Update local ghost values for the new intermediate Ucont.
81 ierr = DMGlobalToLocalBegin(user[bi].fda, user[bi].Ucont, INSERT_VALUES, user[bi].lUcont); CHKERRQ(ierr);
82 ierr = DMGlobalToLocalEnd(user[bi].fda, user[bi].Ucont, INSERT_VALUES, user[bi].lUcont); CHKERRQ(ierr);
83
84 // d. Re-apply boundary conditions for the new intermediate velocity.
85 // This is crucial for the stability and accuracy of the multi-stage scheme.
86 ierr = InflowFlux(&user[bi]); CHKERRQ(ierr);
87 ierr = OutflowFlux(&user[bi]); CHKERRQ(ierr);
88 ierr = FormBCS(&user[bi]); CHKERRQ(ierr);
89
90 } // End of RK stages for one block
91
92 /*
93 // Final IBM Interpolation for the block (if enabled)
94 if (simCtx->immersed) {
95 LOG_ALLOW(LOCAL, LOG_DEBUG, " Performing post-RK IBM interpolation for block %d.\n", bi);
96 for (PetscInt ibi = 0; ibi < simCtx->NumberOfBodies; ibi++) {
97 ierr = ibm_interpolation_advanced(&user[bi], &ibm[ibi], ibi, 1); CHKERRQ(ierr);
98 }
99 }
100 */
101
102 } // End loop over blocks
103
104 // --- 3. Inter-Block Communication (Legacy Logic) ---
105 // This is called after all blocks have completed their RK stages.
106 if (simCtx->block_number > 1) {
107 // LOG_ALLOW(GLOBAL, LOG_DEBUG, "Updating multi-block interfaces after RK stages.\n");
108 // ierr = Block_Interface_U(user); CHKERRQ(ierr);
109 }
110
111 } // End of pseudo-time loop
112
113 // --- 4. Cleanup ---
114 // Destroy the RHS vectors that were created at the start of this function.
115 for (PetscInt bi = 0; bi < simCtx->block_number; bi++) {
116 ierr = VecDestroy(&user[bi].Rhs); CHKERRQ(ierr);
117 }
118
119 LOG_ALLOW(GLOBAL, LOG_INFO, "Runge-Kutta solve completed for all blocks.\n");
120
122
123 PetscFunctionReturn(0);
124}
125
126
127
128
129#undef __FUNCT__
130#define __FUNCT__ "ImpRK"
131/**
132 * @brief Advances the momentum equations using an iterative explicit Runge-Kutta scheme.
133 *
134 * This function uses a pseudo-time-stepping approach. It repeatedly applies an
135 * explicit RK scheme within a while loop until the change between iterations
136 * converges below specified tolerances (imp_atol, imp_rtol) or a maximum
137 * number of iterations (imp_MAX_IT) is reached.
138 *
139 * This refactored version keeps the original internal loop over all blocks. All
140 * former global variables are accessed via the SimCtx. It directly calculates the
141 * full RHS (spatial + temporal terms) instead of calling a separate CalcRHS function.
142 *
143 * @param user The array of UserCtx structs for all blocks.
144 * @param ibm (Optional) Pointer to the full array of IBM data structures. Pass NULL if disabled.
145 * @param fsi (Optional) Pointer to the full array of FSI data structures. Pass NULL if disabled.
146 * @return PetscErrorCode 0 on success.
147 */
148PetscErrorCode ImpRK(UserCtx *user, IBMNodes *ibm, FSInfo *fsi)
149{
150 // --- CONTEXT ACQUISITION BLOCK ---
151 SimCtx *simCtx = user[0].simCtx;
152 const PetscInt block_number = simCtx->block_number;
153 const PetscInt immersed = simCtx->immersed;
154 const PetscInt ti = simCtx->step;
155 const PetscInt tistart = simCtx->StartStep;
156 const PetscInt imp_MAX_IT = simCtx->imp_MAX_IT;
157 const PetscReal imp_atol = simCtx->imp_atol;
158 const PetscReal imp_rtol = simCtx->imp_rtol;
159 const PetscReal dt = simCtx->dt;
160 const PetscReal st = simCtx->st;
161 const PetscReal cfl = simCtx->cfl;
162 // --- END CONTEXT ACQUISITION BLOCK ---
163
164 PetscErrorCode ierr;
165 PetscMPIInt rank;
166 PetscInt istage, pseudot, ibi;
167 PetscReal alfa[] = {0.25, 1.0/3.0, 0.5, 1.0};
168 PetscReal ts, te, cput;
169
170 // --- Convergence Metrics ---
171 PetscReal normdU = 10.0, reldU = 1.0, relF = 1.0;
172 PetscReal *normdU0_bk, *normdU1_bk, *normdU_bk, *reldU_bk;
173 PetscReal *normF0_bk, *normF1_bk, *normF_bk, *relF_bk, *lambda;
174
175 PetscFunctionBeginUser;
177 ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank); CHKERRQ(ierr);
178 LOG_ALLOW(GLOBAL, LOG_INFO, "Executing implicit-like RK solver (ImpRK) for %d block(s).\n", block_number);
179
180 // --- Allocate metric arrays ---
181 ierr = PetscMalloc2(block_number, &normdU0_bk, block_number, &normdU1_bk); CHKERRQ(ierr);
182 ierr = PetscMalloc2(block_number, &normdU_bk, block_number, &reldU_bk); CHKERRQ(ierr);
183 ierr = PetscMalloc2(block_number, &normF0_bk, block_number, &normF1_bk); CHKERRQ(ierr);
184 ierr = PetscMalloc2(block_number, &normF_bk, block_number, &lambda); CHKERRQ(ierr);
185 ierr = PetscMalloc1(block_number, &relF_bk); CHKERRQ(ierr);
186
187 ierr = PetscTime(&ts); CHKERRQ(ierr);
188
189 // --- 1. Pre-Loop Initialization (Legacy Logic) ---
190 if (block_number > 1) {
191 // ierr = Block_Interface_U(user); CHKERRQ(ierr);
192 }
193
194 for (PetscInt bi = 0; bi < block_number; bi++) {
195 // Prepare BCs and workspace vectors
196 ierr = InflowFlux(&user[bi]); CHKERRQ(ierr);
197 ierr = OutflowFlux(&user[bi]); CHKERRQ(ierr);
198 ierr = FormBCS(&user[bi]); CHKERRQ(ierr);
199
200 /*
201 if (immersed) {
202 for (ibi = 0; ibi < simCtx->NumberOfBodies; ibi++) {
203 ierr = ibm_interpolation_advanced(&user[bi], &ibm[ibi], ibi, 1); CHKERRQ(ierr);
204 }
205 }
206 */
207
208 ierr = VecDuplicate(user[bi].Ucont, &user[bi].Rhs); CHKERRQ(ierr);
209 ierr = VecDuplicate(user[bi].Ucont, &user[bi].dUcont); CHKERRQ(ierr);
210 ierr = VecDuplicate(user[bi].Ucont, &user[bi].pUcont); CHKERRQ(ierr);
211 ierr = VecCopy(user[bi].Ucont, user[bi].pUcont); CHKERRQ(ierr);
212
213
214 LOG_ALLOW(GLOBAL,LOG_DEBUG," BCs and workspace vectors prepared for Initial RHS calculation.\n");
215
216 // --- Calculate INITIAL RHS for convergence check ---
217 ierr = ComputeRHS(&user[bi], user[bi].Rhs); CHKERRQ(ierr);
218
219 LOG_ALLOW(GLOBAL,LOG_DEBUG, " Initial RHS calculated for convergence check .\n");
220
221 // Add time derivative part
222 if (COEF_TIME_ACCURACY > 1.1 && ti != tistart && ti != 1) {
223 ierr = VecAXPY(user[bi].Rhs, -COEF_TIME_ACCURACY/dt, user[bi].Ucont); CHKERRQ(ierr);
224 ierr = VecAXPY(user[bi].Rhs, +2.0/dt, user[bi].Ucont_o); CHKERRQ(ierr);
225 ierr = VecAXPY(user[bi].Rhs, -0.5/dt, user[bi].Ucont_rm1); CHKERRQ(ierr);
226 } else {
227 ierr = VecAXPY(user[bi].Rhs, -1.0/dt, user[bi].Ucont); CHKERRQ(ierr);
228 ierr = VecAXPY(user[bi].Rhs, +1.0/dt, user[bi].Ucont_o); CHKERRQ(ierr);
229 }
230
231 LOG_ALLOW(GLOBAL,LOG_DEBUG," Time derivative part added to RHS.\n");
232 ierr = VecNorm(user[bi].Rhs, NORM_INFINITY, &normF0_bk[bi]); CHKERRQ(ierr);
233 normF1_bk[bi] = 1000.0; // Initialize for line search logic
234 normdU1_bk[bi] = 1000.0;
235 lambda[bi] = cfl; // Initialize adaptive step size
236
237 LOG_ALLOW(GLOBAL,LOG_INFO," Block %d | Max RHS = %.6f | CFL = %.4f .\n", bi,normF0_bk[bi],cfl);
238 }
239
240 // --- 2. Main Pseudo-Time Iteration Loop ---
241 pseudot = 0;
242 while (((normdU > imp_atol && reldU > imp_rtol) || pseudot < 1) && pseudot < imp_MAX_IT) {
243 pseudot++;
244
245 for (PetscInt bi = 0; bi < block_number; bi++) {
246 for (istage = 0; istage < 4; istage++) {
247 // Advance in time using RK scheme with adaptive step size `lambda`
248 LOG_ALLOW(GLOBAL,LOG_TRACE," Pseudo-Timestep Solver | RK-Stage : %d | Pseudo-Timestep :%d .\n",istage,pseudot);
249 ierr = VecWAXPY(user[bi].Ucont, lambda[bi] * alfa[istage] * dt * st, user[bi].Rhs, user[bi].pUcont); CHKERRQ(ierr);
250
251 LOG_ALLOW(GLOBAL,LOG_VERBOSE, " Ucont updated as Ucont = [lambda(%.4f)]x[alfa(%.4f)]x[dt(%.4f)]x[st(%.4f)]xRhs + pUcont.\n",lambda[bi],alfa[istage],dt,st);
252
253 // Sync ghosts and re-apply BCs for the intermediate stage
254 ierr = DMGlobalToLocalBegin(user[bi].fda, user[bi].Ucont, INSERT_VALUES, user[bi].lUcont); CHKERRQ(ierr);
255 ierr = DMGlobalToLocalEnd(user[bi].fda, user[bi].Ucont, INSERT_VALUES, user[bi].lUcont); CHKERRQ(ierr);
256 ierr = InflowFlux(&user[bi]); CHKERRQ(ierr);
257 ierr = OutflowFlux(&user[bi]); CHKERRQ(ierr);
258 ierr = FormBCS(&user[bi]); CHKERRQ(ierr);
259
260 LOG_ALLOW(GLOBAL,LOG_TRACE, " Ghosts synced and BCs reapplied.\n");
261
262 // --- Re-calculate the full RHS for the next RK stage ---
263 ierr = ComputeRHS(&user[bi], user[bi].Rhs); CHKERRQ(ierr);
264
265 LOG_ALLOW(GLOBAL,LOG_TRACE, " RHS calculated for the next RK stage.\n");
266
267 if (COEF_TIME_ACCURACY > 1.1 && ti != tistart && ti != 1) {
268 ierr = VecAXPY(user[bi].Rhs, -COEF_TIME_ACCURACY/dt, user[bi].Ucont); CHKERRQ(ierr);
269 ierr = VecAXPY(user[bi].Rhs, +2.0/dt, user[bi].Ucont_o); CHKERRQ(ierr);
270 ierr = VecAXPY(user[bi].Rhs, -0.5/dt, user[bi].Ucont_rm1); CHKERRQ(ierr);
271 } else {
272 ierr = VecAXPY(user[bi].Rhs, -1.0/dt, user[bi].Ucont); CHKERRQ(ierr);
273 ierr = VecAXPY(user[bi].Rhs, +1.0/dt, user[bi].Ucont_o); CHKERRQ(ierr);
274 }
275
276 if(istage <3) LOG_ALLOW(GLOBAL,LOG_TRACE," Time derivative part updated for next RK stage .\n");
277 else LOG_ALLOW(GLOBAL,LOG_DEBUG," All RK stages complete.\n");
278 } // End RK stages
279
280 /*
281 if (immersed) {
282 for (ibi = 0; ibi < simCtx->NumberOfBodies; ibi++) {
283 ierr = ibm_interpolation_advanced(&user[bi], &ibm[ibi], ibi, 1); CHKERRQ(ierr);
284 }
285 }
286 */
287
288 // --- Calculate Convergence Metrics for this block ---
289
290 ierr = VecWAXPY(user[bi].dUcont, -1.0, user[bi].Ucont, user[bi].pUcont); CHKERRQ(ierr);
291
292 LOG_ALLOW(GLOBAL,LOG_TRACE," Block %d | dU calculated .\n",bi);
293
294 normdU1_bk[bi] = normdU_bk[bi];
295
296 ierr = VecNorm(user[bi].dUcont, NORM_INFINITY, &normdU_bk[bi]); CHKERRQ(ierr);
297 ierr = VecNorm(user[bi].Rhs, NORM_INFINITY, &normF_bk[bi]); CHKERRQ(ierr);
298
299 LOG_ALLOW(GLOBAL,LOG_TRACE," Block %d | |dU| = %le | |U| = %le.\n",bi,normdU_bk[bi],normF_bk[bi]);
300
301 if (pseudot == 1) {
302 normdU0_bk[bi] = normdU_bk[bi];
303 reldU_bk[bi] = 1.0;
304 normF0_bk[bi] = normF_bk[bi];
305 relF_bk[bi] = 1.0;
306 } else { //pseudot > 1
307 if (normdU0_bk[bi] > 1.0e-10) {
308 reldU_bk[bi] = normdU_bk[bi] / normdU0_bk[bi];
309 }
310 else {
311 reldU_bk[bi] = 0.0;
312 }
313 if(normF0_bk[bi] > 1.0e-10) {
314 relF_bk[bi] = normF_bk[bi] / normF0_bk[bi];
315 }else {
316 relF_bk[bi] = 0.0;
317 }
318 }
319
320
321 LOG_ALLOW(GLOBAL,LOG_TRACE, " Block %d | |dU|/|dU_0| = %le. \n",bi,reldU_bk[bi]);
322
323 LOG_ALLOW(GLOBAL,LOG_TRACE, " Block %d | |R|/|R_0| = %le. \n",bi,relF_bk[bi]);
324
325 // file logging
326 if (!rank) {
327 FILE *f;
328 char filen[80];
329 sprintf(filen, "%s/Momentum_Solver_Convergence_History_Block_%1d.log", simCtx->log_dir,bi);
330 f = fopen(filen, "a");
331 PetscFPrintf(PETSC_COMM_WORLD, f, "Block %d | Step: %d | Pseudo-Timestep %d | |dU| %le | |dU|/|dU_0| %le | |U| = %le \n",bi,(int)ti, (int)pseudot, normdU_bk[bi], reldU_bk[bi], normF_bk[bi]);
332 fclose(f);
333 }
334
335 LOG_ALLOW(GLOBAL,LOG_DEBUG, " Block %d | Step: %d | Pseudo-Timestep %d | |dU| %le | |dU|/|dU_0| %le | |U| = %le \n",bi,ti,pseudot,normdU_bk[bi], reldU_bk[bi], normF_bk[bi]);
336 } // End loop over blocks
337
338 // --- Update Global Convergence Criteria and Perform Line Search ---
339 normdU = -1.0e20; reldU = -1.0e20; relF = -1.0e20;
340 for (PetscInt bi = 0; bi < block_number; bi++) {
341 normdU = PetscMax(normdU_bk[bi], normdU);
342 reldU = PetscMax(reldU_bk[bi], reldU);
343 relF = PetscMax(relF_bk[bi], relF);
344 }
345
346 ierr = PetscTime(&te); CHKERRQ(ierr);
347 cput = te - ts;
348 LOG_ALLOW(GLOBAL, LOG_INFO, " Iter(Pseudo-Timestep) %d: |dU|=%e, |dU|/|dU_0|=%e, |U|=%e, CPU=%.2fs\n",
349 pseudot, normdU, reldU, relF, cput);
350
351 for (PetscInt bi = 0; bi < block_number; bi++) {
352 // Adaptive step size logic (line search)
353 if (pseudot > 1 && normF_bk[bi] > normF1_bk[bi] && normdU_bk[bi] > normdU1_bk[bi] && lambda[bi] > 0.01) {
354 lambda[bi] *= 0.5;
355 pseudot--; // Retry this iteration with a smaller step
356 LOG_ALLOW(LOCAL, LOG_DEBUG, " Block %d: Line search failed. Reducing lambda to %f and retrying iter.\n", bi, lambda[bi]);
357 ierr = VecCopy(user[bi].pUcont, user[bi].Ucont); CHKERRQ(ierr);
358 ierr = DMGlobalToLocalBegin(user[bi].fda, user[bi].Ucont, INSERT_VALUES, user[bi].lUcont); CHKERRQ(ierr);
359 ierr = DMGlobalToLocalEnd(user[bi].fda, user[bi].Ucont, INSERT_VALUES, user[bi].lUcont); CHKERRQ(ierr);
360 } else {
361 ierr = VecCopy(user[bi].Ucont, user[bi].pUcont); CHKERRQ(ierr);
362 normF1_bk[bi] = normF_bk[bi];
363 normdU1_bk[bi] = normdU_bk[bi];
364 }
365 }
366
367 if (block_number > 1) {
368 LOG_ALLOW(GLOBAL, LOG_DEBUG, " Updating multi-block interfaces after iteration %d.\n", pseudot);
369 // ierr = Block_Interface_U(user); CHKERRQ(ierr);
370 }
371 } // End while loop
372
373 // --- Final Cleanup ---
374 /*
375 if (block_number > 1) {
376 Block_Interface_Correction(user);
377 if (simCtx->blank) {
378 Block_Blank_Correction_adv(&user[0], 0);
379 }
380 }
381 */
382 for (PetscInt bi = 0; bi < block_number; bi++) {
383 ierr = VecDestroy(&user[bi].Rhs); CHKERRQ(ierr);
384 ierr = VecDestroy(&user[bi].dUcont); CHKERRQ(ierr);
385 ierr = VecDestroy(&user[bi].pUcont); CHKERRQ(ierr);
386 }
387
388 ierr = PetscFree2(normdU0_bk, normdU1_bk);CHKERRQ(ierr);
389 ierr = PetscFree2(normdU_bk, reldU_bk);CHKERRQ(ierr);
390 ierr = PetscFree2(normF0_bk, normF1_bk);CHKERRQ(ierr);
391 ierr = PetscFree2(normF_bk, lambda); CHKERRQ(ierr);
392 ierr = PetscFree(relF_bk); CHKERRQ(ierr);
393
394 LOG_ALLOW(GLOBAL, LOG_INFO, "ImpRK solve completed after %d iterations.\n", pseudot);
396 PetscFunctionReturn(0);
397}
PetscErrorCode OutflowFlux(UserCtx *user)
Calculates the total outgoing flux through all OUTLET faces for reporting.
Definition Boundaries.c:932
PetscErrorCode InflowFlux(UserCtx *user)
Applies inlet boundary conditions based on the modern BC handling system.
Definition Boundaries.c:706
PetscErrorCode FormBCS(UserCtx *user)
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 for ALL blocks by one time step using an explicit 4-stage,...
#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
@ LOG_TRACE
Very fine-grained tracing information for in-depth debugging.
Definition logging.h:34
@ LOG_INFO
Informational messages about program execution.
Definition logging.h:32
@ LOG_DEBUG
Detailed debugging information.
Definition logging.h:33
@ LOG_VERBOSE
Extremely detailed logs, typically for development use only.
Definition logging.h:35
#define PROFILE_FUNCTION_BEGIN
Marks the beginning of a profiled code block (typically a function).
Definition logging.h:731
PetscErrorCode ComputeRHS(UserCtx *user, Vec Rhs)
Computes the Right-Hand Side (RHS) of the momentum equations.
Definition rhs.c:1159
PetscInt block_number
Definition variables.h:593
SimCtx * simCtx
Back-pointer to the master simulation context.
Definition variables.h:664
PetscReal imp_rtol
Definition variables.h:576
PetscReal st
Definition variables.h:581
PetscReal imp_atol
Definition variables.h:576
PetscReal dt
Definition variables.h:552
PetscInt StartStep
Definition variables.h:548
char log_dir[PETSC_MAX_PATH_LEN]
Definition variables.h:562
PetscReal cfl
Definition variables.h:581
PetscInt step
Definition variables.h:546
PetscInt imp_MAX_IT
Definition variables.h:575
PetscInt immersed
Definition variables.h:567
#define COEF_TIME_ACCURACY
Coefficient controlling the temporal accuracy scheme (e.g., 1.5 for BDF2).
Definition variables.h:57
Holds all data related to the state and motion of a body in FSI.
Definition variables.h:382
Represents a collection of nodes forming a surface for the IBM.
Definition variables.h:309
The master context for the entire simulation.
Definition variables.h:538
User-defined context containing data specific to a single computational grid level.
Definition variables.h:661