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