xref: /honee/src/setupts.c (revision 0143e3da1c9d35d446c92225da280c31ee8e3678)
1727da7e7SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2727da7e7SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3a515125bSLeila Ghaffari //
4727da7e7SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5a515125bSLeila Ghaffari //
6727da7e7SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
7a515125bSLeila Ghaffari 
8a515125bSLeila Ghaffari /// @file
9a515125bSLeila Ghaffari /// Time-stepping functions for Navier-Stokes example using PETSc
10a515125bSLeila Ghaffari 
11e419654dSJeremy L Thompson #include <ceed.h>
12e419654dSJeremy L Thompson #include <petscdmplex.h>
13e419654dSJeremy L Thompson #include <petscts.h>
14e419654dSJeremy L Thompson 
15a515125bSLeila Ghaffari #include "../navierstokes.h"
16c5e9980aSAdeleke O. Bankole #include "../qfunctions/newtonian_state.h"
17a515125bSLeila Ghaffari 
18a515125bSLeila Ghaffari // Compute mass matrix for explicit scheme
192b916ea7SJeremy L Thompson PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm, CeedData ceed_data, Vec M) {
20a515125bSLeila Ghaffari   CeedQFunction        qf_mass;
21a515125bSLeila Ghaffari   CeedOperator         op_mass;
22*0143e3daSJames Wright   OperatorApplyContext op_mass_ctx;
23*0143e3daSJames Wright   Vec                  Ones_loc;
24a515125bSLeila Ghaffari   CeedInt              num_comp_q, q_data_size;
25a515125bSLeila Ghaffari   PetscFunctionBeginUser;
26a515125bSLeila Ghaffari 
27a515125bSLeila Ghaffari   // CEED Restriction
28a515125bSLeila Ghaffari   CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q);
29a515125bSLeila Ghaffari   CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &q_data_size);
30a515125bSLeila Ghaffari 
31a515125bSLeila Ghaffari   // CEED QFunction
329f59f36eSJames Wright   PetscCall(CreateMassQFunction(ceed, num_comp_q, q_data_size, &qf_mass));
33a515125bSLeila Ghaffari 
34a515125bSLeila Ghaffari   // CEED Operator
35a515125bSLeila Ghaffari   CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass);
369f59f36eSJames Wright   CeedOperatorSetField(op_mass, "u", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
372b916ea7SJeremy L Thompson   CeedOperatorSetField(op_mass, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data);
382b916ea7SJeremy L Thompson   CeedOperatorSetField(op_mass, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
39a515125bSLeila Ghaffari 
40*0143e3daSJames Wright   PetscCall(OperatorApplyContextCreate(NULL, dm, ceed, op_mass, NULL, NULL, NULL, NULL, &op_mass_ctx));
41a515125bSLeila Ghaffari 
42*0143e3daSJames Wright   PetscCall(DMGetLocalVector(dm, &Ones_loc));
43*0143e3daSJames Wright   PetscCall(VecSet(Ones_loc, 1));
44*0143e3daSJames Wright   PetscCall(ApplyCeedOperatorLocalToGlobal(Ones_loc, M, op_mass_ctx));
45a515125bSLeila Ghaffari 
46a515125bSLeila Ghaffari   // Invert diagonally lumped mass vector for RHS function
472b916ea7SJeremy L Thompson   PetscCall(VecReciprocal(M));
48a515125bSLeila Ghaffari 
49a515125bSLeila Ghaffari   // Cleanup
50*0143e3daSJames Wright   PetscCall(OperatorApplyContextDestroy(op_mass_ctx));
51*0143e3daSJames Wright   PetscCall(DMRestoreLocalVector(dm, &Ones_loc));
52a515125bSLeila Ghaffari   CeedQFunctionDestroy(&qf_mass);
53a515125bSLeila Ghaffari   CeedOperatorDestroy(&op_mass);
54a515125bSLeila Ghaffari 
55a515125bSLeila Ghaffari   PetscFunctionReturn(0);
56a515125bSLeila Ghaffari }
57a515125bSLeila Ghaffari 
58c996854bSJames Wright // Insert Boundary values if it's a new time
59c996854bSJames Wright PetscErrorCode UpdateBoundaryValues(User user, Vec Q_loc, PetscReal t) {
60c996854bSJames Wright   PetscFunctionBeginUser;
61c996854bSJames Wright   if (user->time_bc_set != t) {
62c996854bSJames Wright     PetscCall(DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL));
63c996854bSJames Wright     user->time_bc_set = t;
64c996854bSJames Wright   }
65c996854bSJames Wright   PetscFunctionReturn(0);
66c996854bSJames Wright }
67c996854bSJames Wright 
6891f639d2SJames Wright // @brief Update the context label value to new value if necessary.
6991f639d2SJames Wright // @note This only supports labels with scalar label values (ie. not arrays)
7091f639d2SJames Wright PetscErrorCode UpdateContextLabel(MPI_Comm comm, PetscScalar update_value, CeedOperator op, CeedContextFieldLabel label) {
7191f639d2SJames Wright   PetscScalar label_value;
72c2cb7fc8SJames Wright 
7391f639d2SJames Wright   PetscFunctionBeginUser;
7491f639d2SJames Wright   PetscCheck(label, comm, PETSC_ERR_ARG_BADPTR, "Label should be non-NULL");
7591f639d2SJames Wright 
7691f639d2SJames Wright   {
7791f639d2SJames Wright     size_t             num_elements;
7891f639d2SJames Wright     const PetscScalar *label_values;
7991f639d2SJames Wright     CeedOperatorGetContextDoubleRead(op, label, &num_elements, &label_values);
8091f639d2SJames Wright     PetscCheck(num_elements == 1, comm, PETSC_ERR_SUP, "%s does not support labels with more than 1 value. Label has %zu values", __func__,
8191f639d2SJames Wright                num_elements);
8291f639d2SJames Wright     label_value = *label_values;
8391f639d2SJames Wright     CeedOperatorRestoreContextDoubleRead(op, label, &label_values);
84c2cb7fc8SJames Wright   }
8591f639d2SJames Wright 
8691f639d2SJames Wright   if (label_value != update_value) {
8791f639d2SJames Wright     CeedOperatorSetContextDouble(op, label, &update_value);
88c2cb7fc8SJames Wright   }
89c2cb7fc8SJames Wright   PetscFunctionReturn(0);
90c2cb7fc8SJames Wright }
91c2cb7fc8SJames Wright 
92a515125bSLeila Ghaffari // RHS (Explicit time-stepper) function setup
93a515125bSLeila Ghaffari //   This is the RHS of the ODE, given as u_t = G(t,u)
94a515125bSLeila Ghaffari //   This function takes in a state vector Q and writes into G
95a515125bSLeila Ghaffari PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) {
96a515125bSLeila Ghaffari   User         user = *(User *)user_data;
9791f639d2SJames Wright   MPI_Comm     comm = PetscObjectComm((PetscObject)ts);
98fd969b44SJames Wright   PetscScalar  dt;
99e2f84137SJeremy L Thompson   Vec          Q_loc = user->Q_loc, G_loc;
100a515125bSLeila Ghaffari   PetscMemType q_mem_type, g_mem_type;
101a515125bSLeila Ghaffari   PetscFunctionBeginUser;
102a515125bSLeila Ghaffari 
103e2f84137SJeremy L Thompson   // Get local vector
1042b916ea7SJeremy L Thompson   PetscCall(DMGetLocalVector(user->dm, &G_loc));
105e2f84137SJeremy L Thompson 
106e2f84137SJeremy L Thompson   // Update time dependent data
107c996854bSJames Wright   PetscCall(UpdateBoundaryValues(user, Q_loc, t));
10891f639d2SJames Wright   if (user->phys->solution_time_label) PetscCall(UpdateContextLabel(comm, t, user->op_rhs, user->phys->solution_time_label));
1092b916ea7SJeremy L Thompson   PetscCall(TSGetTimeStep(ts, &dt));
11091f639d2SJames Wright   if (user->phys->timestep_size_label) PetscCall(UpdateContextLabel(comm, dt, user->op_rhs, user->phys->timestep_size_label));
111a515125bSLeila Ghaffari 
112a515125bSLeila Ghaffari   // Global-to-local
1132b916ea7SJeremy L Thompson   PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
114a515125bSLeila Ghaffari 
115a515125bSLeila Ghaffari   // Place PETSc vectors in CEED vectors
116fd969b44SJames Wright   PetscCall(VecReadP2C(Q_loc, &q_mem_type, user->q_ceed));
117fd969b44SJames Wright   PetscCall(VecP2C(G_loc, &g_mem_type, user->g_ceed));
118a515125bSLeila Ghaffari 
119a515125bSLeila Ghaffari   // Apply CEED operator
1202b916ea7SJeremy L Thompson   CeedOperatorApply(user->op_rhs, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE);
121a515125bSLeila Ghaffari 
122a515125bSLeila Ghaffari   // Restore vectors
123fd969b44SJames Wright   PetscCall(VecReadC2P(user->q_ceed, q_mem_type, Q_loc));
124fd969b44SJames Wright   PetscCall(VecC2P(user->g_ceed, g_mem_type, G_loc));
125a515125bSLeila Ghaffari 
126a515125bSLeila Ghaffari   // Local-to-Global
1272b916ea7SJeremy L Thompson   PetscCall(VecZeroEntries(G));
1282b916ea7SJeremy L Thompson   PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G));
129a515125bSLeila Ghaffari 
130a515125bSLeila Ghaffari   // Inverse of the lumped mass matrix (M is Minv)
1312b916ea7SJeremy L Thompson   PetscCall(VecPointwiseMult(G, G, user->M));
132a515125bSLeila Ghaffari 
133a515125bSLeila Ghaffari   // Restore vectors
1342b916ea7SJeremy L Thompson   PetscCall(DMRestoreLocalVector(user->dm, &G_loc));
135a515125bSLeila Ghaffari 
136a515125bSLeila Ghaffari   PetscFunctionReturn(0);
137a515125bSLeila Ghaffari }
138a515125bSLeila Ghaffari 
139c5e9980aSAdeleke O. Bankole // Surface forces function setup
140c5e9980aSAdeleke O. Bankole static PetscErrorCode Surface_Forces_NS(DM dm, Vec G_loc, PetscInt num_walls, const PetscInt walls[], PetscScalar *reaction_force) {
141c5e9980aSAdeleke O. Bankole   DMLabel            face_label;
142c5e9980aSAdeleke O. Bankole   const PetscScalar *g;
1432004e3acSAdeleke O. Bankole   PetscInt           dof, dim = 3;
144c5e9980aSAdeleke O. Bankole   MPI_Comm           comm;
1452004e3acSAdeleke O. Bankole   PetscSection       s;
146c5e9980aSAdeleke O. Bankole 
147c5e9980aSAdeleke O. Bankole   PetscFunctionBeginUser;
148c5e9980aSAdeleke O. Bankole   PetscCall(PetscArrayzero(reaction_force, num_walls * dim));
149c5e9980aSAdeleke O. Bankole   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
150c5e9980aSAdeleke O. Bankole   PetscCall(DMGetLabel(dm, "Face Sets", &face_label));
151c5e9980aSAdeleke O. Bankole   PetscCall(VecGetArrayRead(G_loc, &g));
152c5e9980aSAdeleke O. Bankole   for (PetscInt w = 0; w < num_walls; w++) {
153c5e9980aSAdeleke O. Bankole     const PetscInt wall = walls[w];
154c5e9980aSAdeleke O. Bankole     IS             wall_is;
1552004e3acSAdeleke O. Bankole     PetscCall(DMGetLocalSection(dm, &s));
156c5e9980aSAdeleke O. Bankole     PetscCall(DMLabelGetStratumIS(face_label, wall, &wall_is));
157c5e9980aSAdeleke O. Bankole     if (wall_is) {  // There exist such points on this process
158c5e9980aSAdeleke O. Bankole       PetscInt        num_points;
1592004e3acSAdeleke O. Bankole       PetscInt        num_comp = 0;
160c5e9980aSAdeleke O. Bankole       const PetscInt *points;
1612004e3acSAdeleke O. Bankole       PetscCall(PetscSectionGetFieldComponents(s, 0, &num_comp));
162c5e9980aSAdeleke O. Bankole       PetscCall(ISGetSize(wall_is, &num_points));
163c5e9980aSAdeleke O. Bankole       PetscCall(ISGetIndices(wall_is, &points));
164c5e9980aSAdeleke O. Bankole       for (PetscInt i = 0; i < num_points; i++) {
165c5e9980aSAdeleke O. Bankole         const PetscInt           p = points[i];
166c5e9980aSAdeleke O. Bankole         const StateConservative *r;
167c5e9980aSAdeleke O. Bankole         PetscCall(DMPlexPointLocalRead(dm, p, g, &r));
1682004e3acSAdeleke O. Bankole         PetscCall(PetscSectionGetDof(s, p, &dof));
1692004e3acSAdeleke O. Bankole         for (PetscInt node = 0; node < dof / num_comp; node++) {
170c5e9980aSAdeleke O. Bankole           for (PetscInt j = 0; j < 3; j++) {
1712004e3acSAdeleke O. Bankole             reaction_force[w * dim + j] -= r[node].momentum[j];
1722004e3acSAdeleke O. Bankole           }
173c5e9980aSAdeleke O. Bankole         }
174c5e9980aSAdeleke O. Bankole       }
175c5e9980aSAdeleke O. Bankole       PetscCall(ISRestoreIndices(wall_is, &points));
176c5e9980aSAdeleke O. Bankole     }
177c5e9980aSAdeleke O. Bankole     PetscCall(ISDestroy(&wall_is));
178c5e9980aSAdeleke O. Bankole   }
179c5e9980aSAdeleke O. Bankole   PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, reaction_force, dim * num_walls, MPIU_SCALAR, MPI_SUM, comm));
180c5e9980aSAdeleke O. Bankole   //  Restore Vectors
181c5e9980aSAdeleke O. Bankole   PetscCall(VecRestoreArrayRead(G_loc, &g));
182c5e9980aSAdeleke O. Bankole 
183c5e9980aSAdeleke O. Bankole   PetscFunctionReturn(0);
184c5e9980aSAdeleke O. Bankole }
185c5e9980aSAdeleke O. Bankole 
186a515125bSLeila Ghaffari // Implicit time-stepper function setup
1872b916ea7SJeremy L Thompson PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) {
188a515125bSLeila Ghaffari   User         user = *(User *)user_data;
18991f639d2SJames Wright   MPI_Comm     comm = PetscObjectComm((PetscObject)ts);
190fd969b44SJames Wright   PetscScalar  dt;
191e2f84137SJeremy L Thompson   Vec          Q_loc = user->Q_loc, Q_dot_loc = user->Q_dot_loc, G_loc;
192a515125bSLeila Ghaffari   PetscMemType q_mem_type, q_dot_mem_type, g_mem_type;
193a515125bSLeila Ghaffari   PetscFunctionBeginUser;
194a515125bSLeila Ghaffari 
195e2f84137SJeremy L Thompson   // Get local vectors
196c5e9980aSAdeleke O. Bankole   PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
197e2f84137SJeremy L Thompson 
198e2f84137SJeremy L Thompson   // Update time dependent data
199c996854bSJames Wright   PetscCall(UpdateBoundaryValues(user, Q_loc, t));
20091f639d2SJames Wright   if (user->phys->solution_time_label) PetscCall(UpdateContextLabel(comm, t, user->op_ifunction, user->phys->solution_time_label));
2012b916ea7SJeremy L Thompson   PetscCall(TSGetTimeStep(ts, &dt));
20291f639d2SJames Wright   if (user->phys->timestep_size_label) PetscCall(UpdateContextLabel(comm, dt, user->op_ifunction, user->phys->timestep_size_label));
203a515125bSLeila Ghaffari 
204a515125bSLeila Ghaffari   // Global-to-local
20506108310SJames Wright   PetscCall(DMGlobalToLocalBegin(user->dm, Q, INSERT_VALUES, Q_loc));
20606108310SJames Wright   PetscCall(DMGlobalToLocalBegin(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
20706108310SJames Wright   PetscCall(DMGlobalToLocalEnd(user->dm, Q, INSERT_VALUES, Q_loc));
20806108310SJames Wright   PetscCall(DMGlobalToLocalEnd(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
209a515125bSLeila Ghaffari 
210a515125bSLeila Ghaffari   // Place PETSc vectors in CEED vectors
211fd969b44SJames Wright   PetscCall(VecReadP2C(Q_loc, &q_mem_type, user->q_ceed));
212fd969b44SJames Wright   PetscCall(VecReadP2C(Q_dot_loc, &q_dot_mem_type, user->q_dot_ceed));
213fd969b44SJames Wright   PetscCall(VecP2C(G_loc, &g_mem_type, user->g_ceed));
214a515125bSLeila Ghaffari 
215a515125bSLeila Ghaffari   // Apply CEED operator
2162b916ea7SJeremy L Thompson   CeedOperatorApply(user->op_ifunction, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE);
217a515125bSLeila Ghaffari 
218a515125bSLeila Ghaffari   // Restore vectors
219fd969b44SJames Wright   PetscCall(VecReadC2P(user->q_ceed, q_mem_type, Q_loc));
220fd969b44SJames Wright   PetscCall(VecReadC2P(user->q_dot_ceed, q_dot_mem_type, Q_dot_loc));
221fd969b44SJames Wright   PetscCall(VecC2P(user->g_ceed, g_mem_type, G_loc));
222a515125bSLeila Ghaffari 
22301ab89c1SJames Wright   if (user->app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) {
22401ab89c1SJames Wright     PetscCall(SGS_DD_ModelApplyIFunction(user, Q_loc, G_loc));
22501ab89c1SJames Wright   }
2269c678832SJames Wright 
227a515125bSLeila Ghaffari   // Local-to-Global
2282b916ea7SJeremy L Thompson   PetscCall(VecZeroEntries(G));
2292b916ea7SJeremy L Thompson   PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G));
230a515125bSLeila Ghaffari 
231a515125bSLeila Ghaffari   // Restore vectors
232c5e9980aSAdeleke O. Bankole   PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
233a515125bSLeila Ghaffari 
234a515125bSLeila Ghaffari   PetscFunctionReturn(0);
235a515125bSLeila Ghaffari }
236a515125bSLeila Ghaffari 
2372b916ea7SJeremy L Thompson static PetscErrorCode FormPreallocation(User user, PetscBool pbdiagonal, Mat J, CeedVector *coo_values) {
238b107fddaSJed Brown   PetscCount ncoo;
239b107fddaSJed Brown   PetscInt  *rows, *cols;
240b107fddaSJed Brown 
241b107fddaSJed Brown   PetscFunctionBeginUser;
242b107fddaSJed Brown   if (pbdiagonal) {
243b107fddaSJed Brown     CeedSize l_size;
244b107fddaSJed Brown     CeedOperatorGetActiveVectorLengths(user->op_ijacobian, &l_size, NULL);
245b107fddaSJed Brown     ncoo = l_size * 5;
246b107fddaSJed Brown     rows = malloc(ncoo * sizeof(rows[0]));
247b107fddaSJed Brown     cols = malloc(ncoo * sizeof(cols[0]));
248b107fddaSJed Brown     for (PetscCount n = 0; n < l_size / 5; n++) {
249b107fddaSJed Brown       for (PetscInt i = 0; i < 5; i++) {
250b107fddaSJed Brown         for (PetscInt j = 0; j < 5; j++) {
251b107fddaSJed Brown           rows[(n * 5 + i) * 5 + j] = n * 5 + i;
252b107fddaSJed Brown           cols[(n * 5 + i) * 5 + j] = n * 5 + j;
253b107fddaSJed Brown         }
254b107fddaSJed Brown       }
255b107fddaSJed Brown     }
256b107fddaSJed Brown   } else {
2572b916ea7SJeremy L Thompson     PetscCall(CeedOperatorLinearAssembleSymbolic(user->op_ijacobian, &ncoo, &rows, &cols));
258b107fddaSJed Brown   }
259b107fddaSJed Brown   PetscCall(MatSetPreallocationCOOLocal(J, ncoo, rows, cols));
260b107fddaSJed Brown   free(rows);
261b107fddaSJed Brown   free(cols);
262b107fddaSJed Brown   CeedVectorCreate(user->ceed, ncoo, coo_values);
263b107fddaSJed Brown   PetscFunctionReturn(0);
264b107fddaSJed Brown }
265b107fddaSJed Brown 
2662b916ea7SJeremy L Thompson static PetscErrorCode FormSetValues(User user, PetscBool pbdiagonal, Mat J, CeedVector coo_values) {
267b107fddaSJed Brown   CeedMemType        mem_type = CEED_MEM_HOST;
268b107fddaSJed Brown   const PetscScalar *values;
269b107fddaSJed Brown   MatType            mat_type;
270b107fddaSJed Brown 
271b107fddaSJed Brown   PetscFunctionBeginUser;
272b107fddaSJed Brown   PetscCall(MatGetType(J, &mat_type));
2732b916ea7SJeremy L Thompson   if (strstr(mat_type, "kokkos") || strstr(mat_type, "cusparse")) mem_type = CEED_MEM_DEVICE;
274b107fddaSJed Brown   if (user->app_ctx->pmat_pbdiagonal) {
2752b916ea7SJeremy L Thompson     CeedOperatorLinearAssemblePointBlockDiagonal(user->op_ijacobian, coo_values, CEED_REQUEST_IMMEDIATE);
276b107fddaSJed Brown   } else {
277b107fddaSJed Brown     CeedOperatorLinearAssemble(user->op_ijacobian, coo_values);
278b107fddaSJed Brown   }
279b107fddaSJed Brown   CeedVectorGetArrayRead(coo_values, mem_type, &values);
280b107fddaSJed Brown   PetscCall(MatSetValuesCOO(J, values, INSERT_VALUES));
281b107fddaSJed Brown   CeedVectorRestoreArrayRead(coo_values, &values);
282b107fddaSJed Brown   PetscFunctionReturn(0);
283b107fddaSJed Brown }
284b107fddaSJed Brown 
2852b916ea7SJeremy L Thompson PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, PetscReal shift, Mat J, Mat J_pre, void *user_data) {
286f0b65372SJed Brown   User      user = *(User *)user_data;
28704855949SJed Brown   PetscBool J_is_shell, J_is_mffd, J_pre_is_shell;
288f0b65372SJed Brown   PetscFunctionBeginUser;
289a5f46a7bSJeremy L Thompson   if (user->phys->ijacobian_time_shift_label) CeedOperatorSetContextDouble(user->op_ijacobian, user->phys->ijacobian_time_shift_label, &shift);
29004855949SJed Brown   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATMFFD, &J_is_mffd));
291f0b65372SJed Brown   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATSHELL, &J_is_shell));
2922b916ea7SJeremy L Thompson   PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATSHELL, &J_pre_is_shell));
293f0b65372SJed Brown   if (!user->matrices_set_up) {
294f0b65372SJed Brown     if (J_is_shell) {
295f9028c3cSJames Wright       OperatorApplyContext op_ijacobian_ctx;
296f9028c3cSJames Wright       OperatorApplyContextCreate(user->dm, user->dm, user->ceed, user->op_ijacobian, user->q_ceed, user->g_ceed, user->Q_dot_loc, NULL,
297f9028c3cSJames Wright                                  &op_ijacobian_ctx);
298f9028c3cSJames Wright       PetscCall(MatShellSetContext(J, op_ijacobian_ctx));
299f9028c3cSJames Wright       PetscCall(MatShellSetContextDestroy(J, (PetscErrorCode(*)(void *))OperatorApplyContextDestroy));
300f9028c3cSJames Wright       PetscCall(MatShellSetOperation(J, MATOP_MULT, (void (*)(void))MatMult_Ceed));
301f9028c3cSJames Wright       PetscCall(MatShellSetOperation(J, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiag_Ceed));
302f0b65372SJed Brown       PetscCall(MatSetUp(J));
303f0b65372SJed Brown     }
304f0b65372SJed Brown     if (!J_pre_is_shell) {
3052b916ea7SJeremy L Thompson       PetscCall(FormPreallocation(user, user->app_ctx->pmat_pbdiagonal, J_pre, &user->coo_values_pmat));
306b107fddaSJed Brown     }
30704855949SJed Brown     if (J != J_pre && !J_is_shell && !J_is_mffd) {
308b107fddaSJed Brown       PetscCall(FormPreallocation(user, PETSC_FALSE, J, &user->coo_values_amat));
309b107fddaSJed Brown     }
310f0b65372SJed Brown     user->matrices_set_up = true;
311f0b65372SJed Brown   }
312f0b65372SJed Brown   if (!J_pre_is_shell) {
3132b916ea7SJeremy L Thompson     PetscCall(FormSetValues(user, user->app_ctx->pmat_pbdiagonal, J_pre, user->coo_values_pmat));
314f0b65372SJed Brown   }
31504855949SJed Brown   if (user->coo_values_amat) {
31604855949SJed Brown     PetscCall(FormSetValues(user, PETSC_FALSE, J, user->coo_values_amat));
31704855949SJed Brown   } else if (J_is_mffd) {
31804855949SJed Brown     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
31904855949SJed Brown     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
32004855949SJed Brown   }
321f0b65372SJed Brown   PetscFunctionReturn(0);
322f0b65372SJed Brown }
323f0b65372SJed Brown 
3242b916ea7SJeremy L Thompson PetscErrorCode WriteOutput(User user, Vec Q, PetscInt step_no, PetscScalar time) {
325a515125bSLeila Ghaffari   Vec         Q_loc;
326a515125bSLeila Ghaffari   char        file_path[PETSC_MAX_PATH_LEN];
327a515125bSLeila Ghaffari   PetscViewer viewer;
328a515125bSLeila Ghaffari   PetscFunctionBeginUser;
329a515125bSLeila Ghaffari 
330852e5969SJed Brown   if (user->app_ctx->checkpoint_vtk) {
331a515125bSLeila Ghaffari     // Set up output
3327538d537SJames Wright     PetscCall(DMGetLocalVector(user->dm, &Q_loc));
3337538d537SJames Wright     PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec"));
3347538d537SJames Wright     PetscCall(VecZeroEntries(Q_loc));
3357538d537SJames Wright     PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
336a515125bSLeila Ghaffari 
337a515125bSLeila Ghaffari     // Output
338852e5969SJed Brown     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
3397538d537SJames Wright 
3402b916ea7SJeremy L Thompson     PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path, FILE_MODE_WRITE, &viewer));
3417538d537SJames Wright     PetscCall(VecView(Q_loc, viewer));
3427538d537SJames Wright     PetscCall(PetscViewerDestroy(&viewer));
343a515125bSLeila Ghaffari     if (user->dm_viz) {
344a515125bSLeila Ghaffari       Vec         Q_refined, Q_refined_loc;
345a515125bSLeila Ghaffari       char        file_path_refined[PETSC_MAX_PATH_LEN];
346a515125bSLeila Ghaffari       PetscViewer viewer_refined;
347a515125bSLeila Ghaffari 
3487538d537SJames Wright       PetscCall(DMGetGlobalVector(user->dm_viz, &Q_refined));
3497538d537SJames Wright       PetscCall(DMGetLocalVector(user->dm_viz, &Q_refined_loc));
3507538d537SJames Wright       PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined"));
3517538d537SJames Wright 
3527538d537SJames Wright       PetscCall(MatInterpolate(user->interp_viz, Q, Q_refined));
3537538d537SJames Wright       PetscCall(VecZeroEntries(Q_refined_loc));
3542b916ea7SJeremy L Thompson       PetscCall(DMGlobalToLocal(user->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc));
3557538d537SJames Wright 
356852e5969SJed Brown       PetscCall(
357852e5969SJed Brown           PetscSNPrintf(file_path_refined, sizeof file_path_refined, "%s/nsrefined-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
3587538d537SJames Wright 
3592b916ea7SJeremy L Thompson       PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined), file_path_refined, FILE_MODE_WRITE, &viewer_refined));
3607538d537SJames Wright       PetscCall(VecView(Q_refined_loc, viewer_refined));
3617538d537SJames Wright       PetscCall(DMRestoreLocalVector(user->dm_viz, &Q_refined_loc));
3627538d537SJames Wright       PetscCall(DMRestoreGlobalVector(user->dm_viz, &Q_refined));
3637538d537SJames Wright       PetscCall(PetscViewerDestroy(&viewer_refined));
364a515125bSLeila Ghaffari     }
3657538d537SJames Wright     PetscCall(DMRestoreLocalVector(user->dm, &Q_loc));
366852e5969SJed Brown   }
367a515125bSLeila Ghaffari 
368a515125bSLeila Ghaffari   // Save data in a binary file for continuation of simulations
36991a36801SJames Wright   if (user->app_ctx->add_stepnum2bin) {
370852e5969SJed Brown     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", user->app_ctx->output_dir, step_no));
37191a36801SJames Wright   } else {
3722b916ea7SJeremy L Thompson     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", user->app_ctx->output_dir));
37391a36801SJames Wright   }
3742b916ea7SJeremy L Thompson   PetscCall(PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, &viewer));
3757538d537SJames Wright 
3769293eaa1SJed Brown   PetscInt token = FLUIDS_FILE_TOKEN;
3779293eaa1SJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &token, 1, PETSC_INT));
3789293eaa1SJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &step_no, 1, PETSC_INT));
3799293eaa1SJed Brown   time /= user->units->second;  // Dimensionalize time back
3809293eaa1SJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL));
3817538d537SJames Wright   PetscCall(VecView(Q, viewer));
3827538d537SJames Wright   PetscCall(PetscViewerDestroy(&viewer));
3837538d537SJames Wright   PetscFunctionReturn(0);
3847538d537SJames Wright }
3857538d537SJames Wright 
386c5e9980aSAdeleke O. Bankole // CSV Monitor
387c5e9980aSAdeleke O. Bankole PetscErrorCode TSMonitor_WallForce(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
388c5e9980aSAdeleke O. Bankole   User              user = ctx;
389c5e9980aSAdeleke O. Bankole   Vec               G_loc;
390c5e9980aSAdeleke O. Bankole   PetscInt          num_wall = user->app_ctx->wall_forces.num_wall, dim = 3;
391c5e9980aSAdeleke O. Bankole   const PetscInt   *walls  = user->app_ctx->wall_forces.walls;
392c5e9980aSAdeleke O. Bankole   PetscViewer       viewer = user->app_ctx->wall_forces.viewer;
393c5e9980aSAdeleke O. Bankole   PetscViewerFormat format = user->app_ctx->wall_forces.viewer_format;
394c5e9980aSAdeleke O. Bankole   PetscScalar      *reaction_force;
395c5e9980aSAdeleke O. Bankole   PetscBool         iascii;
396c5e9980aSAdeleke O. Bankole 
397c5e9980aSAdeleke O. Bankole   PetscFunctionBeginUser;
398c5e9980aSAdeleke O. Bankole   if (!viewer) PetscFunctionReturn(0);
399c5e9980aSAdeleke O. Bankole   PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
400c5e9980aSAdeleke O. Bankole   PetscCall(PetscMalloc1(num_wall * dim, &reaction_force));
401c5e9980aSAdeleke O. Bankole   PetscCall(Surface_Forces_NS(user->dm, G_loc, num_wall, walls, reaction_force));
402c5e9980aSAdeleke O. Bankole   PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
403c5e9980aSAdeleke O. Bankole 
404c5e9980aSAdeleke O. Bankole   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
405c5e9980aSAdeleke O. Bankole 
406c5e9980aSAdeleke O. Bankole   if (iascii) {
407c5e9980aSAdeleke O. Bankole     if (format == PETSC_VIEWER_ASCII_CSV && !user->app_ctx->wall_forces.header_written) {
408c5e9980aSAdeleke O. Bankole       PetscCall(PetscViewerASCIIPrintf(viewer, "Step,Time,Wall,ForceX,ForceY,ForceZ\n"));
409c5e9980aSAdeleke O. Bankole       user->app_ctx->wall_forces.header_written = PETSC_TRUE;
410c5e9980aSAdeleke O. Bankole     }
411c5e9980aSAdeleke O. Bankole     for (PetscInt w = 0; w < num_wall; w++) {
412c5e9980aSAdeleke O. Bankole       PetscInt wall = walls[w];
413c5e9980aSAdeleke O. Bankole       if (format == PETSC_VIEWER_ASCII_CSV) {
414c5e9980aSAdeleke O. Bankole         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT ",%g,%" PetscInt_FMT ",%g,%g,%g\n", step_no, time, wall,
415c5e9980aSAdeleke O. Bankole                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
416c5e9980aSAdeleke O. Bankole 
417c5e9980aSAdeleke O. Bankole       } else {
418c5e9980aSAdeleke O. Bankole         PetscCall(PetscViewerASCIIPrintf(viewer, "Wall %" PetscInt_FMT " Forces: Force_x = %12g, Force_y = %12g, Force_z = %12g\n", wall,
419c5e9980aSAdeleke O. Bankole                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
420c5e9980aSAdeleke O. Bankole       }
421c5e9980aSAdeleke O. Bankole     }
422c5e9980aSAdeleke O. Bankole   }
423c5e9980aSAdeleke O. Bankole   PetscCall(PetscFree(reaction_force));
424c5e9980aSAdeleke O. Bankole   PetscFunctionReturn(0);
425c5e9980aSAdeleke O. Bankole }
426c5e9980aSAdeleke O. Bankole 
4277538d537SJames Wright // User provided TS Monitor
4282b916ea7SJeremy L Thompson PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
4297538d537SJames Wright   User user = ctx;
4307538d537SJames Wright   PetscFunctionBeginUser;
4317538d537SJames Wright 
432852e5969SJed Brown   // Print every 'checkpoint_interval' steps
433c539088bSJames Wright   if (user->app_ctx->checkpoint_interval <= 0 || step_no % user->app_ctx->checkpoint_interval != 0 ||
434e419654dSJeremy L Thompson       (user->app_ctx->cont_steps == step_no && step_no != 0)) {
435c539088bSJames Wright     PetscFunctionReturn(0);
436e419654dSJeremy L Thompson   }
4377538d537SJames Wright 
4387538d537SJames Wright   PetscCall(WriteOutput(user, Q, step_no, time));
439a515125bSLeila Ghaffari 
440a515125bSLeila Ghaffari   PetscFunctionReturn(0);
441a515125bSLeila Ghaffari }
442a515125bSLeila Ghaffari 
443a515125bSLeila Ghaffari // TS: Create, setup, and solve
4442b916ea7SJeremy L Thompson PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, Vec *Q, PetscScalar *f_time, TS *ts) {
445a515125bSLeila Ghaffari   MPI_Comm    comm = user->comm;
446a515125bSLeila Ghaffari   TSAdapt     adapt;
447a515125bSLeila Ghaffari   PetscScalar final_time;
448a515125bSLeila Ghaffari   PetscFunctionBeginUser;
449a515125bSLeila Ghaffari 
4502b916ea7SJeremy L Thompson   PetscCall(TSCreate(comm, ts));
4512b916ea7SJeremy L Thompson   PetscCall(TSSetDM(*ts, dm));
452a515125bSLeila Ghaffari   if (phys->implicit) {
4532b916ea7SJeremy L Thompson     PetscCall(TSSetType(*ts, TSBDF));
454a515125bSLeila Ghaffari     if (user->op_ifunction) {
4552b916ea7SJeremy L Thompson       PetscCall(TSSetIFunction(*ts, NULL, IFunction_NS, &user));
456a515125bSLeila Ghaffari     } else {  // Implicit integrators can fall back to using an RHSFunction
4572b916ea7SJeremy L Thompson       PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
458a515125bSLeila Ghaffari     }
459f0b65372SJed Brown     if (user->op_ijacobian) {
4602b916ea7SJeremy L Thompson       PetscCall(DMTSSetIJacobian(dm, FormIJacobian_NS, &user));
461b107fddaSJed Brown       if (app_ctx->amat_type) {
462b107fddaSJed Brown         Mat Pmat, Amat;
4632b916ea7SJeremy L Thompson         PetscCall(DMCreateMatrix(dm, &Pmat));
4642b916ea7SJeremy L Thompson         PetscCall(DMSetMatType(dm, app_ctx->amat_type));
4652b916ea7SJeremy L Thompson         PetscCall(DMCreateMatrix(dm, &Amat));
4662b916ea7SJeremy L Thompson         PetscCall(TSSetIJacobian(*ts, Amat, Pmat, NULL, NULL));
4672b916ea7SJeremy L Thompson         PetscCall(MatDestroy(&Amat));
4682b916ea7SJeremy L Thompson         PetscCall(MatDestroy(&Pmat));
469b107fddaSJed Brown       }
470f0b65372SJed Brown     }
471a515125bSLeila Ghaffari   } else {
4725d28dccaSJames Wright     PetscCheck(user->op_rhs, comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction");
4732b916ea7SJeremy L Thompson     PetscCall(TSSetType(*ts, TSRK));
4742b916ea7SJeremy L Thompson     PetscCall(TSRKSetType(*ts, TSRK5F));
4752b916ea7SJeremy L Thompson     PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
476a515125bSLeila Ghaffari   }
4772b916ea7SJeremy L Thompson   PetscCall(TSSetMaxTime(*ts, 500. * user->units->second));
4782b916ea7SJeremy L Thompson   PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER));
47922387d3aSJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(TSSetErrorIfStepFails(*ts, PETSC_FALSE));
4802b916ea7SJeremy L Thompson   PetscCall(TSSetTimeStep(*ts, 1.e-2 * user->units->second));
4810e1e9333SJames Wright   if (app_ctx->test_type != TESTTYPE_NONE) {
4822b916ea7SJeremy L Thompson     PetscCall(TSSetMaxSteps(*ts, 10));
4832b916ea7SJeremy L Thompson   }
4842b916ea7SJeremy L Thompson   PetscCall(TSGetAdapt(*ts, &adapt));
4852b916ea7SJeremy L Thompson   PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * user->units->second, 1.e2 * user->units->second));
4862b916ea7SJeremy L Thompson   PetscCall(TSSetFromOptions(*ts));
48791f639d2SJames Wright   user->time_bc_set = -1.0;    // require all BCs be updated
488a515125bSLeila Ghaffari   if (!app_ctx->cont_steps) {  // print initial condition
4890e1e9333SJames Wright     if (app_ctx->test_type == TESTTYPE_NONE) {
4902b916ea7SJeremy L Thompson       PetscCall(TSMonitor_NS(*ts, 0, 0., *Q, user));
491a515125bSLeila Ghaffari     }
492a515125bSLeila Ghaffari   } else {  // continue from time of last output
493a515125bSLeila Ghaffari     PetscInt    count;
494a515125bSLeila Ghaffari     PetscViewer viewer;
4952b916ea7SJeremy L Thompson 
4969293eaa1SJed Brown     if (app_ctx->cont_time <= 0) {  // Legacy files did not include step number and time
4972b916ea7SJeremy L Thompson       PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_time_file, FILE_MODE_READ, &viewer));
4989293eaa1SJed Brown       PetscCall(PetscViewerBinaryRead(viewer, &app_ctx->cont_time, 1, &count, PETSC_REAL));
4992b916ea7SJeremy L Thompson       PetscCall(PetscViewerDestroy(&viewer));
5009293eaa1SJed Brown       PetscCheck(app_ctx->cont_steps != -1, comm, PETSC_ERR_ARG_INCOMP,
5019293eaa1SJed Brown                  "-continue step number not specified, but checkpoint file does not contain a step number (likely written by older code version)");
5029293eaa1SJed Brown     }
5039293eaa1SJed Brown     PetscCall(TSSetTime(*ts, app_ctx->cont_time * user->units->second));
50474a6f4ddSJed Brown     PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps));
505a515125bSLeila Ghaffari   }
5060e1e9333SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) {
5072b916ea7SJeremy L Thompson     PetscCall(TSMonitorSet(*ts, TSMonitor_NS, user, NULL));
5080e1e9333SJames Wright   }
509c5e9980aSAdeleke O. Bankole   if (app_ctx->wall_forces.viewer) {
510c5e9980aSAdeleke O. Bankole     PetscCall(TSMonitorSet(*ts, TSMonitor_WallForce, user, NULL));
511c5e9980aSAdeleke O. Bankole   }
512c931fa59SJames Wright   if (app_ctx->turb_spanstats_enable) {
513b0488d1fSJames Wright     PetscCall(TSMonitorSet(*ts, TSMonitor_Statistics, user, NULL));
514b8daee98SJames Wright     CeedScalar previous_time = app_ctx->cont_time * user->units->second;
5155cfe26e6SJames Wright     CeedOperatorSetContextDouble(user->spanstats.op_stats_collect_ctx->op, user->spanstats.previous_time_label, &previous_time);
516b0488d1fSJames Wright   }
517a515125bSLeila Ghaffari 
518a515125bSLeila Ghaffari   // Solve
51974a6f4ddSJed Brown   PetscReal start_time;
52074a6f4ddSJed Brown   PetscInt  start_step;
5212b916ea7SJeremy L Thompson   PetscCall(TSGetTime(*ts, &start_time));
52274a6f4ddSJed Brown   PetscCall(TSGetStepNumber(*ts, &start_step));
52391982731SJeremy L Thompson 
52491982731SJeremy L Thompson   PetscPreLoadBegin(PETSC_FALSE, "Fluids Solve");
52591982731SJeremy L Thompson   PetscCall(TSSetTime(*ts, start_time));
52674a6f4ddSJed Brown   PetscCall(TSSetStepNumber(*ts, start_step));
52791982731SJeremy L Thompson   if (PetscPreLoadingOn) {
52891982731SJeremy L Thompson     // LCOV_EXCL_START
52991982731SJeremy L Thompson     SNES      snes;
53091982731SJeremy L Thompson     Vec       Q_preload;
53191982731SJeremy L Thompson     PetscReal rtol;
53291982731SJeremy L Thompson     PetscCall(VecDuplicate(*Q, &Q_preload));
53391982731SJeremy L Thompson     PetscCall(VecCopy(*Q, Q_preload));
53491982731SJeremy L Thompson     PetscCall(TSGetSNES(*ts, &snes));
53591982731SJeremy L Thompson     PetscCall(SNESGetTolerances(snes, NULL, &rtol, NULL, NULL, NULL));
5362b916ea7SJeremy L Thompson     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, .99, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
53791982731SJeremy L Thompson     PetscCall(TSSetSolution(*ts, *Q));
53891982731SJeremy L Thompson     PetscCall(TSStep(*ts));
5392b916ea7SJeremy L Thompson     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, rtol, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
54091982731SJeremy L Thompson     PetscCall(VecDestroy(&Q_preload));
54191982731SJeremy L Thompson     // LCOV_EXCL_STOP
54291982731SJeremy L Thompson   } else {
5432b916ea7SJeremy L Thompson     PetscCall(PetscBarrier((PetscObject)*ts));
5442b916ea7SJeremy L Thompson     PetscCall(TSSolve(*ts, *Q));
54591982731SJeremy L Thompson   }
54691982731SJeremy L Thompson   PetscPreLoadEnd();
54791982731SJeremy L Thompson 
54891982731SJeremy L Thompson   PetscCall(TSGetSolveTime(*ts, &final_time));
549a515125bSLeila Ghaffari   *f_time = final_time;
55091982731SJeremy L Thompson 
5510e1e9333SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) {
5527538d537SJames Wright     PetscInt step_no;
5537538d537SJames Wright     PetscCall(TSGetStepNumber(*ts, &step_no));
554b0488d1fSJames Wright     if (user->app_ctx->checkpoint_interval > 0 || user->app_ctx->checkpoint_interval == -1) {
5557538d537SJames Wright       PetscCall(WriteOutput(user, *Q, step_no, final_time));
5567538d537SJames Wright     }
5577538d537SJames Wright 
55891982731SJeremy L Thompson     PetscLogEvent stage_id;
55991982731SJeremy L Thompson     PetscStageLog stage_log;
56091982731SJeremy L Thompson 
56191982731SJeremy L Thompson     PetscCall(PetscLogStageGetId("Fluids Solve", &stage_id));
56291982731SJeremy L Thompson     PetscCall(PetscLogGetStageLog(&stage_log));
5632b916ea7SJeremy L Thompson     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time taken for solution (sec): %g\n", stage_log->stageInfo[stage_id].perfInfo.time));
564a515125bSLeila Ghaffari   }
565a515125bSLeila Ghaffari   PetscFunctionReturn(0);
566a515125bSLeila Ghaffari }
567