xref: /libCEED/examples/fluids/src/setupts.c (revision ef080ff9ce83a1650979d1b767b88f0d6e3ee6ca)
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 "../navierstokes.h"
12 
13 // Compute mass matrix for explicit scheme
14 PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm, CeedData ceed_data, Vec M) {
15   Vec           M_loc;
16   CeedQFunction qf_mass;
17   CeedOperator  op_mass;
18   CeedVector    m_ceed, ones_vec;
19   CeedInt       num_comp_q, q_data_size;
20   PetscFunctionBeginUser;
21 
22   // CEED Restriction
23   CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q);
24   CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &q_data_size);
25   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &m_ceed, NULL);
26   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &ones_vec, NULL);
27   CeedVectorSetValue(ones_vec, 1.0);
28 
29   // CEED QFunction
30   PetscCall(CreateMassQFunction(ceed, num_comp_q, q_data_size, &qf_mass));
31 
32   // CEED Operator
33   CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass);
34   CeedOperatorSetField(op_mass, "u", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
35   CeedOperatorSetField(op_mass, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data);
36   CeedOperatorSetField(op_mass, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
37 
38   // Place PETSc vector in CEED vector
39   CeedScalar  *m;
40   PetscMemType m_mem_type;
41   PetscCall(DMGetLocalVector(dm, &M_loc));
42   PetscCall(VecGetArrayAndMemType(M_loc, (PetscScalar **)&m, &m_mem_type));
43   CeedVectorSetArray(m_ceed, MemTypeP2C(m_mem_type), CEED_USE_POINTER, m);
44 
45   // Apply CEED Operator
46   CeedOperatorApply(op_mass, ones_vec, m_ceed, CEED_REQUEST_IMMEDIATE);
47 
48   // Restore vectors
49   CeedVectorTakeArray(m_ceed, MemTypeP2C(m_mem_type), NULL);
50   PetscCall(VecRestoreArrayReadAndMemType(M_loc, (const PetscScalar **)&m));
51 
52   // Local-to-Global
53   PetscCall(VecZeroEntries(M));
54   PetscCall(DMLocalToGlobal(dm, M_loc, ADD_VALUES, M));
55   PetscCall(DMRestoreLocalVector(dm, &M_loc));
56 
57   // Invert diagonally lumped mass vector for RHS function
58   PetscCall(VecReciprocal(M));
59 
60   // Cleanup
61   CeedVectorDestroy(&ones_vec);
62   CeedVectorDestroy(&m_ceed);
63   CeedQFunctionDestroy(&qf_mass);
64   CeedOperatorDestroy(&op_mass);
65 
66   PetscFunctionReturn(0);
67 }
68 
69 // RHS (Explicit time-stepper) function setup
70 //   This is the RHS of the ODE, given as u_t = G(t,u)
71 //   This function takes in a state vector Q and writes into G
72 PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) {
73   User         user = *(User *)user_data;
74   PetscScalar *q, *g;
75   Vec          Q_loc = user->Q_loc, G_loc;
76   PetscMemType q_mem_type, g_mem_type;
77   PetscFunctionBeginUser;
78 
79   // Get local vector
80   PetscCall(DMGetLocalVector(user->dm, &G_loc));
81 
82   // Update time dependent data
83   if (user->time != t) {
84     PetscCall(DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL));
85     if (user->phys->solution_time_label) {
86       CeedOperatorContextSetDouble(user->op_rhs, user->phys->solution_time_label, &t);
87     }
88     user->time = t;
89   }
90   if (user->phys->timestep_size_label) {
91     PetscScalar dt;
92     PetscCall(TSGetTimeStep(ts, &dt));
93     if (user->dt != dt) {
94       CeedOperatorContextSetDouble(user->op_rhs, user->phys->timestep_size_label, &dt);
95       user->dt = dt;
96     }
97   }
98 
99   // Global-to-local
100   PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
101 
102   // Place PETSc vectors in CEED vectors
103   PetscCall(VecGetArrayReadAndMemType(Q_loc, (const PetscScalar **)&q, &q_mem_type));
104   PetscCall(VecGetArrayAndMemType(G_loc, &g, &g_mem_type));
105   CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER, q);
106   CeedVectorSetArray(user->g_ceed, MemTypeP2C(g_mem_type), CEED_USE_POINTER, g);
107 
108   // Apply CEED operator
109   CeedOperatorApply(user->op_rhs, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE);
110 
111   // Restore vectors
112   CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL);
113   CeedVectorTakeArray(user->g_ceed, MemTypeP2C(g_mem_type), NULL);
114   PetscCall(VecRestoreArrayReadAndMemType(Q_loc, (const PetscScalar **)&q));
115   PetscCall(VecRestoreArrayAndMemType(G_loc, &g));
116 
117   // Local-to-Global
118   PetscCall(VecZeroEntries(G));
119   PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G));
120 
121   // Inverse of the lumped mass matrix (M is Minv)
122   PetscCall(VecPointwiseMult(G, G, user->M));
123 
124   // Restore vectors
125   PetscCall(DMRestoreLocalVector(user->dm, &G_loc));
126 
127   PetscFunctionReturn(0);
128 }
129 
130 // Implicit time-stepper function setup
131 PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) {
132   User               user = *(User *)user_data;
133   const PetscScalar *q, *q_dot;
134   PetscScalar       *g;
135   Vec                Q_loc = user->Q_loc, Q_dot_loc = user->Q_dot_loc, G_loc;
136   PetscMemType       q_mem_type, q_dot_mem_type, g_mem_type;
137   PetscFunctionBeginUser;
138 
139   // Get local vectors
140   PetscCall(DMGetLocalVector(user->dm, &G_loc));
141 
142   // Update time dependent data
143   if (user->time != t) {
144     PetscCall(DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL));
145     if (user->phys->solution_time_label) {
146       CeedOperatorContextSetDouble(user->op_ifunction, user->phys->solution_time_label, &t);
147     }
148     user->time = t;
149   }
150   if (user->phys->timestep_size_label) {
151     PetscScalar dt;
152     PetscCall(TSGetTimeStep(ts, &dt));
153     if (user->dt != dt) {
154       CeedOperatorContextSetDouble(user->op_ifunction, user->phys->timestep_size_label, &dt);
155       user->dt = dt;
156     }
157   }
158 
159   // Global-to-local
160   PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
161   PetscCall(DMGlobalToLocal(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
162 
163   // Place PETSc vectors in CEED vectors
164   PetscCall(VecGetArrayReadAndMemType(Q_loc, &q, &q_mem_type));
165   PetscCall(VecGetArrayReadAndMemType(Q_dot_loc, &q_dot, &q_dot_mem_type));
166   PetscCall(VecGetArrayAndMemType(G_loc, &g, &g_mem_type));
167   CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER, (PetscScalar *)q);
168   CeedVectorSetArray(user->q_dot_ceed, MemTypeP2C(q_dot_mem_type), CEED_USE_POINTER, (PetscScalar *)q_dot);
169   CeedVectorSetArray(user->g_ceed, MemTypeP2C(g_mem_type), CEED_USE_POINTER, g);
170 
171   // Apply CEED operator
172   CeedOperatorApply(user->op_ifunction, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE);
173 
174   // Restore vectors
175   CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL);
176   CeedVectorTakeArray(user->q_dot_ceed, MemTypeP2C(q_dot_mem_type), NULL);
177   CeedVectorTakeArray(user->g_ceed, MemTypeP2C(g_mem_type), NULL);
178   PetscCall(VecRestoreArrayReadAndMemType(Q_loc, &q));
179   PetscCall(VecRestoreArrayReadAndMemType(Q_dot_loc, &q_dot));
180   PetscCall(VecRestoreArrayAndMemType(G_loc, &g));
181 
182   // Local-to-Global
183   PetscCall(VecZeroEntries(G));
184   PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G));
185 
186   // Restore vectors
187   PetscCall(DMRestoreLocalVector(user->dm, &G_loc));
188 
189   PetscFunctionReturn(0);
190 }
191 
192 static PetscErrorCode MatMult_NS_IJacobian(Mat J, Vec Q, Vec G) {
193   User               user;
194   const PetscScalar *q;
195   PetscScalar       *g;
196   PetscMemType       q_mem_type, g_mem_type;
197   PetscFunctionBeginUser;
198   PetscCall(MatShellGetContext(J, &user));
199   Vec Q_loc = user->Q_dot_loc,  // Note - Q_dot_loc has zero BCs
200       G_loc;
201 
202   // Get local vectors
203   PetscCall(DMGetLocalVector(user->dm, &G_loc));
204 
205   // Global-to-local
206   PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
207 
208   // Place PETSc vectors in CEED vectors
209   PetscCall(VecGetArrayReadAndMemType(Q_loc, &q, &q_mem_type));
210   PetscCall(VecGetArrayAndMemType(G_loc, &g, &g_mem_type));
211   CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER, (PetscScalar *)q);
212   CeedVectorSetArray(user->g_ceed, MemTypeP2C(g_mem_type), CEED_USE_POINTER, g);
213 
214   // Apply CEED operator
215   CeedOperatorApply(user->op_ijacobian, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE);
216 
217   // Restore vectors
218   CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL);
219   CeedVectorTakeArray(user->g_ceed, MemTypeP2C(g_mem_type), NULL);
220   PetscCall(VecRestoreArrayReadAndMemType(Q_loc, &q));
221   PetscCall(VecRestoreArrayAndMemType(G_loc, &g));
222 
223   // Local-to-Global
224   PetscCall(VecZeroEntries(G));
225   PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G));
226 
227   // Restore vectors
228   PetscCall(DMRestoreLocalVector(user->dm, &G_loc));
229   PetscFunctionReturn(0);
230 }
231 
232 PetscErrorCode MatGetDiagonal_NS_IJacobian(Mat A, Vec D) {
233   User         user;
234   Vec          D_loc;
235   PetscScalar *d;
236   PetscMemType mem_type;
237 
238   PetscFunctionBeginUser;
239   MatShellGetContext(A, &user);
240   PetscCall(DMGetLocalVector(user->dm, &D_loc));
241   PetscCall(VecGetArrayAndMemType(D_loc, &d, &mem_type));
242   CeedVectorSetArray(user->g_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, d);
243   CeedOperatorLinearAssembleDiagonal(user->op_ijacobian, user->g_ceed, CEED_REQUEST_IMMEDIATE);
244   CeedVectorTakeArray(user->g_ceed, MemTypeP2C(mem_type), NULL);
245   PetscCall(VecRestoreArrayAndMemType(D_loc, &d));
246   PetscCall(VecZeroEntries(D));
247   PetscCall(DMLocalToGlobal(user->dm, D_loc, ADD_VALUES, D));
248   PetscCall(DMRestoreLocalVector(user->dm, &D_loc));
249   VecViewFromOptions(D, NULL, "-diag_vec_view");
250   PetscFunctionReturn(0);
251 }
252 
253 static PetscErrorCode FormPreallocation(User user, PetscBool pbdiagonal, Mat J, CeedVector *coo_values) {
254   PetscCount ncoo;
255   PetscInt  *rows, *cols;
256 
257   PetscFunctionBeginUser;
258   if (pbdiagonal) {
259     CeedSize l_size;
260     CeedOperatorGetActiveVectorLengths(user->op_ijacobian, &l_size, NULL);
261     ncoo = l_size * 5;
262     rows = malloc(ncoo * sizeof(rows[0]));
263     cols = malloc(ncoo * sizeof(cols[0]));
264     for (PetscCount n = 0; n < l_size / 5; n++) {
265       for (PetscInt i = 0; i < 5; i++) {
266         for (PetscInt j = 0; j < 5; j++) {
267           rows[(n * 5 + i) * 5 + j] = n * 5 + i;
268           cols[(n * 5 + i) * 5 + j] = n * 5 + j;
269         }
270       }
271     }
272   } else {
273     PetscCall(CeedOperatorLinearAssembleSymbolic(user->op_ijacobian, &ncoo, &rows, &cols));
274   }
275   PetscCall(MatSetPreallocationCOOLocal(J, ncoo, rows, cols));
276   free(rows);
277   free(cols);
278   CeedVectorCreate(user->ceed, ncoo, coo_values);
279   PetscFunctionReturn(0);
280 }
281 
282 static PetscErrorCode FormSetValues(User user, PetscBool pbdiagonal, Mat J, CeedVector coo_values) {
283   CeedMemType        mem_type = CEED_MEM_HOST;
284   const PetscScalar *values;
285   MatType            mat_type;
286 
287   PetscFunctionBeginUser;
288   PetscCall(MatGetType(J, &mat_type));
289   if (strstr(mat_type, "kokkos") || strstr(mat_type, "cusparse")) mem_type = CEED_MEM_DEVICE;
290   if (user->app_ctx->pmat_pbdiagonal) {
291     CeedOperatorLinearAssemblePointBlockDiagonal(user->op_ijacobian, coo_values, CEED_REQUEST_IMMEDIATE);
292   } else {
293     CeedOperatorLinearAssemble(user->op_ijacobian, coo_values);
294   }
295   CeedVectorGetArrayRead(coo_values, mem_type, &values);
296   PetscCall(MatSetValuesCOO(J, values, INSERT_VALUES));
297   CeedVectorRestoreArrayRead(coo_values, &values);
298   PetscFunctionReturn(0);
299 }
300 
301 PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, PetscReal shift, Mat J, Mat J_pre, void *user_data) {
302   User      user = *(User *)user_data;
303   PetscBool J_is_shell, J_is_mffd, J_pre_is_shell;
304   PetscFunctionBeginUser;
305   if (user->phys->ijacobian_time_shift_label) CeedOperatorContextSetDouble(user->op_ijacobian, user->phys->ijacobian_time_shift_label, &shift);
306   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATMFFD, &J_is_mffd));
307   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATSHELL, &J_is_shell));
308   PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATSHELL, &J_pre_is_shell));
309   if (!user->matrices_set_up) {
310     if (J_is_shell) {
311       PetscCall(MatShellSetContext(J, user));
312       PetscCall(MatShellSetOperation(J, MATOP_MULT, (void (*)(void))MatMult_NS_IJacobian));
313       PetscCall(MatShellSetOperation(J, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_NS_IJacobian));
314       PetscCall(MatSetUp(J));
315     }
316     if (!J_pre_is_shell) {
317       PetscCall(FormPreallocation(user, user->app_ctx->pmat_pbdiagonal, J_pre, &user->coo_values_pmat));
318     }
319     if (J != J_pre && !J_is_shell && !J_is_mffd) {
320       PetscCall(FormPreallocation(user, PETSC_FALSE, J, &user->coo_values_amat));
321     }
322     user->matrices_set_up = true;
323   }
324   if (!J_pre_is_shell) {
325     PetscCall(FormSetValues(user, user->app_ctx->pmat_pbdiagonal, J_pre, user->coo_values_pmat));
326   }
327   if (user->coo_values_amat) {
328     PetscCall(FormSetValues(user, PETSC_FALSE, J, user->coo_values_amat));
329   } else if (J_is_mffd) {
330     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
331     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
332   }
333   PetscFunctionReturn(0);
334 }
335 
336 PetscErrorCode WriteOutput(User user, Vec Q, PetscInt step_no, PetscScalar time) {
337   Vec         Q_loc;
338   char        file_path[PETSC_MAX_PATH_LEN];
339   PetscViewer viewer;
340   PetscFunctionBeginUser;
341 
342   if (user->app_ctx->checkpoint_vtk) {
343     // Set up output
344     PetscCall(DMGetLocalVector(user->dm, &Q_loc));
345     PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec"));
346     PetscCall(VecZeroEntries(Q_loc));
347     PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
348 
349     // Output
350     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
351 
352     PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path, FILE_MODE_WRITE, &viewer));
353     PetscCall(VecView(Q_loc, viewer));
354     PetscCall(PetscViewerDestroy(&viewer));
355     if (user->dm_viz) {
356       Vec         Q_refined, Q_refined_loc;
357       char        file_path_refined[PETSC_MAX_PATH_LEN];
358       PetscViewer viewer_refined;
359 
360       PetscCall(DMGetGlobalVector(user->dm_viz, &Q_refined));
361       PetscCall(DMGetLocalVector(user->dm_viz, &Q_refined_loc));
362       PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined"));
363 
364       PetscCall(MatInterpolate(user->interp_viz, Q, Q_refined));
365       PetscCall(VecZeroEntries(Q_refined_loc));
366       PetscCall(DMGlobalToLocal(user->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc));
367 
368       PetscCall(
369           PetscSNPrintf(file_path_refined, sizeof file_path_refined, "%s/nsrefined-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
370 
371       PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined), file_path_refined, FILE_MODE_WRITE, &viewer_refined));
372       PetscCall(VecView(Q_refined_loc, viewer_refined));
373       PetscCall(DMRestoreLocalVector(user->dm_viz, &Q_refined_loc));
374       PetscCall(DMRestoreGlobalVector(user->dm_viz, &Q_refined));
375       PetscCall(PetscViewerDestroy(&viewer_refined));
376     }
377     PetscCall(DMRestoreLocalVector(user->dm, &Q_loc));
378   }
379 
380   // Save data in a binary file for continuation of simulations
381   if (user->app_ctx->add_stepnum2bin) {
382     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", user->app_ctx->output_dir, step_no));
383   } else {
384     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", user->app_ctx->output_dir));
385   }
386   PetscCall(PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, &viewer));
387 
388   PetscInt token = FLUIDS_FILE_TOKEN;
389   PetscCall(PetscViewerBinaryWrite(viewer, &token, 1, PETSC_INT));
390   PetscCall(PetscViewerBinaryWrite(viewer, &step_no, 1, PETSC_INT));
391   time /= user->units->second;  // Dimensionalize time back
392   PetscCall(PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL));
393   PetscCall(VecView(Q, viewer));
394   PetscCall(PetscViewerDestroy(&viewer));
395   PetscFunctionReturn(0);
396 }
397 
398 // User provided TS Monitor
399 PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
400   User user = ctx;
401   PetscFunctionBeginUser;
402 
403   // Print every 'checkpoint_interval' steps
404   if (user->app_ctx->checkpoint_interval <= 0 || step_no % user->app_ctx->checkpoint_interval != 0 ||
405       (user->app_ctx->cont_steps == step_no && step_no != 0))
406     PetscFunctionReturn(0);
407 
408   PetscCall(WriteOutput(user, Q, step_no, time));
409 
410   PetscFunctionReturn(0);
411 }
412 
413 // TS: Create, setup, and solve
414 PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, Vec *Q, PetscScalar *f_time, TS *ts) {
415   MPI_Comm    comm = user->comm;
416   TSAdapt     adapt;
417   PetscScalar final_time;
418   PetscFunctionBeginUser;
419 
420   PetscCall(TSCreate(comm, ts));
421   PetscCall(TSSetDM(*ts, dm));
422   if (phys->implicit) {
423     PetscCall(TSSetType(*ts, TSBDF));
424     if (user->op_ifunction) {
425       PetscCall(TSSetIFunction(*ts, NULL, IFunction_NS, &user));
426     } else {  // Implicit integrators can fall back to using an RHSFunction
427       PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
428     }
429     if (user->op_ijacobian) {
430       PetscCall(DMTSSetIJacobian(dm, FormIJacobian_NS, &user));
431       if (app_ctx->amat_type) {
432         Mat Pmat, Amat;
433         PetscCall(DMCreateMatrix(dm, &Pmat));
434         PetscCall(DMSetMatType(dm, app_ctx->amat_type));
435         PetscCall(DMCreateMatrix(dm, &Amat));
436         PetscCall(TSSetIJacobian(*ts, Amat, Pmat, NULL, NULL));
437         PetscCall(MatDestroy(&Amat));
438         PetscCall(MatDestroy(&Pmat));
439       }
440     }
441   } else {
442     if (!user->op_rhs) SETERRQ(comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction");
443     PetscCall(TSSetType(*ts, TSRK));
444     PetscCall(TSRKSetType(*ts, TSRK5F));
445     PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
446   }
447   PetscCall(TSSetMaxTime(*ts, 500. * user->units->second));
448   PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER));
449   PetscCall(TSSetErrorIfStepFails(*ts, PETSC_FALSE));
450   PetscCall(TSSetTimeStep(*ts, 1.e-2 * user->units->second));
451   if (app_ctx->test_mode) {
452     PetscCall(TSSetMaxSteps(*ts, 10));
453   }
454   PetscCall(TSGetAdapt(*ts, &adapt));
455   PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * user->units->second, 1.e2 * user->units->second));
456   PetscCall(TSSetFromOptions(*ts));
457   user->time = -1.0;  // require all BCs and ctx to be updated
458   user->dt   = -1.0;
459   if (!app_ctx->cont_steps) {  // print initial condition
460     if (!app_ctx->test_mode) {
461       PetscCall(TSMonitor_NS(*ts, 0, 0., *Q, user));
462     }
463   } else {  // continue from time of last output
464     PetscInt    count;
465     PetscViewer viewer;
466 
467     if (app_ctx->cont_time <= 0) {  // Legacy files did not include step number and time
468       PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_time_file, FILE_MODE_READ, &viewer));
469       PetscCall(PetscViewerBinaryRead(viewer, &app_ctx->cont_time, 1, &count, PETSC_REAL));
470       PetscCall(PetscViewerDestroy(&viewer));
471       PetscCheck(app_ctx->cont_steps != -1, comm, PETSC_ERR_ARG_INCOMP,
472                  "-continue step number not specified, but checkpoint file does not contain a step number (likely written by older code version)");
473     }
474     PetscCall(TSSetTime(*ts, app_ctx->cont_time * user->units->second));
475     PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps));
476   }
477   if (!app_ctx->test_mode) {
478     PetscCall(TSMonitorSet(*ts, TSMonitor_NS, user, NULL));
479   }
480 
481   // Solve
482   PetscReal start_time;
483   PetscInt  start_step;
484   PetscCall(TSGetTime(*ts, &start_time));
485   PetscCall(TSGetStepNumber(*ts, &start_step));
486 
487   PetscPreLoadBegin(PETSC_FALSE, "Fluids Solve");
488   PetscCall(TSSetTime(*ts, start_time));
489   PetscCall(TSSetStepNumber(*ts, start_step));
490   if (PetscPreLoadingOn) {
491     // LCOV_EXCL_START
492     SNES      snes;
493     Vec       Q_preload;
494     PetscReal rtol;
495     PetscCall(VecDuplicate(*Q, &Q_preload));
496     PetscCall(VecCopy(*Q, Q_preload));
497     PetscCall(TSGetSNES(*ts, &snes));
498     PetscCall(SNESGetTolerances(snes, NULL, &rtol, NULL, NULL, NULL));
499     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, .99, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
500     PetscCall(TSSetSolution(*ts, *Q));
501     PetscCall(TSStep(*ts));
502     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, rtol, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
503     PetscCall(VecDestroy(&Q_preload));
504     // LCOV_EXCL_STOP
505   } else {
506     PetscCall(PetscBarrier((PetscObject)*ts));
507     PetscCall(TSSolve(*ts, *Q));
508   }
509   PetscPreLoadEnd();
510 
511   PetscCall(TSGetSolveTime(*ts, &final_time));
512   *f_time = final_time;
513 
514   if (!app_ctx->test_mode) {
515     if (user->app_ctx->checkpoint_interval > 0 || user->app_ctx->checkpoint_interval == -1) {
516       PetscInt step_no;
517       PetscCall(TSGetStepNumber(*ts, &step_no));
518       PetscCall(WriteOutput(user, *Q, step_no, final_time));
519     }
520 
521     PetscLogEvent stage_id;
522     PetscStageLog stage_log;
523 
524     PetscCall(PetscLogStageGetId("Fluids Solve", &stage_id));
525     PetscCall(PetscLogGetStageLog(&stage_log));
526     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time taken for solution (sec): %g\n", stage_log->stageInfo[stage_id].perfInfo.time));
527   }
528   PetscFunctionReturn(0);
529 }
530