xref: /honee/src/setupts.c (revision da02a6e7947cab8e173fd83bf2b4eda896c06b14)
1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3 
4 /// @file
5 /// Time-stepping functions for Navier-Stokes example using PETSc
6 
7 #include <ceed.h>
8 #include <petscdmplex.h>
9 #include <petscts.h>
10 
11 #include <navierstokes.h>
12 #include "../qfunctions/newtonian_state.h"
13 
14 // @brief Insert Boundary values if it's a new time
15 PetscErrorCode UpdateBoundaryValues(User user, Vec Q_loc, PetscReal t) {
16   PetscFunctionBeginUser;
17   if (user->time_bc_set != t) {
18     PetscCall(DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL));
19     user->time_bc_set = t;
20   }
21   PetscFunctionReturn(PETSC_SUCCESS);
22 }
23 
24 // RHS (Explicit time-stepper) function setup
25 //   This is the RHS of the ODE, given as u_t = G(t,u)
26 //   This function takes in a state vector Q and writes into G
27 PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) {
28   User         user = *(User *)user_data;
29   Ceed         ceed = user->ceed;
30   PetscScalar  dt;
31   Vec          Q_loc = user->Q_loc;
32   PetscMemType q_mem_type;
33 
34   PetscFunctionBeginUser;
35   // Update time dependent data
36   PetscCall(UpdateBoundaryValues(user, Q_loc, t));
37   if (user->phys->solution_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_rhs_ctx->op, user->phys->solution_time_label, &t));
38   PetscCall(TSGetTimeStep(ts, &dt));
39   if (user->phys->timestep_size_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_rhs_ctx->op, user->phys->timestep_size_label, &dt));
40 
41   PetscCall(ApplyCeedOperatorGlobalToGlobal(Q, G, user->op_rhs_ctx));
42 
43   PetscCall(VecReadPetscToCeed(Q_loc, &q_mem_type, user->q_ceed));
44 
45   // Inverse of the mass matrix
46   PetscCall(KSPSolve(user->mass_ksp, G, G));
47 
48   PetscCall(VecReadCeedToPetsc(user->q_ceed, q_mem_type, Q_loc));
49   PetscFunctionReturn(PETSC_SUCCESS);
50 }
51 
52 // Surface forces function setup
53 static PetscErrorCode Surface_Forces_NS(DM dm, Vec G_loc, PetscInt num_walls, const PetscInt walls[], PetscScalar *reaction_force) {
54   DMLabel            face_label;
55   const PetscScalar *g_array;
56   PetscInt           dim  = 3;
57   MPI_Comm           comm = PetscObjectComm((PetscObject)dm);
58   PetscSection       section;
59 
60   PetscFunctionBeginUser;
61   PetscCall(DMGetLabel(dm, "Face Sets", &face_label));
62   PetscCall(VecGetArrayRead(G_loc, &g_array));
63   for (PetscInt w = 0; w < num_walls; w++) {
64     const PetscInt wall = walls[w], *points;
65     IS             wall_is;
66     PetscInt       num_points, num_comp = 0;
67 
68     PetscCall(DMLabelGetStratumIS(face_label, wall, &wall_is));
69     if (!wall_is) continue;  // No wall points on this process, skip
70 
71     PetscCall(DMGetLocalSection(dm, &section));
72     PetscCall(PetscSectionGetFieldComponents(section, 0, &num_comp));
73     PetscCall(ISGetSize(wall_is, &num_points));
74     PetscCall(ISGetIndices(wall_is, &points));
75     for (PetscInt i = 0; i < num_points; i++) {
76       const PetscInt           p = points[i];
77       const StateConservative *r;
78       PetscInt                 dof;
79 
80       PetscCall(DMPlexPointLocalRead(dm, p, g_array, &r));
81       PetscCall(PetscSectionGetDof(section, p, &dof));
82       for (PetscInt node = 0; node < dof / num_comp; node++) {
83         for (PetscInt j = 0; j < dim; j++) {
84           reaction_force[w * dim + j] -= r[node].momentum[j];
85         }
86       }
87     }
88     PetscCall(ISRestoreIndices(wall_is, &points));
89     PetscCall(ISDestroy(&wall_is));
90   }
91   PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, reaction_force, dim * num_walls, MPIU_SCALAR, MPI_SUM, comm));
92   PetscCall(VecRestoreArrayRead(G_loc, &g_array));
93   PetscFunctionReturn(PETSC_SUCCESS);
94 }
95 
96 // Implicit time-stepper function setup
97 PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) {
98   User         user = *(User *)user_data;
99   Ceed         ceed = user->ceed;
100   PetscScalar  dt;
101   Vec          Q_loc = user->Q_loc, Q_dot_loc = user->Q_dot_loc, G_loc;
102   PetscMemType q_mem_type, q_dot_mem_type, g_mem_type;
103 
104   PetscFunctionBeginUser;
105   // Get local vectors
106   PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
107 
108   // Update time dependent data
109   PetscCall(UpdateBoundaryValues(user, Q_loc, t));
110   if (user->phys->solution_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_ifunction, user->phys->solution_time_label, &t));
111   PetscCall(TSGetTimeStep(ts, &dt));
112   if (user->phys->timestep_size_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_ifunction, user->phys->timestep_size_label, &dt));
113 
114   // Global-to-local
115   PetscCall(DMGlobalToLocalBegin(user->dm, Q, INSERT_VALUES, Q_loc));
116   PetscCall(DMGlobalToLocalBegin(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
117   PetscCall(DMGlobalToLocalEnd(user->dm, Q, INSERT_VALUES, Q_loc));
118   if (user->app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) PetscCall(DivDiffFluxProjectionApply(user->diff_flux_proj, Q_loc));
119   PetscCall(DMGlobalToLocalEnd(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
120 
121   // Place PETSc vectors in CEED vectors
122   PetscCall(VecReadPetscToCeed(Q_loc, &q_mem_type, user->q_ceed));
123   PetscCall(VecReadPetscToCeed(Q_dot_loc, &q_dot_mem_type, user->q_dot_ceed));
124   PetscCall(VecPetscToCeed(G_loc, &g_mem_type, user->g_ceed));
125 
126   // Apply CEED operator
127   PetscCall(PetscLogEventBegin(FLUIDS_CeedOperatorApply, Q, G, 0, 0));
128   PetscCall(PetscLogGpuTimeBegin());
129   PetscCallCeed(user->ceed, CeedOperatorApply(user->op_ifunction, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE));
130   PetscCall(PetscLogGpuTimeEnd());
131   PetscCall(PetscLogEventEnd(FLUIDS_CeedOperatorApply, Q, G, 0, 0));
132 
133   // Restore vectors
134   PetscCall(VecReadCeedToPetsc(user->q_ceed, q_mem_type, Q_loc));
135   PetscCall(VecReadCeedToPetsc(user->q_dot_ceed, q_dot_mem_type, Q_dot_loc));
136   PetscCall(VecCeedToPetsc(user->g_ceed, g_mem_type, G_loc));
137 
138   if (user->app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) {
139     PetscCall(SgsDDApplyIFunction(user, Q_loc, G_loc));
140   }
141 
142   // Local-to-Global
143   PetscCall(VecZeroEntries(G));
144   PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G));
145 
146   // Restore vectors
147   PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
148   PetscFunctionReturn(PETSC_SUCCESS);
149 }
150 
151 PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, PetscReal shift, Mat J, Mat J_pre, void *user_data) {
152   User      user = *(User *)user_data;
153   PetscBool J_is_matceed, J_is_mffd, J_pre_is_matceed, J_pre_is_mffd;
154 
155   PetscFunctionBeginUser;
156   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATMFFD, &J_is_mffd));
157   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATCEED, &J_is_matceed));
158   PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATMFFD, &J_pre_is_mffd));
159   PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATCEED, &J_pre_is_matceed));
160 
161   PetscCall(MatCeedSetContextReal(user->mat_ijacobian, "ijacobian time shift", shift));
162 
163   if (J_is_matceed || J_is_mffd) {
164     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
165     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
166   } else PetscCall(MatCeedAssembleCOO(user->mat_ijacobian, J));
167 
168   if (J_pre_is_matceed && J != J_pre) {
169     PetscCall(MatAssemblyBegin(J_pre, MAT_FINAL_ASSEMBLY));
170     PetscCall(MatAssemblyEnd(J_pre, MAT_FINAL_ASSEMBLY));
171   } else if (!J_pre_is_matceed && !J_pre_is_mffd && J != J_pre) {
172     PetscCall(MatCeedAssembleCOO(user->mat_ijacobian, J_pre));
173   }
174   PetscFunctionReturn(PETSC_SUCCESS);
175 }
176 
177 PetscErrorCode WriteOutput(User user, Vec Q, PetscInt step_no, PetscScalar time) {
178   Vec         Q_loc;
179   char        file_path[PETSC_MAX_PATH_LEN];
180   PetscViewer viewer;
181 
182   PetscFunctionBeginUser;
183   if (user->app_ctx->checkpoint_vtk) {
184     // Set up output
185     PetscCall(DMGetLocalVector(user->dm, &Q_loc));
186     PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec"));
187     PetscCall(VecZeroEntries(Q_loc));
188     PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
189 
190     // Output
191     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
192 
193     PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path, FILE_MODE_WRITE, &viewer));
194     PetscCall(VecView(Q_loc, viewer));
195     PetscCall(PetscViewerDestroy(&viewer));
196     if (user->dm_viz) {
197       Vec         Q_refined, Q_refined_loc;
198       char        file_path_refined[PETSC_MAX_PATH_LEN];
199       PetscViewer viewer_refined;
200 
201       PetscCall(DMGetGlobalVector(user->dm_viz, &Q_refined));
202       PetscCall(DMGetLocalVector(user->dm_viz, &Q_refined_loc));
203       PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined"));
204 
205       PetscCall(MatInterpolate(user->interp_viz, Q, Q_refined));
206       PetscCall(VecZeroEntries(Q_refined_loc));
207       PetscCall(DMGlobalToLocal(user->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc));
208 
209       PetscCall(
210           PetscSNPrintf(file_path_refined, sizeof file_path_refined, "%s/nsrefined-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
211 
212       PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined), file_path_refined, FILE_MODE_WRITE, &viewer_refined));
213       PetscCall(VecView(Q_refined_loc, viewer_refined));
214       PetscCall(DMRestoreLocalVector(user->dm_viz, &Q_refined_loc));
215       PetscCall(DMRestoreGlobalVector(user->dm_viz, &Q_refined));
216       PetscCall(PetscViewerDestroy(&viewer_refined));
217     }
218     PetscCall(DMRestoreLocalVector(user->dm, &Q_loc));
219   }
220 
221   // Save data in a binary file for continuation of simulations
222   if (user->app_ctx->add_stepnum2bin) {
223     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", user->app_ctx->output_dir, step_no));
224   } else {
225     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", user->app_ctx->output_dir));
226   }
227   PetscCall(PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, &viewer));
228 
229   time /= user->units->second;  // Dimensionalize time back
230   PetscCall(HoneeWriteBinaryVec(viewer, Q, time, step_no));
231   PetscCall(PetscViewerDestroy(&viewer));
232   PetscFunctionReturn(PETSC_SUCCESS);
233 }
234 
235 // CSV Monitor
236 PetscErrorCode TSMonitor_WallForce(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
237   User              user = ctx;
238   Vec               G_loc;
239   PetscInt          num_wall = user->app_ctx->wall_forces.num_wall, dim = 3;
240   const PetscInt   *walls  = user->app_ctx->wall_forces.walls;
241   PetscViewer       viewer = user->app_ctx->wall_forces.viewer;
242   PetscViewerFormat format = user->app_ctx->wall_forces.viewer_format;
243   PetscScalar      *reaction_force;
244   PetscBool         is_ascii;
245 
246   PetscFunctionBeginUser;
247   if (!viewer) PetscFunctionReturn(PETSC_SUCCESS);
248   PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
249   PetscCall(PetscCalloc1(num_wall * dim, &reaction_force));
250   PetscCall(Surface_Forces_NS(user->dm, G_loc, num_wall, walls, reaction_force));
251   PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
252 
253   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &is_ascii));
254 
255   if (is_ascii) {
256     if (format == PETSC_VIEWER_ASCII_CSV && !user->app_ctx->wall_forces.header_written) {
257       PetscCall(PetscViewerASCIIPrintf(viewer, "Step,Time,Wall,ForceX,ForceY,ForceZ\n"));
258       user->app_ctx->wall_forces.header_written = PETSC_TRUE;
259     }
260     for (PetscInt w = 0; w < num_wall; w++) {
261       PetscInt wall = walls[w];
262       if (format == PETSC_VIEWER_ASCII_CSV) {
263         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT ",%g,%" PetscInt_FMT ",%g,%g,%g\n", step_no, time, wall,
264                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
265 
266       } else {
267         PetscCall(PetscViewerASCIIPrintf(viewer, "Wall %" PetscInt_FMT " Forces: Force_x = %12g, Force_y = %12g, Force_z = %12g\n", wall,
268                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
269       }
270     }
271   }
272   PetscCall(PetscFree(reaction_force));
273   PetscFunctionReturn(PETSC_SUCCESS);
274 }
275 
276 // User provided TS Monitor
277 PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
278   User user = ctx;
279 
280   PetscFunctionBeginUser;
281   // Print every 'checkpoint_interval' steps
282   if (user->app_ctx->checkpoint_interval <= 0 || step_no % user->app_ctx->checkpoint_interval != 0 ||
283       (user->app_ctx->cont_steps == step_no && step_no != 0)) {
284     PetscFunctionReturn(PETSC_SUCCESS);
285   }
286 
287   PetscCall(WriteOutput(user, Q, step_no, time));
288   PetscFunctionReturn(PETSC_SUCCESS);
289 }
290 
291 // TS: Create, setup, and solve
292 PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, ProblemData problem, Vec *Q, PetscScalar *f_time, TS *ts) {
293   MPI_Comm    comm = user->comm;
294   TSAdapt     adapt;
295   PetscScalar final_time;
296 
297   PetscFunctionBeginUser;
298   PetscCall(TSCreate(comm, ts));
299   PetscCall(TSSetDM(*ts, dm));
300   PetscCall(TSSetApplicationContext(*ts, user));
301   if (phys->implicit) {
302     PetscCall(TSSetType(*ts, TSBDF));
303     if (user->op_ifunction) {
304       PetscCall(TSSetIFunction(*ts, NULL, IFunction_NS, &user));
305     } else {  // Implicit integrators can fall back to using an RHSFunction
306       PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
307     }
308     if (user->mat_ijacobian) {
309       PetscCall(DMTSSetIJacobian(dm, FormIJacobian_NS, &user));
310     }
311   } else {
312     PetscCheck(user->op_rhs_ctx, comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction");
313     PetscCall(TSSetType(*ts, TSRK));
314     PetscCall(TSRKSetType(*ts, TSRK5F));
315     PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
316   }
317   PetscCall(TSSetMaxTime(*ts, 500. * user->units->second));
318   PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER));
319   if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(TSSetErrorIfStepFails(*ts, PETSC_FALSE));
320   PetscCall(TSSetTimeStep(*ts, 1.e-2 * user->units->second));
321   PetscCall(TSGetAdapt(*ts, &adapt));
322   PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * user->units->second, 1.e2 * user->units->second));
323   PetscCall(TSSetFromOptions(*ts));
324   if (user->mat_ijacobian) {
325     if (app_ctx->amat_type && !strcmp(app_ctx->amat_type, MATSHELL)) {
326       SNES snes;
327       KSP  ksp;
328       Mat  Pmat, Amat;
329 
330       PetscCall(TSGetSNES(*ts, &snes));
331       PetscCall(SNESGetKSP(snes, &ksp));
332       PetscCall(CreateSolveOperatorsFromMatCeed(ksp, user->mat_ijacobian, PETSC_FALSE, &Amat, &Pmat));
333       PetscCall(TSSetIJacobian(*ts, user->mat_ijacobian, Pmat, NULL, NULL));
334       PetscCall(MatDestroy(&Amat));
335       PetscCall(MatDestroy(&Pmat));
336     }
337   }
338   user->time_bc_set = -1.0;   // require all BCs be updated
339   if (app_ctx->cont_steps) {  // continue from previous timestep data
340     PetscInt    count;
341     PetscViewer viewer;
342 
343     if (app_ctx->cont_time <= 0) {  // Legacy files did not include step number and time
344       PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_time_file, FILE_MODE_READ, &viewer));
345       PetscCall(PetscViewerBinaryRead(viewer, &app_ctx->cont_time, 1, &count, PETSC_REAL));
346       PetscCall(PetscViewerDestroy(&viewer));
347       PetscCheck(app_ctx->cont_steps != -1, comm, PETSC_ERR_ARG_INCOMP,
348                  "-continue step number not specified, but checkpoint file does not contain a step number (likely written by older code version)");
349     }
350     PetscCall(TSSetTime(*ts, app_ctx->cont_time * user->units->second));
351     PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps));
352   }
353   if (app_ctx->test_type == TESTTYPE_NONE) {
354     PetscCall(TSMonitorSet(*ts, TSMonitor_NS, user, NULL));
355   }
356   if (app_ctx->wall_forces.viewer) {
357     PetscCall(TSMonitorSet(*ts, TSMonitor_WallForce, user, NULL));
358   }
359   if (app_ctx->turb_spanstats_enable) {
360     PetscCall(TSMonitorSet(*ts, TSMonitor_TurbulenceStatistics, user, NULL));
361     CeedScalar previous_time = app_ctx->cont_time * user->units->second;
362     PetscCallCeed(user->ceed,
363                   CeedOperatorSetContextDouble(user->spanstats.op_stats_collect_ctx->op, user->spanstats.previous_time_label, &previous_time));
364   }
365   if (app_ctx->diff_filter_monitor) PetscCall(TSMonitorSet(*ts, TSMonitor_DifferentialFilter, user, NULL));
366 
367   if (app_ctx->sgs_train_enable) {
368     PetscCall(TSMonitorSet(*ts, TSMonitor_SGS_DD_Training, user, NULL));
369     PetscCall(TSSetPostStep(*ts, TSPostStep_SGS_DD_Training));
370   }
371 
372   if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(PrintRunInfo(user, user->phys, problem, *ts));
373   // Solve
374   PetscReal start_time;
375   PetscInt  start_step;
376   PetscCall(TSGetTime(*ts, &start_time));
377   PetscCall(TSGetStepNumber(*ts, &start_step));
378 
379   PetscCall(PetscLogDefaultBegin());  // So we can use PetscLogStageGetPerfInfo without -log_view
380   PetscPreLoadBegin(PETSC_FALSE, "Fluids Solve");
381   PetscCall(TSSetTime(*ts, start_time));
382   PetscCall(TSSetStepNumber(*ts, start_step));
383   if (PetscPreLoadingOn) {
384     // LCOV_EXCL_START
385     SNES      snes;
386     Vec       Q_preload;
387     PetscReal rtol;
388     PetscCall(VecDuplicate(*Q, &Q_preload));
389     PetscCall(VecCopy(*Q, Q_preload));
390     PetscCall(TSGetSNES(*ts, &snes));
391     PetscCall(SNESGetTolerances(snes, NULL, &rtol, NULL, NULL, NULL));
392     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, .99, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
393     PetscCall(TSSetSolution(*ts, Q_preload));
394     PetscCall(TSStep(*ts));
395     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, rtol, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
396     PetscCall(VecDestroy(&Q_preload));
397     // LCOV_EXCL_STOP
398   } else {
399     PetscCall(PetscBarrier((PetscObject)*ts));
400     PetscCall(TSSolve(*ts, *Q));
401   }
402   PetscPreLoadEnd();
403 
404   PetscCall(TSGetSolveTime(*ts, &final_time));
405   *f_time = final_time;
406 
407   if (app_ctx->test_type == TESTTYPE_NONE) {
408     PetscInt step_no;
409     PetscCall(TSGetStepNumber(*ts, &step_no));
410     if (user->app_ctx->checkpoint_interval > 0 || user->app_ctx->checkpoint_interval == -1) {
411       PetscCall(WriteOutput(user, *Q, step_no, final_time));
412     }
413 
414     PetscLogStage      stage_id;
415     PetscEventPerfInfo stage_perf;
416 
417     PetscCall(PetscLogStageGetId("Fluids Solve", &stage_id));
418     PetscCall(PetscLogStageGetPerfInfo(stage_id, &stage_perf));
419     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time taken for solution (sec): %g\n", stage_perf.time));
420   }
421   PetscFunctionReturn(PETSC_SUCCESS);
422 }
423