xref: /libCEED/examples/fluids/src/setupts.c (revision 01f0e615c90f4b64214c5f157ea3402fb1010229)
13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
377841947SLeila Ghaffari //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
577841947SLeila Ghaffari //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
777841947SLeila Ghaffari 
877841947SLeila Ghaffari /// @file
977841947SLeila Ghaffari /// Time-stepping functions for Navier-Stokes example using PETSc
1077841947SLeila Ghaffari 
1149aac155SJeremy L Thompson #include <ceed.h>
1249aac155SJeremy L Thompson #include <petscdmplex.h>
1349aac155SJeremy L Thompson #include <petscts.h>
1449aac155SJeremy L Thompson 
1577841947SLeila Ghaffari #include "../navierstokes.h"
16ca69d878SAdeleke O. Bankole #include "../qfunctions/newtonian_state.h"
1777841947SLeila Ghaffari 
1877841947SLeila Ghaffari // Compute mass matrix for explicit scheme
192b730f8bSJeremy L Thompson PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm, CeedData ceed_data, Vec M) {
2077841947SLeila Ghaffari   CeedQFunction        qf_mass;
2177841947SLeila Ghaffari   CeedOperator         op_mass;
223221f4d3SJames Wright   OperatorApplyContext op_mass_ctx;
233221f4d3SJames Wright   Vec                  Ones_loc;
2477841947SLeila Ghaffari   CeedInt              num_comp_q, q_data_size;
2577841947SLeila Ghaffari 
26f17d818dSJames Wright   PetscFunctionBeginUser;
27a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q));
28a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &q_data_size));
2977841947SLeila Ghaffari 
30ef080ff9SJames Wright   PetscCall(CreateMassQFunction(ceed, num_comp_q, q_data_size, &qf_mass));
31a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass));
32a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
33356036faSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
34a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
3577841947SLeila Ghaffari 
363221f4d3SJames Wright   PetscCall(OperatorApplyContextCreate(NULL, dm, ceed, op_mass, NULL, NULL, NULL, NULL, &op_mass_ctx));
3777841947SLeila Ghaffari 
383221f4d3SJames Wright   PetscCall(DMGetLocalVector(dm, &Ones_loc));
393221f4d3SJames Wright   PetscCall(VecSet(Ones_loc, 1));
403221f4d3SJames Wright   PetscCall(ApplyCeedOperatorLocalToGlobal(Ones_loc, M, op_mass_ctx));
4177841947SLeila Ghaffari 
4277841947SLeila Ghaffari   // Invert diagonally lumped mass vector for RHS function
432b730f8bSJeremy L Thompson   PetscCall(VecReciprocal(M));
4477841947SLeila Ghaffari 
4577841947SLeila Ghaffari   // Cleanup
463221f4d3SJames Wright   PetscCall(OperatorApplyContextDestroy(op_mass_ctx));
473221f4d3SJames Wright   PetscCall(DMRestoreLocalVector(dm, &Ones_loc));
48a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
49a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
50ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
5177841947SLeila Ghaffari }
5277841947SLeila Ghaffari 
53a0b9a424SJames Wright // Insert Boundary values if it's a new time
54a0b9a424SJames Wright PetscErrorCode UpdateBoundaryValues(User user, Vec Q_loc, PetscReal t) {
55a0b9a424SJames Wright   PetscFunctionBeginUser;
56a0b9a424SJames Wright   if (user->time_bc_set != t) {
57a0b9a424SJames Wright     PetscCall(DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL));
58a0b9a424SJames Wright     user->time_bc_set = t;
59a0b9a424SJames Wright   }
60ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
61a0b9a424SJames Wright }
62a0b9a424SJames Wright 
6315637395SJames Wright // @brief Update the context label value to new value if necessary.
6415637395SJames Wright // @note This only supports labels with scalar label values (ie. not arrays)
65a424bcd0SJames Wright PetscErrorCode UpdateContextLabel(Ceed ceed, MPI_Comm comm, PetscScalar update_value, CeedOperator op, CeedContextFieldLabel label) {
6615637395SJames Wright   PetscScalar label_value;
67e42c1df6SJames Wright 
6815637395SJames Wright   PetscFunctionBeginUser;
6915637395SJames Wright   PetscCheck(label, comm, PETSC_ERR_ARG_BADPTR, "Label should be non-NULL");
7015637395SJames Wright 
7115637395SJames Wright   {
7215637395SJames Wright     size_t             num_elements;
7315637395SJames Wright     const PetscScalar *label_values;
74a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContextDoubleRead(op, label, &num_elements, &label_values));
7515637395SJames Wright     PetscCheck(num_elements == 1, comm, PETSC_ERR_SUP, "%s does not support labels with more than 1 value. Label has %zu values", __func__,
7615637395SJames Wright                num_elements);
7715637395SJames Wright     label_value = *label_values;
78a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorRestoreContextDoubleRead(op, label, &label_values));
79e42c1df6SJames Wright   }
8015637395SJames Wright 
8115637395SJames Wright   if (label_value != update_value) {
82a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetContextDouble(op, label, &update_value));
83e42c1df6SJames Wright   }
84ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
85e42c1df6SJames Wright }
86e42c1df6SJames Wright 
8777841947SLeila Ghaffari // RHS (Explicit time-stepper) function setup
8877841947SLeila Ghaffari //   This is the RHS of the ODE, given as u_t = G(t,u)
8977841947SLeila Ghaffari //   This function takes in a state vector Q and writes into G
9077841947SLeila Ghaffari PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) {
9177841947SLeila Ghaffari   User        user = *(User *)user_data;
9215637395SJames Wright   MPI_Comm    comm = PetscObjectComm((PetscObject)ts);
93c798d55aSJames Wright   PetscScalar dt;
949ad5e8e4SJames Wright   Vec         Q_loc = user->Q_loc;
9577841947SLeila Ghaffari 
96f17d818dSJames Wright   PetscFunctionBeginUser;
975e82a6e1SJeremy L Thompson   // Update time dependent data
98a0b9a424SJames Wright   PetscCall(UpdateBoundaryValues(user, Q_loc, t));
99a424bcd0SJames Wright   if (user->phys->solution_time_label) PetscCall(UpdateContextLabel(user->ceed, comm, t, user->op_rhs_ctx->op, user->phys->solution_time_label));
1002b730f8bSJeremy L Thompson   PetscCall(TSGetTimeStep(ts, &dt));
101a424bcd0SJames Wright   if (user->phys->timestep_size_label) PetscCall(UpdateContextLabel(user->ceed, comm, dt, user->op_rhs_ctx->op, user->phys->timestep_size_label));
10277841947SLeila Ghaffari 
1039ad5e8e4SJames Wright   PetscCall(ApplyCeedOperatorGlobalToGlobal(Q, G, user->op_rhs_ctx));
10477841947SLeila Ghaffari 
105d6e67e47SJames Wright   // Inverse of the lumped mass matrix
106186595e6SJames Wright   PetscCall(VecPointwiseMult(G, G, user->M_inv));
107ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
10877841947SLeila Ghaffari }
10977841947SLeila Ghaffari 
110ca69d878SAdeleke O. Bankole // Surface forces function setup
111ca69d878SAdeleke O. Bankole static PetscErrorCode Surface_Forces_NS(DM dm, Vec G_loc, PetscInt num_walls, const PetscInt walls[], PetscScalar *reaction_force) {
112ca69d878SAdeleke O. Bankole   DMLabel            face_label;
113ca69d878SAdeleke O. Bankole   const PetscScalar *g;
114d6734f85SAdeleke O. Bankole   PetscInt           dof, dim = 3;
115ca69d878SAdeleke O. Bankole   MPI_Comm           comm;
116d6734f85SAdeleke O. Bankole   PetscSection       s;
117ca69d878SAdeleke O. Bankole 
118ca69d878SAdeleke O. Bankole   PetscFunctionBeginUser;
119ca69d878SAdeleke O. Bankole   PetscCall(PetscArrayzero(reaction_force, num_walls * dim));
120ca69d878SAdeleke O. Bankole   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
121ca69d878SAdeleke O. Bankole   PetscCall(DMGetLabel(dm, "Face Sets", &face_label));
122ca69d878SAdeleke O. Bankole   PetscCall(VecGetArrayRead(G_loc, &g));
123ca69d878SAdeleke O. Bankole   for (PetscInt w = 0; w < num_walls; w++) {
124ca69d878SAdeleke O. Bankole     const PetscInt wall = walls[w];
125ca69d878SAdeleke O. Bankole     IS             wall_is;
126d6734f85SAdeleke O. Bankole     PetscCall(DMGetLocalSection(dm, &s));
127ca69d878SAdeleke O. Bankole     PetscCall(DMLabelGetStratumIS(face_label, wall, &wall_is));
128ca69d878SAdeleke O. Bankole     if (wall_is) {  // There exist such points on this process
129ca69d878SAdeleke O. Bankole       PetscInt        num_points;
130d6734f85SAdeleke O. Bankole       PetscInt        num_comp = 0;
131ca69d878SAdeleke O. Bankole       const PetscInt *points;
132d6734f85SAdeleke O. Bankole       PetscCall(PetscSectionGetFieldComponents(s, 0, &num_comp));
133ca69d878SAdeleke O. Bankole       PetscCall(ISGetSize(wall_is, &num_points));
134ca69d878SAdeleke O. Bankole       PetscCall(ISGetIndices(wall_is, &points));
135ca69d878SAdeleke O. Bankole       for (PetscInt i = 0; i < num_points; i++) {
136ca69d878SAdeleke O. Bankole         const PetscInt           p = points[i];
137ca69d878SAdeleke O. Bankole         const StateConservative *r;
138ca69d878SAdeleke O. Bankole         PetscCall(DMPlexPointLocalRead(dm, p, g, &r));
139d6734f85SAdeleke O. Bankole         PetscCall(PetscSectionGetDof(s, p, &dof));
140d6734f85SAdeleke O. Bankole         for (PetscInt node = 0; node < dof / num_comp; node++) {
141ca69d878SAdeleke O. Bankole           for (PetscInt j = 0; j < 3; j++) {
142d6734f85SAdeleke O. Bankole             reaction_force[w * dim + j] -= r[node].momentum[j];
143d6734f85SAdeleke O. Bankole           }
144ca69d878SAdeleke O. Bankole         }
145ca69d878SAdeleke O. Bankole       }
146ca69d878SAdeleke O. Bankole       PetscCall(ISRestoreIndices(wall_is, &points));
147ca69d878SAdeleke O. Bankole     }
148ca69d878SAdeleke O. Bankole     PetscCall(ISDestroy(&wall_is));
149ca69d878SAdeleke O. Bankole   }
150ca69d878SAdeleke O. Bankole   PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, reaction_force, dim * num_walls, MPIU_SCALAR, MPI_SUM, comm));
151ca69d878SAdeleke O. Bankole   //  Restore Vectors
152ca69d878SAdeleke O. Bankole   PetscCall(VecRestoreArrayRead(G_loc, &g));
153ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
154ca69d878SAdeleke O. Bankole }
155ca69d878SAdeleke O. Bankole 
15677841947SLeila Ghaffari // Implicit time-stepper function setup
1572b730f8bSJeremy L Thompson PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) {
15877841947SLeila Ghaffari   User         user = *(User *)user_data;
15915637395SJames Wright   MPI_Comm     comm = PetscObjectComm((PetscObject)ts);
160c798d55aSJames Wright   PetscScalar  dt;
1615e82a6e1SJeremy L Thompson   Vec          Q_loc = user->Q_loc, Q_dot_loc = user->Q_dot_loc, G_loc;
16277841947SLeila Ghaffari   PetscMemType q_mem_type, q_dot_mem_type, g_mem_type;
16377841947SLeila Ghaffari 
164f17d818dSJames Wright   PetscFunctionBeginUser;
1655e82a6e1SJeremy L Thompson   // Get local vectors
166ca69d878SAdeleke O. Bankole   PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
1675e82a6e1SJeremy L Thompson 
1685e82a6e1SJeremy L Thompson   // Update time dependent data
169a0b9a424SJames Wright   PetscCall(UpdateBoundaryValues(user, Q_loc, t));
170a424bcd0SJames Wright   if (user->phys->solution_time_label) PetscCall(UpdateContextLabel(user->ceed, comm, t, user->op_ifunction, user->phys->solution_time_label));
1712b730f8bSJeremy L Thompson   PetscCall(TSGetTimeStep(ts, &dt));
172a424bcd0SJames Wright   if (user->phys->timestep_size_label) PetscCall(UpdateContextLabel(user->ceed, comm, dt, user->op_ifunction, user->phys->timestep_size_label));
17377841947SLeila Ghaffari 
17477841947SLeila Ghaffari   // Global-to-local
175878eb0e7SJames Wright   PetscCall(DMGlobalToLocalBegin(user->dm, Q, INSERT_VALUES, Q_loc));
176878eb0e7SJames Wright   PetscCall(DMGlobalToLocalBegin(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
177878eb0e7SJames Wright   PetscCall(DMGlobalToLocalEnd(user->dm, Q, INSERT_VALUES, Q_loc));
178878eb0e7SJames Wright   PetscCall(DMGlobalToLocalEnd(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
17977841947SLeila Ghaffari 
18077841947SLeila Ghaffari   // Place PETSc vectors in CEED vectors
181c798d55aSJames Wright   PetscCall(VecReadP2C(Q_loc, &q_mem_type, user->q_ceed));
182c798d55aSJames Wright   PetscCall(VecReadP2C(Q_dot_loc, &q_dot_mem_type, user->q_dot_ceed));
183c798d55aSJames Wright   PetscCall(VecP2C(G_loc, &g_mem_type, user->g_ceed));
18477841947SLeila Ghaffari 
18577841947SLeila Ghaffari   // Apply CEED operator
18675d1798cSJames Wright   PetscCall(PetscLogEventBegin(FLUIDS_CeedOperatorApply, Q, G, 0, 0));
18775d1798cSJames Wright   PetscCall(PetscLogGpuTimeBegin());
188a424bcd0SJames Wright   PetscCallCeed(user->ceed, CeedOperatorApply(user->op_ifunction, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE));
18975d1798cSJames Wright   PetscCall(PetscLogGpuTimeEnd());
19075d1798cSJames Wright   PetscCall(PetscLogEventEnd(FLUIDS_CeedOperatorApply, Q, G, 0, 0));
19177841947SLeila Ghaffari 
19277841947SLeila Ghaffari   // Restore vectors
193c798d55aSJames Wright   PetscCall(VecReadC2P(user->q_ceed, q_mem_type, Q_loc));
194c798d55aSJames Wright   PetscCall(VecReadC2P(user->q_dot_ceed, q_dot_mem_type, Q_dot_loc));
195c798d55aSJames Wright   PetscCall(VecC2P(user->g_ceed, g_mem_type, G_loc));
19677841947SLeila Ghaffari 
19719ffbc25SJames Wright   if (user->app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) {
198cbef7084SJames Wright     PetscCall(SgsDDModelApplyIFunction(user, Q_loc, G_loc));
19919ffbc25SJames Wright   }
200be34b3b7SJames Wright 
20177841947SLeila Ghaffari   // Local-to-Global
2022b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(G));
2032b730f8bSJeremy L Thompson   PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G));
20477841947SLeila Ghaffari 
20577841947SLeila Ghaffari   // Restore vectors
206ca69d878SAdeleke O. Bankole   PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
207ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
20877841947SLeila Ghaffari }
20977841947SLeila Ghaffari 
2102b730f8bSJeremy L Thompson static PetscErrorCode FormPreallocation(User user, PetscBool pbdiagonal, Mat J, CeedVector *coo_values) {
211544be873SJed Brown   PetscCount ncoo;
212457e73b2SJames Wright   PetscInt  *rows_petsc, *cols_petsc;
213*01f0e615SJames Wright   CeedInt   *rows_ceed, *cols_ceed;
214544be873SJed Brown 
215544be873SJed Brown   PetscFunctionBeginUser;
216544be873SJed Brown   if (pbdiagonal) {
217*01f0e615SJames Wright     PetscCallCeed(user->ceed, CeedOperatorLinearAssemblePointBlockDiagonalSymbolic(user->op_ijacobian, &ncoo, &rows_ceed, &cols_ceed));
218544be873SJed Brown   } else {
219a424bcd0SJames Wright     PetscCallCeed(user->ceed, CeedOperatorLinearAssembleSymbolic(user->op_ijacobian, &ncoo, &rows_ceed, &cols_ceed));
220*01f0e615SJames Wright   }
221457e73b2SJames Wright   PetscCall(IntArrayC2P(ncoo, &rows_ceed, &rows_petsc));
222457e73b2SJames Wright   PetscCall(IntArrayC2P(ncoo, &cols_ceed, &cols_petsc));
223457e73b2SJames Wright   PetscCall(MatSetPreallocationCOOLocal(J, ncoo, rows_petsc, cols_petsc));
224457e73b2SJames Wright   free(rows_petsc);
225457e73b2SJames Wright   free(cols_petsc);
226a424bcd0SJames Wright   PetscCallCeed(user->ceed, CeedVectorCreate(user->ceed, ncoo, coo_values));
227ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
228544be873SJed Brown }
229544be873SJed Brown 
2302b730f8bSJeremy L Thompson static PetscErrorCode FormSetValues(User user, PetscBool pbdiagonal, Mat J, CeedVector coo_values) {
231544be873SJed Brown   CeedMemType        mem_type = CEED_MEM_HOST;
232544be873SJed Brown   const PetscScalar *values;
233544be873SJed Brown   MatType            mat_type;
234544be873SJed Brown 
235544be873SJed Brown   PetscFunctionBeginUser;
236544be873SJed Brown   PetscCall(MatGetType(J, &mat_type));
2372b730f8bSJeremy L Thompson   if (strstr(mat_type, "kokkos") || strstr(mat_type, "cusparse")) mem_type = CEED_MEM_DEVICE;
238fbb546ddSJames Wright   if (pbdiagonal) {
23975d1798cSJames Wright     PetscCall(PetscLogEventBegin(FLUIDS_CeedOperatorAssemblePointBlockDiagonal, J, 0, 0, 0));
24075d1798cSJames Wright     PetscCall(PetscLogGpuTimeBegin());
241a424bcd0SJames Wright     PetscCallCeed(user->ceed, CeedOperatorLinearAssemblePointBlockDiagonal(user->op_ijacobian, coo_values, CEED_REQUEST_IMMEDIATE));
24275d1798cSJames Wright     PetscCall(PetscLogGpuTimeEnd());
24375d1798cSJames Wright     PetscCall(PetscLogEventEnd(FLUIDS_CeedOperatorAssemblePointBlockDiagonal, J, 0, 0, 0));
244544be873SJed Brown   } else {
24575d1798cSJames Wright     PetscCall(PetscLogEventBegin(FLUIDS_CeedOperatorAssemble, J, 0, 0, 0));
24675d1798cSJames Wright     PetscCall(PetscLogGpuTimeBegin());
247a424bcd0SJames Wright     PetscCallCeed(user->ceed, CeedOperatorLinearAssemble(user->op_ijacobian, coo_values));
24875d1798cSJames Wright     PetscCall(PetscLogGpuTimeEnd());
24975d1798cSJames Wright     PetscCall(PetscLogEventEnd(FLUIDS_CeedOperatorAssemble, J, 0, 0, 0));
250544be873SJed Brown   }
251a424bcd0SJames Wright   PetscCallCeed(user->ceed, CeedVectorGetArrayRead(coo_values, mem_type, &values));
252544be873SJed Brown   PetscCall(MatSetValuesCOO(J, values, INSERT_VALUES));
253a424bcd0SJames Wright   PetscCallCeed(user->ceed, CeedVectorRestoreArrayRead(coo_values, &values));
254ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
255544be873SJed Brown }
256544be873SJed Brown 
2572b730f8bSJeremy L Thompson PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, PetscReal shift, Mat J, Mat J_pre, void *user_data) {
258e334ad8fSJed Brown   User      user = *(User *)user_data;
259a424bcd0SJames Wright   Ceed      ceed = user->ceed;
260d6bf345cSJed Brown   PetscBool J_is_shell, J_is_mffd, J_pre_is_shell;
261f17d818dSJames Wright 
262e334ad8fSJed Brown   PetscFunctionBeginUser;
263a424bcd0SJames Wright   if (user->phys->ijacobian_time_shift_label)
264a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_ijacobian, user->phys->ijacobian_time_shift_label, &shift));
265d6bf345cSJed Brown   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATMFFD, &J_is_mffd));
266e334ad8fSJed Brown   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATSHELL, &J_is_shell));
2672b730f8bSJeremy L Thompson   PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATSHELL, &J_pre_is_shell));
268e334ad8fSJed Brown   if (!user->matrices_set_up) {
269e334ad8fSJed Brown     if (J_is_shell) {
270b07979f9SJames Wright       OperatorApplyContext op_ijacobian_ctx;
271b07979f9SJames Wright       OperatorApplyContextCreate(user->dm, user->dm, user->ceed, user->op_ijacobian, user->q_ceed, user->g_ceed, user->Q_dot_loc, NULL,
272b07979f9SJames Wright                                  &op_ijacobian_ctx);
273b07979f9SJames Wright       PetscCall(MatShellSetContext(J, op_ijacobian_ctx));
274b07979f9SJames Wright       PetscCall(MatShellSetContextDestroy(J, (PetscErrorCode(*)(void *))OperatorApplyContextDestroy));
275b07979f9SJames Wright       PetscCall(MatShellSetOperation(J, MATOP_MULT, (void (*)(void))MatMult_Ceed));
276b07979f9SJames Wright       PetscCall(MatShellSetOperation(J, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiag_Ceed));
277e334ad8fSJed Brown       PetscCall(MatSetUp(J));
278e334ad8fSJed Brown     }
279e334ad8fSJed Brown     if (!J_pre_is_shell) {
2802b730f8bSJeremy L Thompson       PetscCall(FormPreallocation(user, user->app_ctx->pmat_pbdiagonal, J_pre, &user->coo_values_pmat));
281544be873SJed Brown     }
282d6bf345cSJed Brown     if (J != J_pre && !J_is_shell && !J_is_mffd) {
283544be873SJed Brown       PetscCall(FormPreallocation(user, PETSC_FALSE, J, &user->coo_values_amat));
284544be873SJed Brown     }
285e334ad8fSJed Brown     user->matrices_set_up = true;
286e334ad8fSJed Brown   }
287e334ad8fSJed Brown   if (!J_pre_is_shell) {
2882b730f8bSJeremy L Thompson     PetscCall(FormSetValues(user, user->app_ctx->pmat_pbdiagonal, J_pre, user->coo_values_pmat));
289e334ad8fSJed Brown   }
290d6bf345cSJed Brown   if (user->coo_values_amat) {
291d6bf345cSJed Brown     PetscCall(FormSetValues(user, PETSC_FALSE, J, user->coo_values_amat));
292d6bf345cSJed Brown   } else if (J_is_mffd) {
293d6bf345cSJed Brown     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
294d6bf345cSJed Brown     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
295d6bf345cSJed Brown   }
296ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
297e334ad8fSJed Brown }
298e334ad8fSJed Brown 
2992b730f8bSJeremy L Thompson PetscErrorCode WriteOutput(User user, Vec Q, PetscInt step_no, PetscScalar time) {
30077841947SLeila Ghaffari   Vec         Q_loc;
30177841947SLeila Ghaffari   char        file_path[PETSC_MAX_PATH_LEN];
30277841947SLeila Ghaffari   PetscViewer viewer;
30377841947SLeila Ghaffari 
304f17d818dSJames Wright   PetscFunctionBeginUser;
30537cbb16aSJed Brown   if (user->app_ctx->checkpoint_vtk) {
30677841947SLeila Ghaffari     // Set up output
307f14660a4SJames Wright     PetscCall(DMGetLocalVector(user->dm, &Q_loc));
308f14660a4SJames Wright     PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec"));
309f14660a4SJames Wright     PetscCall(VecZeroEntries(Q_loc));
310f14660a4SJames Wright     PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
31177841947SLeila Ghaffari 
31277841947SLeila Ghaffari     // Output
31337cbb16aSJed Brown     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
314f14660a4SJames Wright 
3152b730f8bSJeremy L Thompson     PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path, FILE_MODE_WRITE, &viewer));
316f14660a4SJames Wright     PetscCall(VecView(Q_loc, viewer));
317f14660a4SJames Wright     PetscCall(PetscViewerDestroy(&viewer));
31877841947SLeila Ghaffari     if (user->dm_viz) {
31977841947SLeila Ghaffari       Vec         Q_refined, Q_refined_loc;
32077841947SLeila Ghaffari       char        file_path_refined[PETSC_MAX_PATH_LEN];
32177841947SLeila Ghaffari       PetscViewer viewer_refined;
32277841947SLeila Ghaffari 
323f14660a4SJames Wright       PetscCall(DMGetGlobalVector(user->dm_viz, &Q_refined));
324f14660a4SJames Wright       PetscCall(DMGetLocalVector(user->dm_viz, &Q_refined_loc));
325f14660a4SJames Wright       PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined"));
326f14660a4SJames Wright 
327f14660a4SJames Wright       PetscCall(MatInterpolate(user->interp_viz, Q, Q_refined));
328f14660a4SJames Wright       PetscCall(VecZeroEntries(Q_refined_loc));
3292b730f8bSJeremy L Thompson       PetscCall(DMGlobalToLocal(user->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc));
330f14660a4SJames Wright 
33137cbb16aSJed Brown       PetscCall(
33237cbb16aSJed Brown           PetscSNPrintf(file_path_refined, sizeof file_path_refined, "%s/nsrefined-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
333f14660a4SJames Wright 
3342b730f8bSJeremy L Thompson       PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined), file_path_refined, FILE_MODE_WRITE, &viewer_refined));
335f14660a4SJames Wright       PetscCall(VecView(Q_refined_loc, viewer_refined));
336f14660a4SJames Wright       PetscCall(DMRestoreLocalVector(user->dm_viz, &Q_refined_loc));
337f14660a4SJames Wright       PetscCall(DMRestoreGlobalVector(user->dm_viz, &Q_refined));
338f14660a4SJames Wright       PetscCall(PetscViewerDestroy(&viewer_refined));
33977841947SLeila Ghaffari     }
340f14660a4SJames Wright     PetscCall(DMRestoreLocalVector(user->dm, &Q_loc));
34137cbb16aSJed Brown   }
34277841947SLeila Ghaffari 
34377841947SLeila Ghaffari   // Save data in a binary file for continuation of simulations
34469293791SJames Wright   if (user->app_ctx->add_stepnum2bin) {
34537cbb16aSJed Brown     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", user->app_ctx->output_dir, step_no));
34669293791SJames Wright   } else {
3472b730f8bSJeremy L Thompson     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", user->app_ctx->output_dir));
34869293791SJames Wright   }
3492b730f8bSJeremy L Thompson   PetscCall(PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, &viewer));
350f14660a4SJames Wright 
3510de6a49fSJames Wright   PetscInt32 token = PetscDefined(USE_64BIT_INDICES) ? FLUIDS_FILE_TOKEN_64 : FLUIDS_FILE_TOKEN_32;
3520de6a49fSJames Wright   PetscCall(PetscViewerBinaryWrite(viewer, &token, 1, PETSC_INT32));
3534de8550aSJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &step_no, 1, PETSC_INT));
3544de8550aSJed Brown   time /= user->units->second;  // Dimensionalize time back
3554de8550aSJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL));
356f14660a4SJames Wright   PetscCall(VecView(Q, viewer));
357f14660a4SJames Wright   PetscCall(PetscViewerDestroy(&viewer));
358ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
359f14660a4SJames Wright }
360f14660a4SJames Wright 
361ca69d878SAdeleke O. Bankole // CSV Monitor
362ca69d878SAdeleke O. Bankole PetscErrorCode TSMonitor_WallForce(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
363ca69d878SAdeleke O. Bankole   User              user = ctx;
364ca69d878SAdeleke O. Bankole   Vec               G_loc;
365ca69d878SAdeleke O. Bankole   PetscInt          num_wall = user->app_ctx->wall_forces.num_wall, dim = 3;
366ca69d878SAdeleke O. Bankole   const PetscInt   *walls  = user->app_ctx->wall_forces.walls;
367ca69d878SAdeleke O. Bankole   PetscViewer       viewer = user->app_ctx->wall_forces.viewer;
368ca69d878SAdeleke O. Bankole   PetscViewerFormat format = user->app_ctx->wall_forces.viewer_format;
369ca69d878SAdeleke O. Bankole   PetscScalar      *reaction_force;
370ca69d878SAdeleke O. Bankole   PetscBool         iascii;
371ca69d878SAdeleke O. Bankole 
372ca69d878SAdeleke O. Bankole   PetscFunctionBeginUser;
373ee4ca9cbSJames Wright   if (!viewer) PetscFunctionReturn(PETSC_SUCCESS);
374ca69d878SAdeleke O. Bankole   PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
375ca69d878SAdeleke O. Bankole   PetscCall(PetscMalloc1(num_wall * dim, &reaction_force));
376ca69d878SAdeleke O. Bankole   PetscCall(Surface_Forces_NS(user->dm, G_loc, num_wall, walls, reaction_force));
377ca69d878SAdeleke O. Bankole   PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
378ca69d878SAdeleke O. Bankole 
379ca69d878SAdeleke O. Bankole   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
380ca69d878SAdeleke O. Bankole 
381ca69d878SAdeleke O. Bankole   if (iascii) {
382ca69d878SAdeleke O. Bankole     if (format == PETSC_VIEWER_ASCII_CSV && !user->app_ctx->wall_forces.header_written) {
383ca69d878SAdeleke O. Bankole       PetscCall(PetscViewerASCIIPrintf(viewer, "Step,Time,Wall,ForceX,ForceY,ForceZ\n"));
384ca69d878SAdeleke O. Bankole       user->app_ctx->wall_forces.header_written = PETSC_TRUE;
385ca69d878SAdeleke O. Bankole     }
386ca69d878SAdeleke O. Bankole     for (PetscInt w = 0; w < num_wall; w++) {
387ca69d878SAdeleke O. Bankole       PetscInt wall = walls[w];
388ca69d878SAdeleke O. Bankole       if (format == PETSC_VIEWER_ASCII_CSV) {
389ca69d878SAdeleke O. Bankole         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT ",%g,%" PetscInt_FMT ",%g,%g,%g\n", step_no, time, wall,
390ca69d878SAdeleke O. Bankole                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
391ca69d878SAdeleke O. Bankole 
392ca69d878SAdeleke O. Bankole       } else {
393ca69d878SAdeleke O. Bankole         PetscCall(PetscViewerASCIIPrintf(viewer, "Wall %" PetscInt_FMT " Forces: Force_x = %12g, Force_y = %12g, Force_z = %12g\n", wall,
394ca69d878SAdeleke O. Bankole                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
395ca69d878SAdeleke O. Bankole       }
396ca69d878SAdeleke O. Bankole     }
397ca69d878SAdeleke O. Bankole   }
398ca69d878SAdeleke O. Bankole   PetscCall(PetscFree(reaction_force));
399ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
400ca69d878SAdeleke O. Bankole }
401ca69d878SAdeleke O. Bankole 
402f14660a4SJames Wright // User provided TS Monitor
4032b730f8bSJeremy L Thompson PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
404f14660a4SJames Wright   User user = ctx;
405f14660a4SJames Wright 
406f17d818dSJames Wright   PetscFunctionBeginUser;
40737cbb16aSJed Brown   // Print every 'checkpoint_interval' steps
408894de27cSJames Wright   if (user->app_ctx->checkpoint_interval <= 0 || step_no % user->app_ctx->checkpoint_interval != 0 ||
40949aac155SJeremy L Thompson       (user->app_ctx->cont_steps == step_no && step_no != 0)) {
410ee4ca9cbSJames Wright     PetscFunctionReturn(PETSC_SUCCESS);
41149aac155SJeremy L Thompson   }
412f14660a4SJames Wright 
413f14660a4SJames Wright   PetscCall(WriteOutput(user, Q, step_no, time));
414ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
41577841947SLeila Ghaffari }
41677841947SLeila Ghaffari 
41777841947SLeila Ghaffari // TS: Create, setup, and solve
4182b730f8bSJeremy L Thompson PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, Vec *Q, PetscScalar *f_time, TS *ts) {
41977841947SLeila Ghaffari   MPI_Comm    comm = user->comm;
42077841947SLeila Ghaffari   TSAdapt     adapt;
42177841947SLeila Ghaffari   PetscScalar final_time;
42277841947SLeila Ghaffari 
423f17d818dSJames Wright   PetscFunctionBeginUser;
4242b730f8bSJeremy L Thompson   PetscCall(TSCreate(comm, ts));
4252b730f8bSJeremy L Thompson   PetscCall(TSSetDM(*ts, dm));
42677841947SLeila Ghaffari   if (phys->implicit) {
4272b730f8bSJeremy L Thompson     PetscCall(TSSetType(*ts, TSBDF));
42877841947SLeila Ghaffari     if (user->op_ifunction) {
4292b730f8bSJeremy L Thompson       PetscCall(TSSetIFunction(*ts, NULL, IFunction_NS, &user));
43077841947SLeila Ghaffari     } else {  // Implicit integrators can fall back to using an RHSFunction
4312b730f8bSJeremy L Thompson       PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
43277841947SLeila Ghaffari     }
433e334ad8fSJed Brown     if (user->op_ijacobian) {
4342b730f8bSJeremy L Thompson       PetscCall(DMTSSetIJacobian(dm, FormIJacobian_NS, &user));
435544be873SJed Brown       if (app_ctx->amat_type) {
436544be873SJed Brown         Mat Pmat, Amat;
4372b730f8bSJeremy L Thompson         PetscCall(DMCreateMatrix(dm, &Pmat));
4382b730f8bSJeremy L Thompson         PetscCall(DMSetMatType(dm, app_ctx->amat_type));
4392b730f8bSJeremy L Thompson         PetscCall(DMCreateMatrix(dm, &Amat));
4402b730f8bSJeremy L Thompson         PetscCall(TSSetIJacobian(*ts, Amat, Pmat, NULL, NULL));
4412b730f8bSJeremy L Thompson         PetscCall(MatDestroy(&Amat));
4422b730f8bSJeremy L Thompson         PetscCall(MatDestroy(&Pmat));
443544be873SJed Brown       }
444e334ad8fSJed Brown     }
44577841947SLeila Ghaffari   } else {
4469ad5e8e4SJames Wright     PetscCheck(user->op_rhs_ctx, comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction");
4472b730f8bSJeremy L Thompson     PetscCall(TSSetType(*ts, TSRK));
4482b730f8bSJeremy L Thompson     PetscCall(TSRKSetType(*ts, TSRK5F));
4492b730f8bSJeremy L Thompson     PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
45077841947SLeila Ghaffari   }
4512b730f8bSJeremy L Thompson   PetscCall(TSSetMaxTime(*ts, 500. * user->units->second));
4522b730f8bSJeremy L Thompson   PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER));
45351b00d91SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(TSSetErrorIfStepFails(*ts, PETSC_FALSE));
4542b730f8bSJeremy L Thompson   PetscCall(TSSetTimeStep(*ts, 1.e-2 * user->units->second));
4558fb33541SJames Wright   if (app_ctx->test_type != TESTTYPE_NONE) {
4562b730f8bSJeremy L Thompson     PetscCall(TSSetMaxSteps(*ts, 10));
4572b730f8bSJeremy L Thompson   }
4582b730f8bSJeremy L Thompson   PetscCall(TSGetAdapt(*ts, &adapt));
4592b730f8bSJeremy L Thompson   PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * user->units->second, 1.e2 * user->units->second));
4602b730f8bSJeremy L Thompson   PetscCall(TSSetFromOptions(*ts));
46115637395SJames Wright   user->time_bc_set = -1.0;   // require all BCs be updated
462d55646a4SJames Wright   if (app_ctx->cont_steps) {  // continue from previous timestep data
46377841947SLeila Ghaffari     PetscInt    count;
46477841947SLeila Ghaffari     PetscViewer viewer;
4652b730f8bSJeremy L Thompson 
4664de8550aSJed Brown     if (app_ctx->cont_time <= 0) {  // Legacy files did not include step number and time
4672b730f8bSJeremy L Thompson       PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_time_file, FILE_MODE_READ, &viewer));
4684de8550aSJed Brown       PetscCall(PetscViewerBinaryRead(viewer, &app_ctx->cont_time, 1, &count, PETSC_REAL));
4692b730f8bSJeremy L Thompson       PetscCall(PetscViewerDestroy(&viewer));
4704de8550aSJed Brown       PetscCheck(app_ctx->cont_steps != -1, comm, PETSC_ERR_ARG_INCOMP,
4714de8550aSJed Brown                  "-continue step number not specified, but checkpoint file does not contain a step number (likely written by older code version)");
4724de8550aSJed Brown     }
4734de8550aSJed Brown     PetscCall(TSSetTime(*ts, app_ctx->cont_time * user->units->second));
474d8e0aecdSJed Brown     PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps));
47577841947SLeila Ghaffari   }
4768fb33541SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) {
4772b730f8bSJeremy L Thompson     PetscCall(TSMonitorSet(*ts, TSMonitor_NS, user, NULL));
4788fb33541SJames Wright   }
479ca69d878SAdeleke O. Bankole   if (app_ctx->wall_forces.viewer) {
480ca69d878SAdeleke O. Bankole     PetscCall(TSMonitorSet(*ts, TSMonitor_WallForce, user, NULL));
481ca69d878SAdeleke O. Bankole   }
482b7d66439SJames Wright   if (app_ctx->turb_spanstats_enable) {
483f5452247SJames Wright     PetscCall(TSMonitorSet(*ts, TSMonitor_TurbulenceStatistics, user, NULL));
484495b9769SJames Wright     CeedScalar previous_time = app_ctx->cont_time * user->units->second;
485a424bcd0SJames Wright     PetscCallCeed(user->ceed,
486a424bcd0SJames Wright                   CeedOperatorSetContextDouble(user->spanstats.op_stats_collect_ctx->op, user->spanstats.previous_time_label, &previous_time));
487a175e481SJames Wright   }
4884e9802d1SJames Wright   if (app_ctx->diff_filter_monitor) PetscCall(TSMonitorSet(*ts, TSMonitor_DifferentialFilter, user, NULL));
48977841947SLeila Ghaffari 
49077841947SLeila Ghaffari   // Solve
491d8e0aecdSJed Brown   PetscReal start_time;
492d8e0aecdSJed Brown   PetscInt  start_step;
4932b730f8bSJeremy L Thompson   PetscCall(TSGetTime(*ts, &start_time));
494d8e0aecdSJed Brown   PetscCall(TSGetStepNumber(*ts, &start_step));
4957b39487dSJeremy L Thompson 
496711a423aSJed Brown   PetscCall(PetscLogDefaultBegin());  // So we can use PetscLogStageGetPerfInfo without -log_view
4977b39487dSJeremy L Thompson   PetscPreLoadBegin(PETSC_FALSE, "Fluids Solve");
4987b39487dSJeremy L Thompson   PetscCall(TSSetTime(*ts, start_time));
499d8e0aecdSJed Brown   PetscCall(TSSetStepNumber(*ts, start_step));
5007b39487dSJeremy L Thompson   if (PetscPreLoadingOn) {
5017b39487dSJeremy L Thompson     // LCOV_EXCL_START
5027b39487dSJeremy L Thompson     SNES      snes;
5037b39487dSJeremy L Thompson     Vec       Q_preload;
5047b39487dSJeremy L Thompson     PetscReal rtol;
5057b39487dSJeremy L Thompson     PetscCall(VecDuplicate(*Q, &Q_preload));
5067b39487dSJeremy L Thompson     PetscCall(VecCopy(*Q, Q_preload));
5077b39487dSJeremy L Thompson     PetscCall(TSGetSNES(*ts, &snes));
5087b39487dSJeremy L Thompson     PetscCall(SNESGetTolerances(snes, NULL, &rtol, NULL, NULL, NULL));
5092b730f8bSJeremy L Thompson     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, .99, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
510fc4c3d69SJames Wright     PetscCall(TSSetSolution(*ts, Q_preload));
5117b39487dSJeremy L Thompson     PetscCall(TSStep(*ts));
5122b730f8bSJeremy L Thompson     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, rtol, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
5137b39487dSJeremy L Thompson     PetscCall(VecDestroy(&Q_preload));
5147b39487dSJeremy L Thompson     // LCOV_EXCL_STOP
5157b39487dSJeremy L Thompson   } else {
5162b730f8bSJeremy L Thompson     PetscCall(PetscBarrier((PetscObject)*ts));
5172b730f8bSJeremy L Thompson     PetscCall(TSSolve(*ts, *Q));
5187b39487dSJeremy L Thompson   }
5197b39487dSJeremy L Thompson   PetscPreLoadEnd();
5207b39487dSJeremy L Thompson 
5217b39487dSJeremy L Thompson   PetscCall(TSGetSolveTime(*ts, &final_time));
52277841947SLeila Ghaffari   *f_time = final_time;
5237b39487dSJeremy L Thompson 
5248fb33541SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) {
525f14660a4SJames Wright     PetscInt step_no;
526f14660a4SJames Wright     PetscCall(TSGetStepNumber(*ts, &step_no));
527a175e481SJames Wright     if (user->app_ctx->checkpoint_interval > 0 || user->app_ctx->checkpoint_interval == -1) {
528f14660a4SJames Wright       PetscCall(WriteOutput(user, *Q, step_no, final_time));
529f14660a4SJames Wright     }
530f14660a4SJames Wright 
53175d1798cSJames Wright     PetscLogStage      stage_id;
532711a423aSJed Brown     PetscEventPerfInfo stage_perf;
5337b39487dSJeremy L Thompson 
5347b39487dSJeremy L Thompson     PetscCall(PetscLogStageGetId("Fluids Solve", &stage_id));
535711a423aSJed Brown     PetscCall(PetscLogStageGetPerfInfo(stage_id, &stage_perf));
536711a423aSJed Brown     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time taken for solution (sec): %g\n", stage_perf.time));
53777841947SLeila Ghaffari   }
538ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
53977841947SLeila Ghaffari }
540