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

Go to the source code of this file.

Macros

#define __FUNCT__   "RungeKutta"
 
#define __FUNCT__   "ImpRK"
 

Functions

PetscErrorCode RungeKutta (UserCtx *user, IBMNodes *ibm, FSInfo *fsi)
 Advances the momentum equations for ALL blocks by one time step using an explicit 4-stage, 4th-order Runge-Kutta scheme.
 
PetscErrorCode ImpRK (UserCtx *user, IBMNodes *ibm, FSInfo *fsi)
 Advances the momentum equations using an iterative explicit Runge-Kutta scheme.
 

Macro Definition Documentation

◆ __FUNCT__ [1/2]

#define __FUNCT__   "RungeKutta"

Definition at line 4 of file implicitsolvers.c.

◆ __FUNCT__ [2/2]

#define __FUNCT__   "ImpRK"

Definition at line 4 of file implicitsolvers.c.

Function Documentation

◆ RungeKutta()

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

Advances the momentum equations for ALL blocks by one time step using an explicit 4-stage, 4th-order Runge-Kutta scheme.

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:

◆ ImpRK()

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

Advances the momentum equations using an iterative explicit Runge-Kutta scheme.

This function uses a pseudo-time-stepping approach. It repeatedly applies an explicit RK scheme within a while loop until the change between iterations converges below specified tolerances (imp_atol, imp_rtol) or a maximum number of iterations (imp_MAX_IT) is reached.

This refactored version keeps the original internal loop over all blocks. All former global variables are accessed via the SimCtx. It directly calculates the full RHS (spatial + temporal terms) instead of calling a separate CalcRHS function.

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 142 of file implicitsolvers.c.

143{
144 // --- CONTEXT ACQUISITION BLOCK ---
145 SimCtx *simCtx = user[0].simCtx;
146 const PetscInt block_number = simCtx->block_number;
147 const PetscInt immersed = simCtx->immersed;
148 const PetscInt ti = simCtx->step;
149 const PetscInt tistart = simCtx->StartStep;
150 const PetscInt imp_MAX_IT = simCtx->imp_MAX_IT;
151 const PetscReal imp_atol = simCtx->imp_atol;
152 const PetscReal imp_rtol = simCtx->imp_rtol;
153 const PetscReal dt = simCtx->dt;
154 const PetscReal st = simCtx->st;
155 const PetscReal cfl = simCtx->cfl;
156 // --- END CONTEXT ACQUISITION BLOCK ---
157
158 PetscErrorCode ierr;
159 PetscMPIInt rank;
160 PetscInt istage, pseudot, ibi;
161 PetscReal alfa[] = {0.25, 1.0/3.0, 0.5, 1.0};
162 PetscReal ts, te, cput;
163
164 // --- Convergence Metrics ---
165 PetscReal normdU = 10.0, reldU = 1.0, relF = 1.0;
166 PetscReal *normdU0_bk, *normdU1_bk, *normdU_bk, *reldU_bk;
167 PetscReal *normF0_bk, *normF1_bk, *normF_bk, *relF_bk, *lambda;
168
169 PetscFunctionBeginUser;
171 ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank); CHKERRQ(ierr);
172 LOG_ALLOW(GLOBAL, LOG_INFO, "Executing implicit-like RK solver (ImpRK) for %d block(s).\n", block_number);
173
174 // --- Allocate metric arrays ---
175 ierr = PetscMalloc2(block_number, &normdU0_bk, block_number, &normdU1_bk); CHKERRQ(ierr);
176 ierr = PetscMalloc2(block_number, &normdU_bk, block_number, &reldU_bk); CHKERRQ(ierr);
177 ierr = PetscMalloc2(block_number, &normF0_bk, block_number, &normF1_bk); CHKERRQ(ierr);
178 ierr = PetscMalloc2(block_number, &normF_bk, block_number, &lambda); CHKERRQ(ierr);
179 ierr = PetscMalloc1(block_number, &relF_bk); CHKERRQ(ierr);
180
181 ierr = PetscTime(&ts); CHKERRQ(ierr);
182
183 // --- 1. Pre-Loop Initialization (Legacy Logic) ---
184 if (block_number > 1) {
185 // ierr = Block_Interface_U(user); CHKERRQ(ierr);
186 }
187
188 for (PetscInt bi = 0; bi < block_number; bi++) {
189 // Prepare BCs and workspace vectors
190 ierr = InflowFlux(&user[bi]); CHKERRQ(ierr);
191 ierr = OutflowFlux(&user[bi]); CHKERRQ(ierr);
192 ierr = FormBCS(&user[bi]); CHKERRQ(ierr);
193
194 /*
195 if (immersed) {
196 for (ibi = 0; ibi < simCtx->NumberOfBodies; ibi++) {
197 ierr = ibm_interpolation_advanced(&user[bi], &ibm[ibi], ibi, 1); CHKERRQ(ierr);
198 }
199 }
200 */
201
202 ierr = VecDuplicate(user[bi].Ucont, &user[bi].Rhs); CHKERRQ(ierr);
203 ierr = VecDuplicate(user[bi].Ucont, &user[bi].dUcont); CHKERRQ(ierr);
204 ierr = VecDuplicate(user[bi].Ucont, &user[bi].pUcont); CHKERRQ(ierr);
205 ierr = VecCopy(user[bi].Ucont, user[bi].pUcont); CHKERRQ(ierr);
206
207
208 LOG_ALLOW(GLOBAL,LOG_DEBUG," BCs and workspace vectors prepared for Initial RHS calculation.\n");
209
210 // --- Calculate INITIAL RHS for convergence check ---
211 ierr = FormFunction1(&user[bi], user[bi].Rhs); CHKERRQ(ierr);
212
213 LOG_ALLOW(GLOBAL,LOG_DEBUG, " Initial RHS calculated for convergence check .\n");
214
215 // Add time derivative part
216 if (COEF_TIME_ACCURACY > 1.1 && ti != tistart && ti != 1) {
217 ierr = VecAXPY(user[bi].Rhs, -COEF_TIME_ACCURACY/dt, user[bi].Ucont); CHKERRQ(ierr);
218 ierr = VecAXPY(user[bi].Rhs, +2.0/dt, user[bi].Ucont_o); CHKERRQ(ierr);
219 ierr = VecAXPY(user[bi].Rhs, -0.5/dt, user[bi].Ucont_rm1); CHKERRQ(ierr);
220 } else {
221 ierr = VecAXPY(user[bi].Rhs, -1.0/dt, user[bi].Ucont); CHKERRQ(ierr);
222 ierr = VecAXPY(user[bi].Rhs, +1.0/dt, user[bi].Ucont_o); CHKERRQ(ierr);
223 }
224
225 LOG_ALLOW(GLOBAL,LOG_DEBUG," Time derivative part added to RHS.\n");
226 ierr = VecNorm(user[bi].Rhs, NORM_INFINITY, &normF0_bk[bi]); CHKERRQ(ierr);
227 normF1_bk[bi] = 1000.0; // Initialize for line search logic
228 normdU1_bk[bi] = 1000.0;
229 lambda[bi] = cfl; // Initialize adaptive step size
230
231 LOG_ALLOW(GLOBAL,LOG_INFO," Block %d | Max RHS = %.6f | CFL = %.4f .\n", bi,normF0_bk[bi],cfl);
232 }
233
234 // --- 2. Main Pseudo-Time Iteration Loop ---
235 pseudot = 0;
236 while (((normdU > imp_atol && reldU > imp_rtol) || pseudot < 1) && pseudot < imp_MAX_IT) {
237 pseudot++;
238
239 for (PetscInt bi = 0; bi < block_number; bi++) {
240 for (istage = 0; istage < 4; istage++) {
241 // Advance in time using RK scheme with adaptive step size `lambda`
242 LOG_ALLOW(GLOBAL,LOG_INFO," Pseudo-Timestep Solver | RK-Stage : %d | Pseudo-Timestep :%d .\n",istage,pseudot);
243 ierr = VecWAXPY(user[bi].Ucont, lambda[bi] * alfa[istage] * dt * st, user[bi].Rhs, user[bi].pUcont); CHKERRQ(ierr);
244
245 LOG_ALLOW(GLOBAL,LOG_INFO, " Ucont updated as Ucont = [lambda(%.4f)]x[alfa(%.4f)]x[dt(%.4f)]x[st(%.4f)]xRhs + pUcont.\n",lambda[bi],alfa[istage],dt,st);
246
247 // Sync ghosts and re-apply BCs for the intermediate stage
248 ierr = DMGlobalToLocalBegin(user[bi].fda, user[bi].Ucont, INSERT_VALUES, user[bi].lUcont); CHKERRQ(ierr);
249 ierr = DMGlobalToLocalEnd(user[bi].fda, user[bi].Ucont, INSERT_VALUES, user[bi].lUcont); CHKERRQ(ierr);
250 ierr = InflowFlux(&user[bi]); CHKERRQ(ierr);
251 ierr = OutflowFlux(&user[bi]); CHKERRQ(ierr);
252 ierr = FormBCS(&user[bi]); CHKERRQ(ierr);
253
254 LOG_ALLOW(GLOBAL,LOG_INFO, " Ghosts synced and BCs reapplied.\n");
255
256 // --- Re-calculate the full RHS for the next RK stage ---
257 ierr = FormFunction1(&user[bi], user[bi].Rhs); CHKERRQ(ierr);
258
259 LOG_ALLOW(GLOBAL,LOG_INFO, " RHS calculated for the next RK stage.\n");
260
261 if (COEF_TIME_ACCURACY > 1.1 && ti != tistart && ti != 1) {
262 ierr = VecAXPY(user[bi].Rhs, -COEF_TIME_ACCURACY/dt, user[bi].Ucont); CHKERRQ(ierr);
263 ierr = VecAXPY(user[bi].Rhs, +2.0/dt, user[bi].Ucont_o); CHKERRQ(ierr);
264 ierr = VecAXPY(user[bi].Rhs, -0.5/dt, user[bi].Ucont_rm1); CHKERRQ(ierr);
265 } else {
266 ierr = VecAXPY(user[bi].Rhs, -1.0/dt, user[bi].Ucont); CHKERRQ(ierr);
267 ierr = VecAXPY(user[bi].Rhs, +1.0/dt, user[bi].Ucont_o); CHKERRQ(ierr);
268 }
269
270 if(istage <3) LOG_ALLOW(GLOBAL,LOG_INFO," Time derivative part updated for next RK stage .\n");
271 else LOG_ALLOW(GLOBAL,LOG_INFO," All RK stages complete.\n");
272 } // End RK stages
273
274 /*
275 if (immersed) {
276 for (ibi = 0; ibi < simCtx->NumberOfBodies; ibi++) {
277 ierr = ibm_interpolation_advanced(&user[bi], &ibm[ibi], ibi, 1); CHKERRQ(ierr);
278 }
279 }
280 */
281
282 // --- Calculate Convergence Metrics for this block ---
283
284 ierr = VecWAXPY(user[bi].dUcont, -1.0, user[bi].Ucont, user[bi].pUcont); CHKERRQ(ierr);
285
286 LOG_ALLOW(GLOBAL,LOG_DEBUG," Block %d | dU calculated .\n",bi);
287
288 normdU1_bk[bi] = normdU_bk[bi];
289
290 ierr = VecNorm(user[bi].dUcont, NORM_INFINITY, &normdU_bk[bi]); CHKERRQ(ierr);
291 ierr = VecNorm(user[bi].Rhs, NORM_INFINITY, &normF_bk[bi]); CHKERRQ(ierr);
292
293 LOG_ALLOW(GLOBAL,LOG_DEBUG," Block %d | |dU| = %le | |U| = %le.\n",bi,normdU_bk[bi],normF_bk[bi]);
294
295 if (pseudot == 1) {
296 normdU0_bk[bi] = normdU_bk[bi];
297 reldU_bk[bi] = 1.0;
298 normF0_bk[bi] = normF_bk[bi];
299 relF_bk[bi] = 1.0;
300 } else { //pseudot > 1
301 if (normdU0_bk[bi] > 1.0e-10) {
302 reldU_bk[bi] = normdU_bk[bi] / normdU0_bk[bi];
303 }
304 else {
305 reldU_bk[bi] = 0.0;
306 }
307 if(normF0_bk[bi] > 1.0e-10) {
308 relF_bk[bi] = normF_bk[bi] / normF0_bk[bi];
309 }else {
310 relF_bk[bi] = 0.0;
311 }
312 }
313
314
315 LOG_ALLOW(GLOBAL,LOG_DEBUG, " Block %d | |dU|/|dU_0| = %le. \n",bi,reldU_bk[bi]);
316
317 LOG_ALLOW(GLOBAL,LOG_DEBUG, " Block %d | |R|/|R_0| = %le. \n",bi,relF_bk[bi]);
318
319 // file logging
320 if (!rank) {
321 FILE *f;
322 char filen[80];
323 sprintf(filen, "logs/Momentum_Solver_Convergence_History_Block_%1d.log", bi);
324 f = fopen(filen, "a");
325 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]);
326 fclose(f);
327 }
328
329 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]);
330 } // End loop over blocks
331
332 // --- Update Global Convergence Criteria and Perform Line Search ---
333 normdU = -1.0e20; reldU = -1.0e20; relF = -1.0e20;
334 for (PetscInt bi = 0; bi < block_number; bi++) {
335 normdU = PetscMax(normdU_bk[bi], normdU);
336 reldU = PetscMax(reldU_bk[bi], reldU);
337 relF = PetscMax(relF_bk[bi], relF);
338 }
339
340 ierr = PetscTime(&te); CHKERRQ(ierr);
341 cput = te - ts;
342 LOG_ALLOW(GLOBAL, LOG_INFO, " Iter(Pseudo-Timestep) %d: |dU|=%e, |dU|/|dU_0|=%e, |U|=%e, CPU=%.2fs\n",
343 pseudot, normdU, reldU, relF, cput);
344
345 for (PetscInt bi = 0; bi < block_number; bi++) {
346 // Adaptive step size logic (line search)
347 if (pseudot > 1 && normF_bk[bi] > normF1_bk[bi] && normdU_bk[bi] > normdU1_bk[bi] && lambda[bi] > 0.01) {
348 lambda[bi] *= 0.5;
349 pseudot--; // Retry this iteration with a smaller step
350 LOG_ALLOW(LOCAL, LOG_DEBUG, " Block %d: Line search failed. Reducing lambda to %f and retrying iter.\n", bi, lambda[bi]);
351 ierr = VecCopy(user[bi].pUcont, user[bi].Ucont); CHKERRQ(ierr);
352 ierr = DMGlobalToLocalBegin(user[bi].fda, user[bi].Ucont, INSERT_VALUES, user[bi].lUcont); CHKERRQ(ierr);
353 ierr = DMGlobalToLocalEnd(user[bi].fda, user[bi].Ucont, INSERT_VALUES, user[bi].lUcont); CHKERRQ(ierr);
354 } else {
355 ierr = VecCopy(user[bi].Ucont, user[bi].pUcont); CHKERRQ(ierr);
356 normF1_bk[bi] = normF_bk[bi];
357 normdU1_bk[bi] = normdU_bk[bi];
358 }
359 }
360
361 if (block_number > 1) {
362 LOG_ALLOW(GLOBAL, LOG_DEBUG, " Updating multi-block interfaces after iteration %d.\n", pseudot);
363 // ierr = Block_Interface_U(user); CHKERRQ(ierr);
364 }
365 } // End while loop
366
367 // --- Final Cleanup ---
368 /*
369 if (block_number > 1) {
370 Block_Interface_Correction(user);
371 if (simCtx->blank) {
372 Block_Blank_Correction_adv(&user[0], 0);
373 }
374 }
375 */
376 for (PetscInt bi = 0; bi < block_number; bi++) {
377 ierr = VecDestroy(&user[bi].Rhs); CHKERRQ(ierr);
378 ierr = VecDestroy(&user[bi].dUcont); CHKERRQ(ierr);
379 ierr = VecDestroy(&user[bi].pUcont); CHKERRQ(ierr);
380 }
381
382 ierr = PetscFree2(normdU0_bk, normdU1_bk);CHKERRQ(ierr);
383 ierr = PetscFree2(normdU_bk, reldU_bk);CHKERRQ(ierr);
384 ierr = PetscFree2(normF0_bk, normF1_bk);CHKERRQ(ierr);
385 ierr = PetscFree2(normF_bk, lambda); CHKERRQ(ierr);
386 ierr = PetscFree(relF_bk); CHKERRQ(ierr);
387
388 LOG_ALLOW(GLOBAL, LOG_INFO, "ImpRK solve completed after %d iterations.\n", pseudot);
390 PetscFunctionReturn(0);
391}
#define PROFILE_FUNCTION_END
Marks the end of a profiled code block.
Definition logging.h:767
#define PROFILE_FUNCTION_BEGIN
Marks the beginning of a profiled code block (typically a function).
Definition logging.h:758
PetscReal imp_rtol
Definition variables.h:544
PetscReal imp_atol
Definition variables.h:544
PetscInt StartStep
Definition variables.h:523
PetscReal cfl
Definition variables.h:549
PetscInt step
Definition variables.h:521
PetscInt imp_MAX_IT
Definition variables.h:543
PetscInt immersed
Definition variables.h:535
#define COEF_TIME_ACCURACY
Coefficient controlling the temporal accuracy scheme (e.g., 1.5 for BDF2).
Definition variables.h:57
Here is the call graph for this function:
Here is the caller graph for this function: