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