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