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