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