xref: /libCEED/examples/fluids/src/setupts.c (revision 04292f7dad43a72a08e1eaa490df719405885952)
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;
39     Vec      Ones_loc;
40     MPI_Comm comm = PetscObjectComm((PetscObject)dm);
41 
42     PetscCall(DMCreateLocalVector(dm, &Ones_loc));
43     PetscCall(VecSet(Ones_loc, 1));
44     PetscCall(OperatorApplyContextCreate(dm, dm, ceed, op_mass, NULL, NULL, Ones_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(&Ones_loc));
59   }
60 
61   // Cleanup
62   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
63   PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
64   PetscFunctionReturn(PETSC_SUCCESS);
65 }
66 
67 // Insert Boundary values if it's a new time
68 PetscErrorCode UpdateBoundaryValues(User user, Vec Q_loc, PetscReal t) {
69   PetscFunctionBeginUser;
70   if (user->time_bc_set != t) {
71     PetscCall(DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL));
72     user->time_bc_set = t;
73   }
74   PetscFunctionReturn(PETSC_SUCCESS);
75 }
76 
77 // RHS (Explicit time-stepper) function setup
78 //   This is the RHS of the ODE, given as u_t = G(t,u)
79 //   This function takes in a state vector Q and writes into G
80 PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) {
81   User        user = *(User *)user_data;
82   Ceed        ceed = user->ceed;
83   PetscScalar dt;
84   Vec         Q_loc = user->Q_loc;
85 
86   PetscFunctionBeginUser;
87   // Update time dependent data
88   PetscCall(UpdateBoundaryValues(user, Q_loc, t));
89   if (user->phys->solution_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_rhs_ctx->op, user->phys->solution_time_label, &t));
90   PetscCall(TSGetTimeStep(ts, &dt));
91   if (user->phys->timestep_size_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_rhs_ctx->op, user->phys->timestep_size_label, &dt));
92 
93   PetscCall(ApplyCeedOperatorGlobalToGlobal(Q, G, user->op_rhs_ctx));
94 
95   // Inverse of the lumped mass matrix
96   PetscCall(KSPSolve(user->mass_ksp, G, G));
97   PetscFunctionReturn(PETSC_SUCCESS);
98 }
99 
100 // Surface forces function setup
101 static PetscErrorCode Surface_Forces_NS(DM dm, Vec G_loc, PetscInt num_walls, const PetscInt walls[], PetscScalar *reaction_force) {
102   DMLabel            face_label;
103   const PetscScalar *g;
104   PetscInt           dof, dim = 3;
105   MPI_Comm           comm;
106   PetscSection       s;
107 
108   PetscFunctionBeginUser;
109   PetscCall(PetscArrayzero(reaction_force, num_walls * dim));
110   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
111   PetscCall(DMGetLabel(dm, "Face Sets", &face_label));
112   PetscCall(VecGetArrayRead(G_loc, &g));
113   for (PetscInt w = 0; w < num_walls; w++) {
114     const PetscInt wall = walls[w];
115     IS             wall_is;
116     PetscCall(DMGetLocalSection(dm, &s));
117     PetscCall(DMLabelGetStratumIS(face_label, wall, &wall_is));
118     if (wall_is) {  // There exist such points on this process
119       PetscInt        num_points;
120       PetscInt        num_comp = 0;
121       const PetscInt *points;
122       PetscCall(PetscSectionGetFieldComponents(s, 0, &num_comp));
123       PetscCall(ISGetSize(wall_is, &num_points));
124       PetscCall(ISGetIndices(wall_is, &points));
125       for (PetscInt i = 0; i < num_points; i++) {
126         const PetscInt           p = points[i];
127         const StateConservative *r;
128         PetscCall(DMPlexPointLocalRead(dm, p, g, &r));
129         PetscCall(PetscSectionGetDof(s, p, &dof));
130         for (PetscInt node = 0; node < dof / num_comp; node++) {
131           for (PetscInt j = 0; j < 3; j++) {
132             reaction_force[w * dim + j] -= r[node].momentum[j];
133           }
134         }
135       }
136       PetscCall(ISRestoreIndices(wall_is, &points));
137     }
138     PetscCall(ISDestroy(&wall_is));
139   }
140   PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, reaction_force, dim * num_walls, MPIU_SCALAR, MPI_SUM, comm));
141   //  Restore Vectors
142   PetscCall(VecRestoreArrayRead(G_loc, &g));
143   PetscFunctionReturn(PETSC_SUCCESS);
144 }
145 
146 // Implicit time-stepper function setup
147 PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) {
148   User         user = *(User *)user_data;
149   Ceed         ceed = user->ceed;
150   PetscScalar  dt;
151   Vec          Q_loc = user->Q_loc, Q_dot_loc = user->Q_dot_loc, G_loc;
152   PetscMemType q_mem_type, q_dot_mem_type, g_mem_type;
153 
154   PetscFunctionBeginUser;
155   // Get local vectors
156   PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
157 
158   // Update time dependent data
159   PetscCall(UpdateBoundaryValues(user, Q_loc, t));
160   if (user->phys->solution_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_ifunction, user->phys->solution_time_label, &t));
161   PetscCall(TSGetTimeStep(ts, &dt));
162   if (user->phys->timestep_size_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_ifunction, user->phys->timestep_size_label, &dt));
163 
164   // Global-to-local
165   PetscCall(DMGlobalToLocalBegin(user->dm, Q, INSERT_VALUES, Q_loc));
166   PetscCall(DMGlobalToLocalBegin(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
167   PetscCall(DMGlobalToLocalEnd(user->dm, Q, INSERT_VALUES, Q_loc));
168   PetscCall(DMGlobalToLocalEnd(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
169 
170   // Place PETSc vectors in CEED vectors
171   PetscCall(VecReadP2C(Q_loc, &q_mem_type, user->q_ceed));
172   PetscCall(VecReadP2C(Q_dot_loc, &q_dot_mem_type, user->q_dot_ceed));
173   PetscCall(VecP2C(G_loc, &g_mem_type, user->g_ceed));
174 
175   // Apply CEED operator
176   PetscCall(PetscLogEventBegin(FLUIDS_CeedOperatorApply, Q, G, 0, 0));
177   PetscCall(PetscLogGpuTimeBegin());
178   PetscCallCeed(user->ceed, CeedOperatorApply(user->op_ifunction, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE));
179   PetscCall(PetscLogGpuTimeEnd());
180   PetscCall(PetscLogEventEnd(FLUIDS_CeedOperatorApply, Q, G, 0, 0));
181 
182   // Restore vectors
183   PetscCall(VecReadC2P(user->q_ceed, q_mem_type, Q_loc));
184   PetscCall(VecReadC2P(user->q_dot_ceed, q_dot_mem_type, Q_dot_loc));
185   PetscCall(VecC2P(user->g_ceed, g_mem_type, G_loc));
186 
187   if (user->app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) {
188     PetscCall(SgsDDModelApplyIFunction(user, Q_loc, G_loc));
189   }
190 
191   // Local-to-Global
192   PetscCall(VecZeroEntries(G));
193   PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G));
194 
195   // Restore vectors
196   PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
197   PetscFunctionReturn(PETSC_SUCCESS);
198 }
199 
200 static PetscErrorCode FormPreallocation(User user, PetscBool pbdiagonal, Mat J, CeedVector *coo_values) {
201   PetscCount ncoo;
202   PetscInt  *rows_petsc, *cols_petsc;
203   CeedInt   *rows_ceed, *cols_ceed;
204 
205   PetscFunctionBeginUser;
206   if (pbdiagonal) {
207     PetscCallCeed(user->ceed, CeedOperatorLinearAssemblePointBlockDiagonalSymbolic(user->op_ijacobian, &ncoo, &rows_ceed, &cols_ceed));
208   } else {
209     PetscCallCeed(user->ceed, CeedOperatorLinearAssembleSymbolic(user->op_ijacobian, &ncoo, &rows_ceed, &cols_ceed));
210   }
211   PetscCall(IntArrayC2P(ncoo, &rows_ceed, &rows_petsc));
212   PetscCall(IntArrayC2P(ncoo, &cols_ceed, &cols_petsc));
213   PetscCall(MatSetPreallocationCOOLocal(J, ncoo, rows_petsc, cols_petsc));
214   free(rows_petsc);
215   free(cols_petsc);
216   PetscCallCeed(user->ceed, CeedVectorCreate(user->ceed, ncoo, coo_values));
217   PetscFunctionReturn(PETSC_SUCCESS);
218 }
219 
220 static PetscErrorCode FormSetValues(User user, PetscBool pbdiagonal, Mat J, CeedVector coo_values) {
221   CeedMemType        mem_type = CEED_MEM_HOST;
222   const PetscScalar *values;
223   MatType            mat_type;
224 
225   PetscFunctionBeginUser;
226   PetscCall(MatGetType(J, &mat_type));
227   if (strstr(mat_type, "kokkos") || strstr(mat_type, "cusparse")) mem_type = CEED_MEM_DEVICE;
228   if (pbdiagonal) {
229     PetscCall(PetscLogEventBegin(FLUIDS_CeedOperatorAssemblePointBlockDiagonal, J, 0, 0, 0));
230     PetscCall(PetscLogGpuTimeBegin());
231     PetscCallCeed(user->ceed, CeedOperatorLinearAssemblePointBlockDiagonal(user->op_ijacobian, coo_values, CEED_REQUEST_IMMEDIATE));
232     PetscCall(PetscLogGpuTimeEnd());
233     PetscCall(PetscLogEventEnd(FLUIDS_CeedOperatorAssemblePointBlockDiagonal, J, 0, 0, 0));
234   } else {
235     PetscCall(PetscLogEventBegin(FLUIDS_CeedOperatorAssemble, J, 0, 0, 0));
236     PetscCall(PetscLogGpuTimeBegin());
237     PetscCallCeed(user->ceed, CeedOperatorLinearAssemble(user->op_ijacobian, coo_values));
238     PetscCall(PetscLogGpuTimeEnd());
239     PetscCall(PetscLogEventEnd(FLUIDS_CeedOperatorAssemble, J, 0, 0, 0));
240   }
241   PetscCallCeed(user->ceed, CeedVectorGetArrayRead(coo_values, mem_type, &values));
242   PetscCall(MatSetValuesCOO(J, values, INSERT_VALUES));
243   PetscCallCeed(user->ceed, CeedVectorRestoreArrayRead(coo_values, &values));
244   PetscFunctionReturn(PETSC_SUCCESS);
245 }
246 
247 PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, PetscReal shift, Mat J, Mat J_pre, void *user_data) {
248   User      user = *(User *)user_data;
249   Ceed      ceed = user->ceed;
250   PetscBool J_is_shell, J_is_mffd, J_pre_is_shell;
251 
252   PetscFunctionBeginUser;
253   if (user->phys->ijacobian_time_shift_label)
254     PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_ijacobian, user->phys->ijacobian_time_shift_label, &shift));
255   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATMFFD, &J_is_mffd));
256   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATSHELL, &J_is_shell));
257   PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATSHELL, &J_pre_is_shell));
258   if (!user->matrices_set_up) {
259     if (J_is_shell) {
260       OperatorApplyContext op_ijacobian_ctx;
261       OperatorApplyContextCreate(user->dm, user->dm, user->ceed, user->op_ijacobian, user->q_ceed, user->g_ceed, user->Q_dot_loc, NULL,
262                                  &op_ijacobian_ctx);
263       PetscCall(MatShellSetContext(J, op_ijacobian_ctx));
264       PetscCall(MatShellSetContextDestroy(J, (PetscErrorCode(*)(void *))OperatorApplyContextDestroy));
265       PetscCall(MatShellSetOperation(J, MATOP_MULT, (void (*)(void))MatMult_Ceed));
266       PetscCall(MatShellSetOperation(J, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiag_Ceed));
267       PetscCall(MatSetUp(J));
268     }
269     if (!J_pre_is_shell) {
270       PetscCall(FormPreallocation(user, user->app_ctx->pmat_pbdiagonal, J_pre, &user->coo_values_pmat));
271     }
272     if (J != J_pre && !J_is_shell && !J_is_mffd) {
273       PetscCall(FormPreallocation(user, PETSC_FALSE, J, &user->coo_values_amat));
274     }
275     user->matrices_set_up = true;
276   }
277   if (!J_pre_is_shell) {
278     PetscCall(FormSetValues(user, user->app_ctx->pmat_pbdiagonal, J_pre, user->coo_values_pmat));
279   }
280   if (user->coo_values_amat) {
281     PetscCall(FormSetValues(user, PETSC_FALSE, J, user->coo_values_amat));
282   } else if (J_is_mffd) {
283     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
284     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
285   }
286   PetscFunctionReturn(PETSC_SUCCESS);
287 }
288 
289 PetscErrorCode WriteOutput(User user, Vec Q, PetscInt step_no, PetscScalar time) {
290   Vec         Q_loc;
291   char        file_path[PETSC_MAX_PATH_LEN];
292   PetscViewer viewer;
293 
294   PetscFunctionBeginUser;
295   if (user->app_ctx->checkpoint_vtk) {
296     // Set up output
297     PetscCall(DMGetLocalVector(user->dm, &Q_loc));
298     PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec"));
299     PetscCall(VecZeroEntries(Q_loc));
300     PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
301 
302     // Output
303     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
304 
305     PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path, FILE_MODE_WRITE, &viewer));
306     PetscCall(VecView(Q_loc, viewer));
307     PetscCall(PetscViewerDestroy(&viewer));
308     if (user->dm_viz) {
309       Vec         Q_refined, Q_refined_loc;
310       char        file_path_refined[PETSC_MAX_PATH_LEN];
311       PetscViewer viewer_refined;
312 
313       PetscCall(DMGetGlobalVector(user->dm_viz, &Q_refined));
314       PetscCall(DMGetLocalVector(user->dm_viz, &Q_refined_loc));
315       PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined"));
316 
317       PetscCall(MatInterpolate(user->interp_viz, Q, Q_refined));
318       PetscCall(VecZeroEntries(Q_refined_loc));
319       PetscCall(DMGlobalToLocal(user->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc));
320 
321       PetscCall(
322           PetscSNPrintf(file_path_refined, sizeof file_path_refined, "%s/nsrefined-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
323 
324       PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined), file_path_refined, FILE_MODE_WRITE, &viewer_refined));
325       PetscCall(VecView(Q_refined_loc, viewer_refined));
326       PetscCall(DMRestoreLocalVector(user->dm_viz, &Q_refined_loc));
327       PetscCall(DMRestoreGlobalVector(user->dm_viz, &Q_refined));
328       PetscCall(PetscViewerDestroy(&viewer_refined));
329     }
330     PetscCall(DMRestoreLocalVector(user->dm, &Q_loc));
331   }
332 
333   // Save data in a binary file for continuation of simulations
334   if (user->app_ctx->add_stepnum2bin) {
335     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", user->app_ctx->output_dir, step_no));
336   } else {
337     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", user->app_ctx->output_dir));
338   }
339   PetscCall(PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, &viewer));
340 
341   PetscInt32 token = PetscDefined(USE_64BIT_INDICES) ? FLUIDS_FILE_TOKEN_64 : FLUIDS_FILE_TOKEN_32;
342   PetscCall(PetscViewerBinaryWrite(viewer, &token, 1, PETSC_INT32));
343   PetscCall(PetscViewerBinaryWrite(viewer, &step_no, 1, PETSC_INT));
344   time /= user->units->second;  // Dimensionalize time back
345   PetscCall(PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL));
346   PetscCall(VecView(Q, viewer));
347   PetscCall(PetscViewerDestroy(&viewer));
348   PetscFunctionReturn(PETSC_SUCCESS);
349 }
350 
351 // CSV Monitor
352 PetscErrorCode TSMonitor_WallForce(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
353   User              user = ctx;
354   Vec               G_loc;
355   PetscInt          num_wall = user->app_ctx->wall_forces.num_wall, dim = 3;
356   const PetscInt   *walls  = user->app_ctx->wall_forces.walls;
357   PetscViewer       viewer = user->app_ctx->wall_forces.viewer;
358   PetscViewerFormat format = user->app_ctx->wall_forces.viewer_format;
359   PetscScalar      *reaction_force;
360   PetscBool         iascii;
361 
362   PetscFunctionBeginUser;
363   if (!viewer) PetscFunctionReturn(PETSC_SUCCESS);
364   PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
365   PetscCall(PetscMalloc1(num_wall * dim, &reaction_force));
366   PetscCall(Surface_Forces_NS(user->dm, G_loc, num_wall, walls, reaction_force));
367   PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
368 
369   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
370 
371   if (iascii) {
372     if (format == PETSC_VIEWER_ASCII_CSV && !user->app_ctx->wall_forces.header_written) {
373       PetscCall(PetscViewerASCIIPrintf(viewer, "Step,Time,Wall,ForceX,ForceY,ForceZ\n"));
374       user->app_ctx->wall_forces.header_written = PETSC_TRUE;
375     }
376     for (PetscInt w = 0; w < num_wall; w++) {
377       PetscInt wall = walls[w];
378       if (format == PETSC_VIEWER_ASCII_CSV) {
379         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT ",%g,%" PetscInt_FMT ",%g,%g,%g\n", step_no, time, wall,
380                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
381 
382       } else {
383         PetscCall(PetscViewerASCIIPrintf(viewer, "Wall %" PetscInt_FMT " Forces: Force_x = %12g, Force_y = %12g, Force_z = %12g\n", wall,
384                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
385       }
386     }
387   }
388   PetscCall(PetscFree(reaction_force));
389   PetscFunctionReturn(PETSC_SUCCESS);
390 }
391 
392 // User provided TS Monitor
393 PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
394   User user = ctx;
395 
396   PetscFunctionBeginUser;
397   // Print every 'checkpoint_interval' steps
398   if (user->app_ctx->checkpoint_interval <= 0 || step_no % user->app_ctx->checkpoint_interval != 0 ||
399       (user->app_ctx->cont_steps == step_no && step_no != 0)) {
400     PetscFunctionReturn(PETSC_SUCCESS);
401   }
402 
403   PetscCall(WriteOutput(user, Q, step_no, time));
404   PetscFunctionReturn(PETSC_SUCCESS);
405 }
406 
407 // TS: Create, setup, and solve
408 PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, Vec *Q, PetscScalar *f_time, TS *ts) {
409   MPI_Comm    comm = user->comm;
410   TSAdapt     adapt;
411   PetscScalar final_time;
412 
413   PetscFunctionBeginUser;
414   PetscCall(TSCreate(comm, ts));
415   PetscCall(TSSetDM(*ts, dm));
416   PetscCall(TSSetApplicationContext(*ts, user));
417   if (phys->implicit) {
418     PetscCall(TSSetType(*ts, TSBDF));
419     if (user->op_ifunction) {
420       PetscCall(TSSetIFunction(*ts, NULL, IFunction_NS, &user));
421     } else {  // Implicit integrators can fall back to using an RHSFunction
422       PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
423     }
424     if (user->op_ijacobian) {
425       PetscCall(DMTSSetIJacobian(dm, FormIJacobian_NS, &user));
426       if (app_ctx->amat_type) {
427         Mat Pmat, Amat;
428         PetscCall(DMCreateMatrix(dm, &Pmat));
429         PetscCall(DMSetMatType(dm, app_ctx->amat_type));
430         PetscCall(DMCreateMatrix(dm, &Amat));
431         PetscCall(TSSetIJacobian(*ts, Amat, Pmat, NULL, NULL));
432         PetscCall(MatDestroy(&Amat));
433         PetscCall(MatDestroy(&Pmat));
434       }
435     }
436   } else {
437     PetscCheck(user->op_rhs_ctx, comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction");
438     PetscCall(TSSetType(*ts, TSRK));
439     PetscCall(TSRKSetType(*ts, TSRK5F));
440     PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
441   }
442   PetscCall(TSSetMaxTime(*ts, 500. * user->units->second));
443   PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER));
444   if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(TSSetErrorIfStepFails(*ts, PETSC_FALSE));
445   PetscCall(TSSetTimeStep(*ts, 1.e-2 * user->units->second));
446   PetscCall(TSGetAdapt(*ts, &adapt));
447   PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * user->units->second, 1.e2 * user->units->second));
448   PetscCall(TSSetFromOptions(*ts));
449   user->time_bc_set = -1.0;   // require all BCs be updated
450   if (app_ctx->cont_steps) {  // continue from previous timestep data
451     PetscInt    count;
452     PetscViewer viewer;
453 
454     if (app_ctx->cont_time <= 0) {  // Legacy files did not include step number and time
455       PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_time_file, FILE_MODE_READ, &viewer));
456       PetscCall(PetscViewerBinaryRead(viewer, &app_ctx->cont_time, 1, &count, PETSC_REAL));
457       PetscCall(PetscViewerDestroy(&viewer));
458       PetscCheck(app_ctx->cont_steps != -1, comm, PETSC_ERR_ARG_INCOMP,
459                  "-continue step number not specified, but checkpoint file does not contain a step number (likely written by older code version)");
460     }
461     PetscCall(TSSetTime(*ts, app_ctx->cont_time * user->units->second));
462     PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps));
463   }
464   if (app_ctx->test_type == TESTTYPE_NONE) {
465     PetscCall(TSMonitorSet(*ts, TSMonitor_NS, user, NULL));
466   }
467   if (app_ctx->wall_forces.viewer) {
468     PetscCall(TSMonitorSet(*ts, TSMonitor_WallForce, user, NULL));
469   }
470   if (app_ctx->turb_spanstats_enable) {
471     PetscCall(TSMonitorSet(*ts, TSMonitor_TurbulenceStatistics, user, NULL));
472     CeedScalar previous_time = app_ctx->cont_time * user->units->second;
473     PetscCallCeed(user->ceed,
474                   CeedOperatorSetContextDouble(user->spanstats.op_stats_collect_ctx->op, user->spanstats.previous_time_label, &previous_time));
475   }
476   if (app_ctx->diff_filter_monitor) PetscCall(TSMonitorSet(*ts, TSMonitor_DifferentialFilter, user, NULL));
477 
478   if (app_ctx->sgs_train_enable) {
479     PetscCall(TSMonitorSet(*ts, TSMonitor_SGS_DD_Training, user, NULL));
480     PetscCall(TSSetPostStep(*ts, TSPostStep_SGS_DD_Training));
481   }
482   // Solve
483   PetscReal start_time;
484   PetscInt  start_step;
485   PetscCall(TSGetTime(*ts, &start_time));
486   PetscCall(TSGetStepNumber(*ts, &start_step));
487 
488   PetscCall(PetscLogDefaultBegin());  // So we can use PetscLogStageGetPerfInfo without -log_view
489   PetscPreLoadBegin(PETSC_FALSE, "Fluids Solve");
490   PetscCall(TSSetTime(*ts, start_time));
491   PetscCall(TSSetStepNumber(*ts, start_step));
492   if (PetscPreLoadingOn) {
493     // LCOV_EXCL_START
494     SNES      snes;
495     Vec       Q_preload;
496     PetscReal rtol;
497     PetscCall(VecDuplicate(*Q, &Q_preload));
498     PetscCall(VecCopy(*Q, Q_preload));
499     PetscCall(TSGetSNES(*ts, &snes));
500     PetscCall(SNESGetTolerances(snes, NULL, &rtol, NULL, NULL, NULL));
501     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, .99, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
502     PetscCall(TSSetSolution(*ts, Q_preload));
503     PetscCall(TSStep(*ts));
504     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, rtol, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
505     PetscCall(VecDestroy(&Q_preload));
506     // LCOV_EXCL_STOP
507   } else {
508     PetscCall(PetscBarrier((PetscObject)*ts));
509     PetscCall(TSSolve(*ts, *Q));
510   }
511   PetscPreLoadEnd();
512 
513   PetscCall(TSGetSolveTime(*ts, &final_time));
514   *f_time = final_time;
515 
516   if (app_ctx->test_type == TESTTYPE_NONE) {
517     PetscInt step_no;
518     PetscCall(TSGetStepNumber(*ts, &step_no));
519     if (user->app_ctx->checkpoint_interval > 0 || user->app_ctx->checkpoint_interval == -1) {
520       PetscCall(WriteOutput(user, *Q, step_no, final_time));
521     }
522 
523     PetscLogStage      stage_id;
524     PetscEventPerfInfo stage_perf;
525 
526     PetscCall(PetscLogStageGetId("Fluids Solve", &stage_id));
527     PetscCall(PetscLogStageGetPerfInfo(stage_id, &stage_perf));
528     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time taken for solution (sec): %g\n", stage_perf.time));
529   }
530   PetscFunctionReturn(PETSC_SUCCESS);
531 }
532