xref: /libCEED/examples/fluids/src/setupts.c (revision 457e73b2156816fd50927125e96139f6cbacb95c)
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   PetscFunctionBeginUser;
2677841947SLeila Ghaffari 
2777841947SLeila Ghaffari   // CEED Restriction
2877841947SLeila Ghaffari   CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q);
2977841947SLeila Ghaffari   CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &q_data_size);
3077841947SLeila Ghaffari 
3177841947SLeila Ghaffari   // CEED QFunction
32ef080ff9SJames Wright   PetscCall(CreateMassQFunction(ceed, num_comp_q, q_data_size, &qf_mass));
3377841947SLeila Ghaffari 
3477841947SLeila Ghaffari   // CEED Operator
3577841947SLeila Ghaffari   CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass);
36ef080ff9SJames Wright   CeedOperatorSetField(op_mass, "u", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
372b730f8bSJeremy L Thompson   CeedOperatorSetField(op_mass, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data);
382b730f8bSJeremy L Thompson   CeedOperatorSetField(op_mass, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
3977841947SLeila Ghaffari 
403221f4d3SJames Wright   PetscCall(OperatorApplyContextCreate(NULL, dm, ceed, op_mass, NULL, NULL, NULL, NULL, &op_mass_ctx));
4177841947SLeila Ghaffari 
423221f4d3SJames Wright   PetscCall(DMGetLocalVector(dm, &Ones_loc));
433221f4d3SJames Wright   PetscCall(VecSet(Ones_loc, 1));
443221f4d3SJames Wright   PetscCall(ApplyCeedOperatorLocalToGlobal(Ones_loc, M, op_mass_ctx));
4577841947SLeila Ghaffari 
4677841947SLeila Ghaffari   // Invert diagonally lumped mass vector for RHS function
472b730f8bSJeremy L Thompson   PetscCall(VecReciprocal(M));
4877841947SLeila Ghaffari 
4977841947SLeila Ghaffari   // Cleanup
503221f4d3SJames Wright   PetscCall(OperatorApplyContextDestroy(op_mass_ctx));
513221f4d3SJames Wright   PetscCall(DMRestoreLocalVector(dm, &Ones_loc));
5277841947SLeila Ghaffari   CeedQFunctionDestroy(&qf_mass);
5377841947SLeila Ghaffari   CeedOperatorDestroy(&op_mass);
5477841947SLeila Ghaffari 
55ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
5677841947SLeila Ghaffari }
5777841947SLeila Ghaffari 
58a0b9a424SJames Wright // Insert Boundary values if it's a new time
59a0b9a424SJames Wright PetscErrorCode UpdateBoundaryValues(User user, Vec Q_loc, PetscReal t) {
60a0b9a424SJames Wright   PetscFunctionBeginUser;
61a0b9a424SJames Wright   if (user->time_bc_set != t) {
62a0b9a424SJames Wright     PetscCall(DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL));
63a0b9a424SJames Wright     user->time_bc_set = t;
64a0b9a424SJames Wright   }
65ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
66a0b9a424SJames Wright }
67a0b9a424SJames Wright 
6815637395SJames Wright // @brief Update the context label value to new value if necessary.
6915637395SJames Wright // @note This only supports labels with scalar label values (ie. not arrays)
7015637395SJames Wright PetscErrorCode UpdateContextLabel(MPI_Comm comm, PetscScalar update_value, CeedOperator op, CeedContextFieldLabel label) {
7115637395SJames Wright   PetscScalar label_value;
72e42c1df6SJames Wright 
7315637395SJames Wright   PetscFunctionBeginUser;
7415637395SJames Wright   PetscCheck(label, comm, PETSC_ERR_ARG_BADPTR, "Label should be non-NULL");
7515637395SJames Wright 
7615637395SJames Wright   {
7715637395SJames Wright     size_t             num_elements;
7815637395SJames Wright     const PetscScalar *label_values;
7915637395SJames Wright     CeedOperatorGetContextDoubleRead(op, label, &num_elements, &label_values);
8015637395SJames Wright     PetscCheck(num_elements == 1, comm, PETSC_ERR_SUP, "%s does not support labels with more than 1 value. Label has %zu values", __func__,
8115637395SJames Wright                num_elements);
8215637395SJames Wright     label_value = *label_values;
8315637395SJames Wright     CeedOperatorRestoreContextDoubleRead(op, label, &label_values);
84e42c1df6SJames Wright   }
8515637395SJames Wright 
8615637395SJames Wright   if (label_value != update_value) {
8715637395SJames Wright     CeedOperatorSetContextDouble(op, label, &update_value);
88e42c1df6SJames Wright   }
89ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
90e42c1df6SJames Wright }
91e42c1df6SJames Wright 
9277841947SLeila Ghaffari // RHS (Explicit time-stepper) function setup
9377841947SLeila Ghaffari //   This is the RHS of the ODE, given as u_t = G(t,u)
9477841947SLeila Ghaffari //   This function takes in a state vector Q and writes into G
9577841947SLeila Ghaffari PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) {
9677841947SLeila Ghaffari   User        user = *(User *)user_data;
9715637395SJames Wright   MPI_Comm    comm = PetscObjectComm((PetscObject)ts);
98c798d55aSJames Wright   PetscScalar dt;
999ad5e8e4SJames Wright   Vec         Q_loc = user->Q_loc;
10077841947SLeila Ghaffari   PetscFunctionBeginUser;
10177841947SLeila Ghaffari 
1025e82a6e1SJeremy L Thompson   // Update time dependent data
103a0b9a424SJames Wright   PetscCall(UpdateBoundaryValues(user, Q_loc, t));
1049ad5e8e4SJames Wright   if (user->phys->solution_time_label) PetscCall(UpdateContextLabel(comm, t, user->op_rhs_ctx->op, user->phys->solution_time_label));
1052b730f8bSJeremy L Thompson   PetscCall(TSGetTimeStep(ts, &dt));
1069ad5e8e4SJames Wright   if (user->phys->timestep_size_label) PetscCall(UpdateContextLabel(comm, dt, user->op_rhs_ctx->op, user->phys->timestep_size_label));
10777841947SLeila Ghaffari 
1089ad5e8e4SJames Wright   PetscCall(ApplyCeedOperatorGlobalToGlobal(Q, G, user->op_rhs_ctx));
10977841947SLeila Ghaffari 
11077841947SLeila Ghaffari   // Inverse of the lumped mass matrix (M is Minv)
111186595e6SJames Wright   PetscCall(VecPointwiseMult(G, G, user->M_inv));
11277841947SLeila Ghaffari 
113ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
11477841947SLeila Ghaffari }
11577841947SLeila Ghaffari 
116ca69d878SAdeleke O. Bankole // Surface forces function setup
117ca69d878SAdeleke O. Bankole static PetscErrorCode Surface_Forces_NS(DM dm, Vec G_loc, PetscInt num_walls, const PetscInt walls[], PetscScalar *reaction_force) {
118ca69d878SAdeleke O. Bankole   DMLabel            face_label;
119ca69d878SAdeleke O. Bankole   const PetscScalar *g;
120d6734f85SAdeleke O. Bankole   PetscInt           dof, dim = 3;
121ca69d878SAdeleke O. Bankole   MPI_Comm           comm;
122d6734f85SAdeleke O. Bankole   PetscSection       s;
123ca69d878SAdeleke O. Bankole 
124ca69d878SAdeleke O. Bankole   PetscFunctionBeginUser;
125ca69d878SAdeleke O. Bankole   PetscCall(PetscArrayzero(reaction_force, num_walls * dim));
126ca69d878SAdeleke O. Bankole   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
127ca69d878SAdeleke O. Bankole   PetscCall(DMGetLabel(dm, "Face Sets", &face_label));
128ca69d878SAdeleke O. Bankole   PetscCall(VecGetArrayRead(G_loc, &g));
129ca69d878SAdeleke O. Bankole   for (PetscInt w = 0; w < num_walls; w++) {
130ca69d878SAdeleke O. Bankole     const PetscInt wall = walls[w];
131ca69d878SAdeleke O. Bankole     IS             wall_is;
132d6734f85SAdeleke O. Bankole     PetscCall(DMGetLocalSection(dm, &s));
133ca69d878SAdeleke O. Bankole     PetscCall(DMLabelGetStratumIS(face_label, wall, &wall_is));
134ca69d878SAdeleke O. Bankole     if (wall_is) {  // There exist such points on this process
135ca69d878SAdeleke O. Bankole       PetscInt        num_points;
136d6734f85SAdeleke O. Bankole       PetscInt        num_comp = 0;
137ca69d878SAdeleke O. Bankole       const PetscInt *points;
138d6734f85SAdeleke O. Bankole       PetscCall(PetscSectionGetFieldComponents(s, 0, &num_comp));
139ca69d878SAdeleke O. Bankole       PetscCall(ISGetSize(wall_is, &num_points));
140ca69d878SAdeleke O. Bankole       PetscCall(ISGetIndices(wall_is, &points));
141ca69d878SAdeleke O. Bankole       for (PetscInt i = 0; i < num_points; i++) {
142ca69d878SAdeleke O. Bankole         const PetscInt           p = points[i];
143ca69d878SAdeleke O. Bankole         const StateConservative *r;
144ca69d878SAdeleke O. Bankole         PetscCall(DMPlexPointLocalRead(dm, p, g, &r));
145d6734f85SAdeleke O. Bankole         PetscCall(PetscSectionGetDof(s, p, &dof));
146d6734f85SAdeleke O. Bankole         for (PetscInt node = 0; node < dof / num_comp; node++) {
147ca69d878SAdeleke O. Bankole           for (PetscInt j = 0; j < 3; j++) {
148d6734f85SAdeleke O. Bankole             reaction_force[w * dim + j] -= r[node].momentum[j];
149d6734f85SAdeleke O. Bankole           }
150ca69d878SAdeleke O. Bankole         }
151ca69d878SAdeleke O. Bankole       }
152ca69d878SAdeleke O. Bankole       PetscCall(ISRestoreIndices(wall_is, &points));
153ca69d878SAdeleke O. Bankole     }
154ca69d878SAdeleke O. Bankole     PetscCall(ISDestroy(&wall_is));
155ca69d878SAdeleke O. Bankole   }
156ca69d878SAdeleke O. Bankole   PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, reaction_force, dim * num_walls, MPIU_SCALAR, MPI_SUM, comm));
157ca69d878SAdeleke O. Bankole   //  Restore Vectors
158ca69d878SAdeleke O. Bankole   PetscCall(VecRestoreArrayRead(G_loc, &g));
159ca69d878SAdeleke O. Bankole 
160ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
161ca69d878SAdeleke O. Bankole }
162ca69d878SAdeleke O. Bankole 
16377841947SLeila Ghaffari // Implicit time-stepper function setup
1642b730f8bSJeremy L Thompson PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) {
16577841947SLeila Ghaffari   User         user = *(User *)user_data;
16615637395SJames Wright   MPI_Comm     comm = PetscObjectComm((PetscObject)ts);
167c798d55aSJames Wright   PetscScalar  dt;
1685e82a6e1SJeremy L Thompson   Vec          Q_loc = user->Q_loc, Q_dot_loc = user->Q_dot_loc, G_loc;
16977841947SLeila Ghaffari   PetscMemType q_mem_type, q_dot_mem_type, g_mem_type;
17077841947SLeila Ghaffari   PetscFunctionBeginUser;
17177841947SLeila Ghaffari 
1725e82a6e1SJeremy L Thompson   // Get local vectors
173ca69d878SAdeleke O. Bankole   PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
1745e82a6e1SJeremy L Thompson 
1755e82a6e1SJeremy L Thompson   // Update time dependent data
176a0b9a424SJames Wright   PetscCall(UpdateBoundaryValues(user, Q_loc, t));
17715637395SJames Wright   if (user->phys->solution_time_label) PetscCall(UpdateContextLabel(comm, t, user->op_ifunction, user->phys->solution_time_label));
1782b730f8bSJeremy L Thompson   PetscCall(TSGetTimeStep(ts, &dt));
17915637395SJames Wright   if (user->phys->timestep_size_label) PetscCall(UpdateContextLabel(comm, dt, user->op_ifunction, user->phys->timestep_size_label));
18077841947SLeila Ghaffari 
18177841947SLeila Ghaffari   // Global-to-local
182878eb0e7SJames Wright   PetscCall(DMGlobalToLocalBegin(user->dm, Q, INSERT_VALUES, Q_loc));
183878eb0e7SJames Wright   PetscCall(DMGlobalToLocalBegin(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
184878eb0e7SJames Wright   PetscCall(DMGlobalToLocalEnd(user->dm, Q, INSERT_VALUES, Q_loc));
185878eb0e7SJames Wright   PetscCall(DMGlobalToLocalEnd(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
18677841947SLeila Ghaffari 
18777841947SLeila Ghaffari   // Place PETSc vectors in CEED vectors
188c798d55aSJames Wright   PetscCall(VecReadP2C(Q_loc, &q_mem_type, user->q_ceed));
189c798d55aSJames Wright   PetscCall(VecReadP2C(Q_dot_loc, &q_dot_mem_type, user->q_dot_ceed));
190c798d55aSJames Wright   PetscCall(VecP2C(G_loc, &g_mem_type, user->g_ceed));
19177841947SLeila Ghaffari 
19277841947SLeila Ghaffari   // Apply CEED operator
19375d1798cSJames Wright   PetscCall(PetscLogEventBegin(FLUIDS_CeedOperatorApply, Q, G, 0, 0));
19475d1798cSJames Wright   PetscCall(PetscLogGpuTimeBegin());
1952b730f8bSJeremy L Thompson   CeedOperatorApply(user->op_ifunction, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE);
19675d1798cSJames Wright   PetscCall(PetscLogGpuTimeEnd());
19775d1798cSJames Wright   PetscCall(PetscLogEventEnd(FLUIDS_CeedOperatorApply, Q, G, 0, 0));
19877841947SLeila Ghaffari 
19977841947SLeila Ghaffari   // Restore vectors
200c798d55aSJames Wright   PetscCall(VecReadC2P(user->q_ceed, q_mem_type, Q_loc));
201c798d55aSJames Wright   PetscCall(VecReadC2P(user->q_dot_ceed, q_dot_mem_type, Q_dot_loc));
202c798d55aSJames Wright   PetscCall(VecC2P(user->g_ceed, g_mem_type, G_loc));
20377841947SLeila Ghaffari 
20419ffbc25SJames Wright   if (user->app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) {
20519ffbc25SJames Wright     PetscCall(SGS_DD_ModelApplyIFunction(user, Q_loc, G_loc));
20619ffbc25SJames Wright   }
207be34b3b7SJames Wright 
20877841947SLeila Ghaffari   // Local-to-Global
2092b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(G));
2102b730f8bSJeremy L Thompson   PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G));
21177841947SLeila Ghaffari 
21277841947SLeila Ghaffari   // Restore vectors
213ca69d878SAdeleke O. Bankole   PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
21477841947SLeila Ghaffari 
215ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
21677841947SLeila Ghaffari }
21777841947SLeila Ghaffari 
2182b730f8bSJeremy L Thompson static PetscErrorCode FormPreallocation(User user, PetscBool pbdiagonal, Mat J, CeedVector *coo_values) {
219544be873SJed Brown   PetscCount ncoo;
220*457e73b2SJames Wright   PetscInt  *rows_petsc, *cols_petsc;
221544be873SJed Brown 
222544be873SJed Brown   PetscFunctionBeginUser;
223544be873SJed Brown   if (pbdiagonal) {
224544be873SJed Brown     CeedSize l_size;
225544be873SJed Brown     CeedOperatorGetActiveVectorLengths(user->op_ijacobian, &l_size, NULL);
226544be873SJed Brown     ncoo       = l_size * 5;
227*457e73b2SJames Wright     rows_petsc = malloc(ncoo * sizeof(rows_petsc));
228*457e73b2SJames Wright     cols_petsc = malloc(ncoo * sizeof(cols_petsc));
229544be873SJed Brown     for (PetscCount n = 0; n < l_size / 5; n++) {
230544be873SJed Brown       for (PetscInt i = 0; i < 5; i++) {
231544be873SJed Brown         for (PetscInt j = 0; j < 5; j++) {
232*457e73b2SJames Wright           rows_petsc[(n * 5 + i) * 5 + j] = n * 5 + i;
233*457e73b2SJames Wright           cols_petsc[(n * 5 + i) * 5 + j] = n * 5 + j;
234544be873SJed Brown         }
235544be873SJed Brown       }
236544be873SJed Brown     }
237544be873SJed Brown   } else {
238*457e73b2SJames Wright     CeedInt *rows_ceed, *cols_ceed;
239*457e73b2SJames Wright     PetscCall(CeedOperatorLinearAssembleSymbolic(user->op_ijacobian, &ncoo, &rows_ceed, &cols_ceed));
240*457e73b2SJames Wright     PetscCall(IntArrayC2P(ncoo, &rows_ceed, &rows_petsc));
241*457e73b2SJames Wright     PetscCall(IntArrayC2P(ncoo, &cols_ceed, &cols_petsc));
242544be873SJed Brown   }
243*457e73b2SJames Wright   PetscCall(MatSetPreallocationCOOLocal(J, ncoo, rows_petsc, cols_petsc));
244*457e73b2SJames Wright   free(rows_petsc);
245*457e73b2SJames Wright   free(cols_petsc);
246544be873SJed Brown   CeedVectorCreate(user->ceed, ncoo, coo_values);
247ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
248544be873SJed Brown }
249544be873SJed Brown 
2502b730f8bSJeremy L Thompson static PetscErrorCode FormSetValues(User user, PetscBool pbdiagonal, Mat J, CeedVector coo_values) {
251544be873SJed Brown   CeedMemType        mem_type = CEED_MEM_HOST;
252544be873SJed Brown   const PetscScalar *values;
253544be873SJed Brown   MatType            mat_type;
254544be873SJed Brown 
255544be873SJed Brown   PetscFunctionBeginUser;
256544be873SJed Brown   PetscCall(MatGetType(J, &mat_type));
2572b730f8bSJeremy L Thompson   if (strstr(mat_type, "kokkos") || strstr(mat_type, "cusparse")) mem_type = CEED_MEM_DEVICE;
258fbb546ddSJames Wright   if (pbdiagonal) {
25975d1798cSJames Wright     PetscCall(PetscLogEventBegin(FLUIDS_CeedOperatorAssemblePointBlockDiagonal, J, 0, 0, 0));
26075d1798cSJames Wright     PetscCall(PetscLogGpuTimeBegin());
2612b730f8bSJeremy L Thompson     CeedOperatorLinearAssemblePointBlockDiagonal(user->op_ijacobian, coo_values, CEED_REQUEST_IMMEDIATE);
26275d1798cSJames Wright     PetscCall(PetscLogGpuTimeEnd());
26375d1798cSJames Wright     PetscCall(PetscLogEventEnd(FLUIDS_CeedOperatorAssemblePointBlockDiagonal, J, 0, 0, 0));
264544be873SJed Brown   } else {
26575d1798cSJames Wright     PetscCall(PetscLogEventBegin(FLUIDS_CeedOperatorAssemble, J, 0, 0, 0));
26675d1798cSJames Wright     PetscCall(PetscLogGpuTimeBegin());
267544be873SJed Brown     CeedOperatorLinearAssemble(user->op_ijacobian, coo_values);
26875d1798cSJames Wright     PetscCall(PetscLogGpuTimeEnd());
26975d1798cSJames Wright     PetscCall(PetscLogEventEnd(FLUIDS_CeedOperatorAssemble, J, 0, 0, 0));
270544be873SJed Brown   }
271544be873SJed Brown   CeedVectorGetArrayRead(coo_values, mem_type, &values);
272544be873SJed Brown   PetscCall(MatSetValuesCOO(J, values, INSERT_VALUES));
273544be873SJed Brown   CeedVectorRestoreArrayRead(coo_values, &values);
274ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
275544be873SJed Brown }
276544be873SJed Brown 
2772b730f8bSJeremy L Thompson PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, PetscReal shift, Mat J, Mat J_pre, void *user_data) {
278e334ad8fSJed Brown   User      user = *(User *)user_data;
279d6bf345cSJed Brown   PetscBool J_is_shell, J_is_mffd, J_pre_is_shell;
280e334ad8fSJed Brown   PetscFunctionBeginUser;
28117b0d5c6SJeremy L Thompson   if (user->phys->ijacobian_time_shift_label) CeedOperatorSetContextDouble(user->op_ijacobian, user->phys->ijacobian_time_shift_label, &shift);
282d6bf345cSJed Brown   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATMFFD, &J_is_mffd));
283e334ad8fSJed Brown   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATSHELL, &J_is_shell));
2842b730f8bSJeremy L Thompson   PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATSHELL, &J_pre_is_shell));
285e334ad8fSJed Brown   if (!user->matrices_set_up) {
286e334ad8fSJed Brown     if (J_is_shell) {
287b07979f9SJames Wright       OperatorApplyContext op_ijacobian_ctx;
288b07979f9SJames Wright       OperatorApplyContextCreate(user->dm, user->dm, user->ceed, user->op_ijacobian, user->q_ceed, user->g_ceed, user->Q_dot_loc, NULL,
289b07979f9SJames Wright                                  &op_ijacobian_ctx);
290b07979f9SJames Wright       PetscCall(MatShellSetContext(J, op_ijacobian_ctx));
291b07979f9SJames Wright       PetscCall(MatShellSetContextDestroy(J, (PetscErrorCode(*)(void *))OperatorApplyContextDestroy));
292b07979f9SJames Wright       PetscCall(MatShellSetOperation(J, MATOP_MULT, (void (*)(void))MatMult_Ceed));
293b07979f9SJames Wright       PetscCall(MatShellSetOperation(J, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiag_Ceed));
294e334ad8fSJed Brown       PetscCall(MatSetUp(J));
295e334ad8fSJed Brown     }
296e334ad8fSJed Brown     if (!J_pre_is_shell) {
2972b730f8bSJeremy L Thompson       PetscCall(FormPreallocation(user, user->app_ctx->pmat_pbdiagonal, J_pre, &user->coo_values_pmat));
298544be873SJed Brown     }
299d6bf345cSJed Brown     if (J != J_pre && !J_is_shell && !J_is_mffd) {
300544be873SJed Brown       PetscCall(FormPreallocation(user, PETSC_FALSE, J, &user->coo_values_amat));
301544be873SJed Brown     }
302e334ad8fSJed Brown     user->matrices_set_up = true;
303e334ad8fSJed Brown   }
304e334ad8fSJed Brown   if (!J_pre_is_shell) {
3052b730f8bSJeremy L Thompson     PetscCall(FormSetValues(user, user->app_ctx->pmat_pbdiagonal, J_pre, user->coo_values_pmat));
306e334ad8fSJed Brown   }
307d6bf345cSJed Brown   if (user->coo_values_amat) {
308d6bf345cSJed Brown     PetscCall(FormSetValues(user, PETSC_FALSE, J, user->coo_values_amat));
309d6bf345cSJed Brown   } else if (J_is_mffd) {
310d6bf345cSJed Brown     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
311d6bf345cSJed Brown     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
312d6bf345cSJed Brown   }
313ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
314e334ad8fSJed Brown }
315e334ad8fSJed Brown 
3162b730f8bSJeremy L Thompson PetscErrorCode WriteOutput(User user, Vec Q, PetscInt step_no, PetscScalar time) {
31777841947SLeila Ghaffari   Vec         Q_loc;
31877841947SLeila Ghaffari   char        file_path[PETSC_MAX_PATH_LEN];
31977841947SLeila Ghaffari   PetscViewer viewer;
32077841947SLeila Ghaffari   PetscFunctionBeginUser;
32177841947SLeila Ghaffari 
32237cbb16aSJed Brown   if (user->app_ctx->checkpoint_vtk) {
32377841947SLeila Ghaffari     // Set up output
324f14660a4SJames Wright     PetscCall(DMGetLocalVector(user->dm, &Q_loc));
325f14660a4SJames Wright     PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec"));
326f14660a4SJames Wright     PetscCall(VecZeroEntries(Q_loc));
327f14660a4SJames Wright     PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
32877841947SLeila Ghaffari 
32977841947SLeila Ghaffari     // Output
33037cbb16aSJed Brown     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
331f14660a4SJames Wright 
3322b730f8bSJeremy L Thompson     PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path, FILE_MODE_WRITE, &viewer));
333f14660a4SJames Wright     PetscCall(VecView(Q_loc, viewer));
334f14660a4SJames Wright     PetscCall(PetscViewerDestroy(&viewer));
33577841947SLeila Ghaffari     if (user->dm_viz) {
33677841947SLeila Ghaffari       Vec         Q_refined, Q_refined_loc;
33777841947SLeila Ghaffari       char        file_path_refined[PETSC_MAX_PATH_LEN];
33877841947SLeila Ghaffari       PetscViewer viewer_refined;
33977841947SLeila Ghaffari 
340f14660a4SJames Wright       PetscCall(DMGetGlobalVector(user->dm_viz, &Q_refined));
341f14660a4SJames Wright       PetscCall(DMGetLocalVector(user->dm_viz, &Q_refined_loc));
342f14660a4SJames Wright       PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined"));
343f14660a4SJames Wright 
344f14660a4SJames Wright       PetscCall(MatInterpolate(user->interp_viz, Q, Q_refined));
345f14660a4SJames Wright       PetscCall(VecZeroEntries(Q_refined_loc));
3462b730f8bSJeremy L Thompson       PetscCall(DMGlobalToLocal(user->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc));
347f14660a4SJames Wright 
34837cbb16aSJed Brown       PetscCall(
34937cbb16aSJed Brown           PetscSNPrintf(file_path_refined, sizeof file_path_refined, "%s/nsrefined-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
350f14660a4SJames Wright 
3512b730f8bSJeremy L Thompson       PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined), file_path_refined, FILE_MODE_WRITE, &viewer_refined));
352f14660a4SJames Wright       PetscCall(VecView(Q_refined_loc, viewer_refined));
353f14660a4SJames Wright       PetscCall(DMRestoreLocalVector(user->dm_viz, &Q_refined_loc));
354f14660a4SJames Wright       PetscCall(DMRestoreGlobalVector(user->dm_viz, &Q_refined));
355f14660a4SJames Wright       PetscCall(PetscViewerDestroy(&viewer_refined));
35677841947SLeila Ghaffari     }
357f14660a4SJames Wright     PetscCall(DMRestoreLocalVector(user->dm, &Q_loc));
35837cbb16aSJed Brown   }
35977841947SLeila Ghaffari 
36077841947SLeila Ghaffari   // Save data in a binary file for continuation of simulations
36169293791SJames Wright   if (user->app_ctx->add_stepnum2bin) {
36237cbb16aSJed Brown     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", user->app_ctx->output_dir, step_no));
36369293791SJames Wright   } else {
3642b730f8bSJeremy L Thompson     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", user->app_ctx->output_dir));
36569293791SJames Wright   }
3662b730f8bSJeremy L Thompson   PetscCall(PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, &viewer));
367f14660a4SJames Wright 
3684de8550aSJed Brown   PetscInt token = FLUIDS_FILE_TOKEN;
3694de8550aSJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &token, 1, PETSC_INT));
3704de8550aSJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &step_no, 1, PETSC_INT));
3714de8550aSJed Brown   time /= user->units->second;  // Dimensionalize time back
3724de8550aSJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL));
373f14660a4SJames Wright   PetscCall(VecView(Q, viewer));
374f14660a4SJames Wright   PetscCall(PetscViewerDestroy(&viewer));
375ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
376f14660a4SJames Wright }
377f14660a4SJames Wright 
378ca69d878SAdeleke O. Bankole // CSV Monitor
379ca69d878SAdeleke O. Bankole PetscErrorCode TSMonitor_WallForce(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
380ca69d878SAdeleke O. Bankole   User              user = ctx;
381ca69d878SAdeleke O. Bankole   Vec               G_loc;
382ca69d878SAdeleke O. Bankole   PetscInt          num_wall = user->app_ctx->wall_forces.num_wall, dim = 3;
383ca69d878SAdeleke O. Bankole   const PetscInt   *walls  = user->app_ctx->wall_forces.walls;
384ca69d878SAdeleke O. Bankole   PetscViewer       viewer = user->app_ctx->wall_forces.viewer;
385ca69d878SAdeleke O. Bankole   PetscViewerFormat format = user->app_ctx->wall_forces.viewer_format;
386ca69d878SAdeleke O. Bankole   PetscScalar      *reaction_force;
387ca69d878SAdeleke O. Bankole   PetscBool         iascii;
388ca69d878SAdeleke O. Bankole 
389ca69d878SAdeleke O. Bankole   PetscFunctionBeginUser;
390ee4ca9cbSJames Wright   if (!viewer) PetscFunctionReturn(PETSC_SUCCESS);
391ca69d878SAdeleke O. Bankole   PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
392ca69d878SAdeleke O. Bankole   PetscCall(PetscMalloc1(num_wall * dim, &reaction_force));
393ca69d878SAdeleke O. Bankole   PetscCall(Surface_Forces_NS(user->dm, G_loc, num_wall, walls, reaction_force));
394ca69d878SAdeleke O. Bankole   PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
395ca69d878SAdeleke O. Bankole 
396ca69d878SAdeleke O. Bankole   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
397ca69d878SAdeleke O. Bankole 
398ca69d878SAdeleke O. Bankole   if (iascii) {
399ca69d878SAdeleke O. Bankole     if (format == PETSC_VIEWER_ASCII_CSV && !user->app_ctx->wall_forces.header_written) {
400ca69d878SAdeleke O. Bankole       PetscCall(PetscViewerASCIIPrintf(viewer, "Step,Time,Wall,ForceX,ForceY,ForceZ\n"));
401ca69d878SAdeleke O. Bankole       user->app_ctx->wall_forces.header_written = PETSC_TRUE;
402ca69d878SAdeleke O. Bankole     }
403ca69d878SAdeleke O. Bankole     for (PetscInt w = 0; w < num_wall; w++) {
404ca69d878SAdeleke O. Bankole       PetscInt wall = walls[w];
405ca69d878SAdeleke O. Bankole       if (format == PETSC_VIEWER_ASCII_CSV) {
406ca69d878SAdeleke O. Bankole         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT ",%g,%" PetscInt_FMT ",%g,%g,%g\n", step_no, time, wall,
407ca69d878SAdeleke O. Bankole                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
408ca69d878SAdeleke O. Bankole 
409ca69d878SAdeleke O. Bankole       } else {
410ca69d878SAdeleke O. Bankole         PetscCall(PetscViewerASCIIPrintf(viewer, "Wall %" PetscInt_FMT " Forces: Force_x = %12g, Force_y = %12g, Force_z = %12g\n", wall,
411ca69d878SAdeleke O. Bankole                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
412ca69d878SAdeleke O. Bankole       }
413ca69d878SAdeleke O. Bankole     }
414ca69d878SAdeleke O. Bankole   }
415ca69d878SAdeleke O. Bankole   PetscCall(PetscFree(reaction_force));
416ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
417ca69d878SAdeleke O. Bankole }
418ca69d878SAdeleke O. Bankole 
419f14660a4SJames Wright // User provided TS Monitor
4202b730f8bSJeremy L Thompson PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
421f14660a4SJames Wright   User user = ctx;
422f14660a4SJames Wright   PetscFunctionBeginUser;
423f14660a4SJames Wright 
42437cbb16aSJed Brown   // Print every 'checkpoint_interval' steps
425894de27cSJames Wright   if (user->app_ctx->checkpoint_interval <= 0 || step_no % user->app_ctx->checkpoint_interval != 0 ||
42649aac155SJeremy L Thompson       (user->app_ctx->cont_steps == step_no && step_no != 0)) {
427ee4ca9cbSJames Wright     PetscFunctionReturn(PETSC_SUCCESS);
42849aac155SJeremy L Thompson   }
429f14660a4SJames Wright 
430f14660a4SJames Wright   PetscCall(WriteOutput(user, Q, step_no, time));
43177841947SLeila Ghaffari 
432ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
43377841947SLeila Ghaffari }
43477841947SLeila Ghaffari 
43577841947SLeila Ghaffari // TS: Create, setup, and solve
4362b730f8bSJeremy L Thompson PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, Vec *Q, PetscScalar *f_time, TS *ts) {
43777841947SLeila Ghaffari   MPI_Comm    comm = user->comm;
43877841947SLeila Ghaffari   TSAdapt     adapt;
43977841947SLeila Ghaffari   PetscScalar final_time;
44077841947SLeila Ghaffari   PetscFunctionBeginUser;
44177841947SLeila Ghaffari 
4422b730f8bSJeremy L Thompson   PetscCall(TSCreate(comm, ts));
4432b730f8bSJeremy L Thompson   PetscCall(TSSetDM(*ts, dm));
44477841947SLeila Ghaffari   if (phys->implicit) {
4452b730f8bSJeremy L Thompson     PetscCall(TSSetType(*ts, TSBDF));
44677841947SLeila Ghaffari     if (user->op_ifunction) {
4472b730f8bSJeremy L Thompson       PetscCall(TSSetIFunction(*ts, NULL, IFunction_NS, &user));
44877841947SLeila Ghaffari     } else {  // Implicit integrators can fall back to using an RHSFunction
4492b730f8bSJeremy L Thompson       PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
45077841947SLeila Ghaffari     }
451e334ad8fSJed Brown     if (user->op_ijacobian) {
4522b730f8bSJeremy L Thompson       PetscCall(DMTSSetIJacobian(dm, FormIJacobian_NS, &user));
453544be873SJed Brown       if (app_ctx->amat_type) {
454544be873SJed Brown         Mat Pmat, Amat;
4552b730f8bSJeremy L Thompson         PetscCall(DMCreateMatrix(dm, &Pmat));
4562b730f8bSJeremy L Thompson         PetscCall(DMSetMatType(dm, app_ctx->amat_type));
4572b730f8bSJeremy L Thompson         PetscCall(DMCreateMatrix(dm, &Amat));
4582b730f8bSJeremy L Thompson         PetscCall(TSSetIJacobian(*ts, Amat, Pmat, NULL, NULL));
4592b730f8bSJeremy L Thompson         PetscCall(MatDestroy(&Amat));
4602b730f8bSJeremy L Thompson         PetscCall(MatDestroy(&Pmat));
461544be873SJed Brown       }
462e334ad8fSJed Brown     }
46377841947SLeila Ghaffari   } else {
4649ad5e8e4SJames Wright     PetscCheck(user->op_rhs_ctx, comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction");
4652b730f8bSJeremy L Thompson     PetscCall(TSSetType(*ts, TSRK));
4662b730f8bSJeremy L Thompson     PetscCall(TSRKSetType(*ts, TSRK5F));
4672b730f8bSJeremy L Thompson     PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
46877841947SLeila Ghaffari   }
4692b730f8bSJeremy L Thompson   PetscCall(TSSetMaxTime(*ts, 500. * user->units->second));
4702b730f8bSJeremy L Thompson   PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER));
47151b00d91SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(TSSetErrorIfStepFails(*ts, PETSC_FALSE));
4722b730f8bSJeremy L Thompson   PetscCall(TSSetTimeStep(*ts, 1.e-2 * user->units->second));
4738fb33541SJames Wright   if (app_ctx->test_type != TESTTYPE_NONE) {
4742b730f8bSJeremy L Thompson     PetscCall(TSSetMaxSteps(*ts, 10));
4752b730f8bSJeremy L Thompson   }
4762b730f8bSJeremy L Thompson   PetscCall(TSGetAdapt(*ts, &adapt));
4772b730f8bSJeremy L Thompson   PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * user->units->second, 1.e2 * user->units->second));
4782b730f8bSJeremy L Thompson   PetscCall(TSSetFromOptions(*ts));
47915637395SJames Wright   user->time_bc_set = -1.0;   // require all BCs be updated
480d55646a4SJames Wright   if (app_ctx->cont_steps) {  // continue from previous timestep data
48177841947SLeila Ghaffari     PetscInt    count;
48277841947SLeila Ghaffari     PetscViewer viewer;
4832b730f8bSJeremy L Thompson 
4844de8550aSJed Brown     if (app_ctx->cont_time <= 0) {  // Legacy files did not include step number and time
4852b730f8bSJeremy L Thompson       PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_time_file, FILE_MODE_READ, &viewer));
4864de8550aSJed Brown       PetscCall(PetscViewerBinaryRead(viewer, &app_ctx->cont_time, 1, &count, PETSC_REAL));
4872b730f8bSJeremy L Thompson       PetscCall(PetscViewerDestroy(&viewer));
4884de8550aSJed Brown       PetscCheck(app_ctx->cont_steps != -1, comm, PETSC_ERR_ARG_INCOMP,
4894de8550aSJed Brown                  "-continue step number not specified, but checkpoint file does not contain a step number (likely written by older code version)");
4904de8550aSJed Brown     }
4914de8550aSJed Brown     PetscCall(TSSetTime(*ts, app_ctx->cont_time * user->units->second));
492d8e0aecdSJed Brown     PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps));
49377841947SLeila Ghaffari   }
4948fb33541SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) {
4952b730f8bSJeremy L Thompson     PetscCall(TSMonitorSet(*ts, TSMonitor_NS, user, NULL));
4968fb33541SJames Wright   }
497ca69d878SAdeleke O. Bankole   if (app_ctx->wall_forces.viewer) {
498ca69d878SAdeleke O. Bankole     PetscCall(TSMonitorSet(*ts, TSMonitor_WallForce, user, NULL));
499ca69d878SAdeleke O. Bankole   }
500b7d66439SJames Wright   if (app_ctx->turb_spanstats_enable) {
501f5452247SJames Wright     PetscCall(TSMonitorSet(*ts, TSMonitor_TurbulenceStatistics, user, NULL));
502495b9769SJames Wright     CeedScalar previous_time = app_ctx->cont_time * user->units->second;
503ed331efdSJames Wright     CeedOperatorSetContextDouble(user->spanstats.op_stats_collect_ctx->op, user->spanstats.previous_time_label, &previous_time);
504a175e481SJames Wright   }
5054e9802d1SJames Wright   if (app_ctx->diff_filter_monitor) PetscCall(TSMonitorSet(*ts, TSMonitor_DifferentialFilter, user, NULL));
50677841947SLeila Ghaffari 
50777841947SLeila Ghaffari   // Solve
508d8e0aecdSJed Brown   PetscReal start_time;
509d8e0aecdSJed Brown   PetscInt  start_step;
5102b730f8bSJeremy L Thompson   PetscCall(TSGetTime(*ts, &start_time));
511d8e0aecdSJed Brown   PetscCall(TSGetStepNumber(*ts, &start_step));
5127b39487dSJeremy L Thompson 
5137b39487dSJeremy L Thompson   PetscPreLoadBegin(PETSC_FALSE, "Fluids Solve");
5147b39487dSJeremy L Thompson   PetscCall(TSSetTime(*ts, start_time));
515d8e0aecdSJed Brown   PetscCall(TSSetStepNumber(*ts, start_step));
5167b39487dSJeremy L Thompson   if (PetscPreLoadingOn) {
5177b39487dSJeremy L Thompson     // LCOV_EXCL_START
5187b39487dSJeremy L Thompson     SNES      snes;
5197b39487dSJeremy L Thompson     Vec       Q_preload;
5207b39487dSJeremy L Thompson     PetscReal rtol;
5217b39487dSJeremy L Thompson     PetscCall(VecDuplicate(*Q, &Q_preload));
5227b39487dSJeremy L Thompson     PetscCall(VecCopy(*Q, Q_preload));
5237b39487dSJeremy L Thompson     PetscCall(TSGetSNES(*ts, &snes));
5247b39487dSJeremy L Thompson     PetscCall(SNESGetTolerances(snes, NULL, &rtol, NULL, NULL, NULL));
5252b730f8bSJeremy L Thompson     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, .99, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
5267b39487dSJeremy L Thompson     PetscCall(TSSetSolution(*ts, *Q));
5277b39487dSJeremy L Thompson     PetscCall(TSStep(*ts));
5282b730f8bSJeremy L Thompson     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, rtol, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
5297b39487dSJeremy L Thompson     PetscCall(VecDestroy(&Q_preload));
5307b39487dSJeremy L Thompson     // LCOV_EXCL_STOP
5317b39487dSJeremy L Thompson   } else {
5322b730f8bSJeremy L Thompson     PetscCall(PetscBarrier((PetscObject)*ts));
5332b730f8bSJeremy L Thompson     PetscCall(TSSolve(*ts, *Q));
5347b39487dSJeremy L Thompson   }
5357b39487dSJeremy L Thompson   PetscPreLoadEnd();
5367b39487dSJeremy L Thompson 
5377b39487dSJeremy L Thompson   PetscCall(TSGetSolveTime(*ts, &final_time));
53877841947SLeila Ghaffari   *f_time = final_time;
5397b39487dSJeremy L Thompson 
5408fb33541SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) {
541f14660a4SJames Wright     PetscInt step_no;
542f14660a4SJames Wright     PetscCall(TSGetStepNumber(*ts, &step_no));
543a175e481SJames Wright     if (user->app_ctx->checkpoint_interval > 0 || user->app_ctx->checkpoint_interval == -1) {
544f14660a4SJames Wright       PetscCall(WriteOutput(user, *Q, step_no, final_time));
545f14660a4SJames Wright     }
546f14660a4SJames Wright 
54775d1798cSJames Wright     PetscLogStage stage_id;
5487b39487dSJeremy L Thompson     PetscStageLog stage_log;
5497b39487dSJeremy L Thompson 
5507b39487dSJeremy L Thompson     PetscCall(PetscLogStageGetId("Fluids Solve", &stage_id));
5517b39487dSJeremy L Thompson     PetscCall(PetscLogGetStageLog(&stage_log));
5522b730f8bSJeremy L Thompson     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time taken for solution (sec): %g\n", stage_log->stageInfo[stage_id].perfInfo.time));
55377841947SLeila Ghaffari   }
554ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
55577841947SLeila Ghaffari }
556