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