xref: /honee/src/setupts.c (revision d6ca8aaa1fd88ecd7e943e5b4f181b6026914de2)
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 
18f8e2d240SJames Wright // @brief Create KSP to solve the inverse mass operator for explicit time stepping schemes
19*d6ca8aaaSJames Wright PetscErrorCode CreateKSPMassOperator(User user, CeedData ceed_data) {
20f8e2d240SJames Wright   Ceed                 ceed = user->ceed;
21f8e2d240SJames Wright   DM                   dm   = user->dm;
22a515125bSLeila Ghaffari   CeedQFunction        qf_mass;
23a515125bSLeila Ghaffari   CeedOperator         op_mass;
24f8e2d240SJames Wright   OperatorApplyContext mass_matop_ctx;
25a515125bSLeila Ghaffari   CeedInt              num_comp_q, q_data_size;
26a515125bSLeila Ghaffari 
2706f41313SJames Wright   PetscFunctionBeginUser;
28b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q));
29b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &q_data_size));
30a515125bSLeila Ghaffari 
319f59f36eSJames Wright   PetscCall(CreateMassQFunction(ceed, num_comp_q, q_data_size, &qf_mass));
32b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass));
33b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
3458e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
35b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
36a515125bSLeila Ghaffari 
37f8e2d240SJames Wright   {  // -- Setup KSP for mass operator
38f8e2d240SJames Wright     Mat      mat_mass;
39614823f4SJames Wright     Vec      Zeros_loc;
40f8e2d240SJames Wright     MPI_Comm comm = PetscObjectComm((PetscObject)dm);
41a515125bSLeila Ghaffari 
42614823f4SJames Wright     PetscCall(DMCreateLocalVector(dm, &Zeros_loc));
43614823f4SJames Wright     PetscCall(VecZeroEntries(Zeros_loc));
44614823f4SJames Wright     PetscCall(OperatorApplyContextCreate(dm, dm, ceed, op_mass, NULL, NULL, Zeros_loc, NULL, &mass_matop_ctx));
45f8e2d240SJames Wright     PetscCall(CreateMatShell_Ceed(mass_matop_ctx, &mat_mass));
46a515125bSLeila Ghaffari 
47f8e2d240SJames Wright     PetscCall(KSPCreate(comm, &user->mass_ksp));
48f8e2d240SJames Wright     PetscCall(KSPSetOptionsPrefix(user->mass_ksp, "mass_"));
49f8e2d240SJames Wright     {  // lumped by default
50f8e2d240SJames Wright       PC pc;
51f8e2d240SJames Wright       PetscCall(KSPGetPC(user->mass_ksp, &pc));
52f8e2d240SJames Wright       PetscCall(PCSetType(pc, PCJACOBI));
53f8e2d240SJames Wright       PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
54f8e2d240SJames Wright       PetscCall(KSPSetType(user->mass_ksp, KSPPREONLY));
55f8e2d240SJames Wright     }
56f8e2d240SJames Wright     PetscCall(KSPSetOperators(user->mass_ksp, mat_mass, mat_mass));
57f8e2d240SJames Wright     PetscCall(KSPSetFromOptions(user->mass_ksp));
58614823f4SJames Wright     PetscCall(VecDestroy(&Zeros_loc));
59614823f4SJames Wright     PetscCall(MatDestroy(&mat_mass));
60f8e2d240SJames Wright   }
61a515125bSLeila Ghaffari 
62a515125bSLeila Ghaffari   // Cleanup
63b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
64b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
65d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
66a515125bSLeila Ghaffari }
67a515125bSLeila Ghaffari 
68c996854bSJames Wright // Insert Boundary values if it's a new time
69c996854bSJames Wright PetscErrorCode UpdateBoundaryValues(User user, Vec Q_loc, PetscReal t) {
70c996854bSJames Wright   PetscFunctionBeginUser;
71c996854bSJames Wright   if (user->time_bc_set != t) {
72c996854bSJames Wright     PetscCall(DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL));
73c996854bSJames Wright     user->time_bc_set = t;
74c996854bSJames Wright   }
75d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
76c996854bSJames Wright }
77c996854bSJames Wright 
78a515125bSLeila Ghaffari // RHS (Explicit time-stepper) function setup
79a515125bSLeila Ghaffari //   This is the RHS of the ODE, given as u_t = G(t,u)
80a515125bSLeila Ghaffari //   This function takes in a state vector Q and writes into G
81a515125bSLeila Ghaffari PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) {
82a515125bSLeila Ghaffari   User        user = *(User *)user_data;
83701e5830SJames Wright   Ceed        ceed = user->ceed;
84fd969b44SJames Wright   PetscScalar dt;
85da5fe0e4SJames Wright   Vec         Q_loc = user->Q_loc;
86a515125bSLeila Ghaffari 
8706f41313SJames Wright   PetscFunctionBeginUser;
88e2f84137SJeremy L Thompson   // Update time dependent data
89c996854bSJames Wright   PetscCall(UpdateBoundaryValues(user, Q_loc, t));
90701e5830SJames Wright   if (user->phys->solution_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_rhs_ctx->op, user->phys->solution_time_label, &t));
912b916ea7SJeremy L Thompson   PetscCall(TSGetTimeStep(ts, &dt));
92701e5830SJames Wright   if (user->phys->timestep_size_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_rhs_ctx->op, user->phys->timestep_size_label, &dt));
93a515125bSLeila Ghaffari 
94da5fe0e4SJames Wright   PetscCall(ApplyCeedOperatorGlobalToGlobal(Q, G, user->op_rhs_ctx));
95a515125bSLeila Ghaffari 
96f3fcf8f4SJames Wright   // Inverse of the lumped mass matrix
97f8e2d240SJames Wright   PetscCall(KSPSolve(user->mass_ksp, G, G));
98d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
99a515125bSLeila Ghaffari }
100a515125bSLeila Ghaffari 
101c5e9980aSAdeleke O. Bankole // Surface forces function setup
102c5e9980aSAdeleke O. Bankole static PetscErrorCode Surface_Forces_NS(DM dm, Vec G_loc, PetscInt num_walls, const PetscInt walls[], PetscScalar *reaction_force) {
103c5e9980aSAdeleke O. Bankole   DMLabel            face_label;
104c5e9980aSAdeleke O. Bankole   const PetscScalar *g;
1052004e3acSAdeleke O. Bankole   PetscInt           dof, dim = 3;
106c5e9980aSAdeleke O. Bankole   MPI_Comm           comm;
1072004e3acSAdeleke O. Bankole   PetscSection       s;
108c5e9980aSAdeleke O. Bankole 
109c5e9980aSAdeleke O. Bankole   PetscFunctionBeginUser;
110c5e9980aSAdeleke O. Bankole   PetscCall(PetscArrayzero(reaction_force, num_walls * dim));
111c5e9980aSAdeleke O. Bankole   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
112c5e9980aSAdeleke O. Bankole   PetscCall(DMGetLabel(dm, "Face Sets", &face_label));
113c5e9980aSAdeleke O. Bankole   PetscCall(VecGetArrayRead(G_loc, &g));
114c5e9980aSAdeleke O. Bankole   for (PetscInt w = 0; w < num_walls; w++) {
115c5e9980aSAdeleke O. Bankole     const PetscInt wall = walls[w];
116c5e9980aSAdeleke O. Bankole     IS             wall_is;
1172004e3acSAdeleke O. Bankole     PetscCall(DMGetLocalSection(dm, &s));
118c5e9980aSAdeleke O. Bankole     PetscCall(DMLabelGetStratumIS(face_label, wall, &wall_is));
119c5e9980aSAdeleke O. Bankole     if (wall_is) {  // There exist such points on this process
120c5e9980aSAdeleke O. Bankole       PetscInt        num_points;
1212004e3acSAdeleke O. Bankole       PetscInt        num_comp = 0;
122c5e9980aSAdeleke O. Bankole       const PetscInt *points;
1232004e3acSAdeleke O. Bankole       PetscCall(PetscSectionGetFieldComponents(s, 0, &num_comp));
124c5e9980aSAdeleke O. Bankole       PetscCall(ISGetSize(wall_is, &num_points));
125c5e9980aSAdeleke O. Bankole       PetscCall(ISGetIndices(wall_is, &points));
126c5e9980aSAdeleke O. Bankole       for (PetscInt i = 0; i < num_points; i++) {
127c5e9980aSAdeleke O. Bankole         const PetscInt           p = points[i];
128c5e9980aSAdeleke O. Bankole         const StateConservative *r;
129c5e9980aSAdeleke O. Bankole         PetscCall(DMPlexPointLocalRead(dm, p, g, &r));
1302004e3acSAdeleke O. Bankole         PetscCall(PetscSectionGetDof(s, p, &dof));
1312004e3acSAdeleke O. Bankole         for (PetscInt node = 0; node < dof / num_comp; node++) {
132c5e9980aSAdeleke O. Bankole           for (PetscInt j = 0; j < 3; j++) {
1332004e3acSAdeleke O. Bankole             reaction_force[w * dim + j] -= r[node].momentum[j];
1342004e3acSAdeleke O. Bankole           }
135c5e9980aSAdeleke O. Bankole         }
136c5e9980aSAdeleke O. Bankole       }
137c5e9980aSAdeleke O. Bankole       PetscCall(ISRestoreIndices(wall_is, &points));
138c5e9980aSAdeleke O. Bankole     }
139c5e9980aSAdeleke O. Bankole     PetscCall(ISDestroy(&wall_is));
140c5e9980aSAdeleke O. Bankole   }
141c5e9980aSAdeleke O. Bankole   PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, reaction_force, dim * num_walls, MPIU_SCALAR, MPI_SUM, comm));
142c5e9980aSAdeleke O. Bankole   //  Restore Vectors
143c5e9980aSAdeleke O. Bankole   PetscCall(VecRestoreArrayRead(G_loc, &g));
144d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
145c5e9980aSAdeleke O. Bankole }
146c5e9980aSAdeleke O. Bankole 
147a515125bSLeila Ghaffari // Implicit time-stepper function setup
1482b916ea7SJeremy L Thompson PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) {
149a515125bSLeila Ghaffari   User         user = *(User *)user_data;
150701e5830SJames Wright   Ceed         ceed = user->ceed;
151fd969b44SJames Wright   PetscScalar  dt;
152e2f84137SJeremy L Thompson   Vec          Q_loc = user->Q_loc, Q_dot_loc = user->Q_dot_loc, G_loc;
153a515125bSLeila Ghaffari   PetscMemType q_mem_type, q_dot_mem_type, g_mem_type;
154a515125bSLeila Ghaffari 
15506f41313SJames Wright   PetscFunctionBeginUser;
156e2f84137SJeremy L Thompson   // Get local vectors
157c5e9980aSAdeleke O. Bankole   PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
158e2f84137SJeremy L Thompson 
159e2f84137SJeremy L Thompson   // Update time dependent data
160c996854bSJames Wright   PetscCall(UpdateBoundaryValues(user, Q_loc, t));
161701e5830SJames Wright   if (user->phys->solution_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_ifunction, user->phys->solution_time_label, &t));
1622b916ea7SJeremy L Thompson   PetscCall(TSGetTimeStep(ts, &dt));
163701e5830SJames Wright   if (user->phys->timestep_size_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_ifunction, user->phys->timestep_size_label, &dt));
164a515125bSLeila Ghaffari 
165a515125bSLeila Ghaffari   // Global-to-local
16606108310SJames Wright   PetscCall(DMGlobalToLocalBegin(user->dm, Q, INSERT_VALUES, Q_loc));
16706108310SJames Wright   PetscCall(DMGlobalToLocalBegin(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
16806108310SJames Wright   PetscCall(DMGlobalToLocalEnd(user->dm, Q, INSERT_VALUES, Q_loc));
16906108310SJames Wright   PetscCall(DMGlobalToLocalEnd(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
170a515125bSLeila Ghaffari 
171a515125bSLeila Ghaffari   // Place PETSc vectors in CEED vectors
172fd969b44SJames Wright   PetscCall(VecReadP2C(Q_loc, &q_mem_type, user->q_ceed));
173fd969b44SJames Wright   PetscCall(VecReadP2C(Q_dot_loc, &q_dot_mem_type, user->q_dot_ceed));
174fd969b44SJames Wright   PetscCall(VecP2C(G_loc, &g_mem_type, user->g_ceed));
175a515125bSLeila Ghaffari 
176a515125bSLeila Ghaffari   // Apply CEED operator
1777eedc94cSJames Wright   PetscCall(PetscLogEventBegin(FLUIDS_CeedOperatorApply, Q, G, 0, 0));
1787eedc94cSJames Wright   PetscCall(PetscLogGpuTimeBegin());
179b4c37c5cSJames Wright   PetscCallCeed(user->ceed, CeedOperatorApply(user->op_ifunction, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE));
1807eedc94cSJames Wright   PetscCall(PetscLogGpuTimeEnd());
1817eedc94cSJames Wright   PetscCall(PetscLogEventEnd(FLUIDS_CeedOperatorApply, Q, G, 0, 0));
182a515125bSLeila Ghaffari 
183a515125bSLeila Ghaffari   // Restore vectors
184fd969b44SJames Wright   PetscCall(VecReadC2P(user->q_ceed, q_mem_type, Q_loc));
185fd969b44SJames Wright   PetscCall(VecReadC2P(user->q_dot_ceed, q_dot_mem_type, Q_dot_loc));
186fd969b44SJames Wright   PetscCall(VecC2P(user->g_ceed, g_mem_type, G_loc));
187a515125bSLeila Ghaffari 
18801ab89c1SJames Wright   if (user->app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) {
18942454adaSJames Wright     PetscCall(SgsDDModelApplyIFunction(user, Q_loc, G_loc));
19001ab89c1SJames Wright   }
1919c678832SJames Wright 
192a515125bSLeila Ghaffari   // Local-to-Global
1932b916ea7SJeremy L Thompson   PetscCall(VecZeroEntries(G));
1942b916ea7SJeremy L Thompson   PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G));
195a515125bSLeila Ghaffari 
196a515125bSLeila Ghaffari   // Restore vectors
197c5e9980aSAdeleke O. Bankole   PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
198d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
199a515125bSLeila Ghaffari }
200a515125bSLeila Ghaffari 
2012b916ea7SJeremy L Thompson static PetscErrorCode FormPreallocation(User user, PetscBool pbdiagonal, Mat J, CeedVector *coo_values) {
202b107fddaSJed Brown   PetscCount ncoo;
203defe8520SJames Wright   PetscInt  *rows_petsc, *cols_petsc;
20471c848e3SJames Wright   CeedInt   *rows_ceed, *cols_ceed;
205b107fddaSJed Brown 
206b107fddaSJed Brown   PetscFunctionBeginUser;
207b107fddaSJed Brown   if (pbdiagonal) {
20871c848e3SJames Wright     PetscCallCeed(user->ceed, CeedOperatorLinearAssemblePointBlockDiagonalSymbolic(user->op_ijacobian, &ncoo, &rows_ceed, &cols_ceed));
209b107fddaSJed Brown   } else {
210b4c37c5cSJames Wright     PetscCallCeed(user->ceed, CeedOperatorLinearAssembleSymbolic(user->op_ijacobian, &ncoo, &rows_ceed, &cols_ceed));
21171c848e3SJames Wright   }
212defe8520SJames Wright   PetscCall(IntArrayC2P(ncoo, &rows_ceed, &rows_petsc));
213defe8520SJames Wright   PetscCall(IntArrayC2P(ncoo, &cols_ceed, &cols_petsc));
214defe8520SJames Wright   PetscCall(MatSetPreallocationCOOLocal(J, ncoo, rows_petsc, cols_petsc));
215defe8520SJames Wright   free(rows_petsc);
216defe8520SJames Wright   free(cols_petsc);
217b4c37c5cSJames Wright   PetscCallCeed(user->ceed, CeedVectorCreate(user->ceed, ncoo, coo_values));
218d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
219b107fddaSJed Brown }
220b107fddaSJed Brown 
2212b916ea7SJeremy L Thompson static PetscErrorCode FormSetValues(User user, PetscBool pbdiagonal, Mat J, CeedVector coo_values) {
222b107fddaSJed Brown   CeedMemType        mem_type = CEED_MEM_HOST;
223b107fddaSJed Brown   const PetscScalar *values;
224b107fddaSJed Brown   MatType            mat_type;
225b107fddaSJed Brown 
226b107fddaSJed Brown   PetscFunctionBeginUser;
227b107fddaSJed Brown   PetscCall(MatGetType(J, &mat_type));
2282b916ea7SJeremy L Thompson   if (strstr(mat_type, "kokkos") || strstr(mat_type, "cusparse")) mem_type = CEED_MEM_DEVICE;
229cb315d14SJames Wright   if (pbdiagonal) {
2307eedc94cSJames Wright     PetscCall(PetscLogEventBegin(FLUIDS_CeedOperatorAssemblePointBlockDiagonal, J, 0, 0, 0));
2317eedc94cSJames Wright     PetscCall(PetscLogGpuTimeBegin());
232b4c37c5cSJames Wright     PetscCallCeed(user->ceed, CeedOperatorLinearAssemblePointBlockDiagonal(user->op_ijacobian, coo_values, CEED_REQUEST_IMMEDIATE));
2337eedc94cSJames Wright     PetscCall(PetscLogGpuTimeEnd());
2347eedc94cSJames Wright     PetscCall(PetscLogEventEnd(FLUIDS_CeedOperatorAssemblePointBlockDiagonal, J, 0, 0, 0));
235b107fddaSJed Brown   } else {
2367eedc94cSJames Wright     PetscCall(PetscLogEventBegin(FLUIDS_CeedOperatorAssemble, J, 0, 0, 0));
2377eedc94cSJames Wright     PetscCall(PetscLogGpuTimeBegin());
238b4c37c5cSJames Wright     PetscCallCeed(user->ceed, CeedOperatorLinearAssemble(user->op_ijacobian, coo_values));
2397eedc94cSJames Wright     PetscCall(PetscLogGpuTimeEnd());
2407eedc94cSJames Wright     PetscCall(PetscLogEventEnd(FLUIDS_CeedOperatorAssemble, J, 0, 0, 0));
241b107fddaSJed Brown   }
242b4c37c5cSJames Wright   PetscCallCeed(user->ceed, CeedVectorGetArrayRead(coo_values, mem_type, &values));
243b107fddaSJed Brown   PetscCall(MatSetValuesCOO(J, values, INSERT_VALUES));
244b4c37c5cSJames Wright   PetscCallCeed(user->ceed, CeedVectorRestoreArrayRead(coo_values, &values));
245d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
246b107fddaSJed Brown }
247b107fddaSJed Brown 
2482b916ea7SJeremy L Thompson PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, PetscReal shift, Mat J, Mat J_pre, void *user_data) {
249f0b65372SJed Brown   User      user = *(User *)user_data;
250b4c37c5cSJames Wright   Ceed      ceed = user->ceed;
25104855949SJed Brown   PetscBool J_is_shell, J_is_mffd, J_pre_is_shell;
25206f41313SJames Wright 
253f0b65372SJed Brown   PetscFunctionBeginUser;
254b4c37c5cSJames Wright   if (user->phys->ijacobian_time_shift_label)
255b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_ijacobian, user->phys->ijacobian_time_shift_label, &shift));
25604855949SJed Brown   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATMFFD, &J_is_mffd));
257f0b65372SJed Brown   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATSHELL, &J_is_shell));
2582b916ea7SJeremy L Thompson   PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATSHELL, &J_pre_is_shell));
259f0b65372SJed Brown   if (!user->matrices_set_up) {
260f0b65372SJed Brown     if (J_is_shell) {
261f9028c3cSJames Wright       OperatorApplyContext op_ijacobian_ctx;
262f9028c3cSJames Wright       OperatorApplyContextCreate(user->dm, user->dm, user->ceed, user->op_ijacobian, user->q_ceed, user->g_ceed, user->Q_dot_loc, NULL,
263f9028c3cSJames Wright                                  &op_ijacobian_ctx);
264f9028c3cSJames Wright       PetscCall(MatShellSetContext(J, op_ijacobian_ctx));
265f9028c3cSJames Wright       PetscCall(MatShellSetContextDestroy(J, (PetscErrorCode(*)(void *))OperatorApplyContextDestroy));
266f9028c3cSJames Wright       PetscCall(MatShellSetOperation(J, MATOP_MULT, (void (*)(void))MatMult_Ceed));
267f9028c3cSJames Wright       PetscCall(MatShellSetOperation(J, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiag_Ceed));
268f0b65372SJed Brown       PetscCall(MatSetUp(J));
269f0b65372SJed Brown     }
270f0b65372SJed Brown     if (!J_pre_is_shell) {
2712b916ea7SJeremy L Thompson       PetscCall(FormPreallocation(user, user->app_ctx->pmat_pbdiagonal, J_pre, &user->coo_values_pmat));
272b107fddaSJed Brown     }
27304855949SJed Brown     if (J != J_pre && !J_is_shell && !J_is_mffd) {
274b107fddaSJed Brown       PetscCall(FormPreallocation(user, PETSC_FALSE, J, &user->coo_values_amat));
275b107fddaSJed Brown     }
276f0b65372SJed Brown     user->matrices_set_up = true;
277f0b65372SJed Brown   }
278f0b65372SJed Brown   if (!J_pre_is_shell) {
2792b916ea7SJeremy L Thompson     PetscCall(FormSetValues(user, user->app_ctx->pmat_pbdiagonal, J_pre, user->coo_values_pmat));
280f0b65372SJed Brown   }
28104855949SJed Brown   if (user->coo_values_amat) {
28204855949SJed Brown     PetscCall(FormSetValues(user, PETSC_FALSE, J, user->coo_values_amat));
28304855949SJed Brown   } else if (J_is_mffd) {
28404855949SJed Brown     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
28504855949SJed Brown     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
28604855949SJed Brown   }
287d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
288f0b65372SJed Brown }
289f0b65372SJed Brown 
2902b916ea7SJeremy L Thompson PetscErrorCode WriteOutput(User user, Vec Q, PetscInt step_no, PetscScalar time) {
291a515125bSLeila Ghaffari   Vec         Q_loc;
292a515125bSLeila Ghaffari   char        file_path[PETSC_MAX_PATH_LEN];
293a515125bSLeila Ghaffari   PetscViewer viewer;
294a515125bSLeila Ghaffari 
29506f41313SJames Wright   PetscFunctionBeginUser;
296852e5969SJed Brown   if (user->app_ctx->checkpoint_vtk) {
297a515125bSLeila Ghaffari     // Set up output
2987538d537SJames Wright     PetscCall(DMGetLocalVector(user->dm, &Q_loc));
2997538d537SJames Wright     PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec"));
3007538d537SJames Wright     PetscCall(VecZeroEntries(Q_loc));
3017538d537SJames Wright     PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
302a515125bSLeila Ghaffari 
303a515125bSLeila Ghaffari     // Output
304852e5969SJed Brown     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
3057538d537SJames Wright 
3062b916ea7SJeremy L Thompson     PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path, FILE_MODE_WRITE, &viewer));
3077538d537SJames Wright     PetscCall(VecView(Q_loc, viewer));
3087538d537SJames Wright     PetscCall(PetscViewerDestroy(&viewer));
309a515125bSLeila Ghaffari     if (user->dm_viz) {
310a515125bSLeila Ghaffari       Vec         Q_refined, Q_refined_loc;
311a515125bSLeila Ghaffari       char        file_path_refined[PETSC_MAX_PATH_LEN];
312a515125bSLeila Ghaffari       PetscViewer viewer_refined;
313a515125bSLeila Ghaffari 
3147538d537SJames Wright       PetscCall(DMGetGlobalVector(user->dm_viz, &Q_refined));
3157538d537SJames Wright       PetscCall(DMGetLocalVector(user->dm_viz, &Q_refined_loc));
3167538d537SJames Wright       PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined"));
3177538d537SJames Wright 
3187538d537SJames Wright       PetscCall(MatInterpolate(user->interp_viz, Q, Q_refined));
3197538d537SJames Wright       PetscCall(VecZeroEntries(Q_refined_loc));
3202b916ea7SJeremy L Thompson       PetscCall(DMGlobalToLocal(user->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc));
3217538d537SJames Wright 
322852e5969SJed Brown       PetscCall(
323852e5969SJed Brown           PetscSNPrintf(file_path_refined, sizeof file_path_refined, "%s/nsrefined-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
3247538d537SJames Wright 
3252b916ea7SJeremy L Thompson       PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined), file_path_refined, FILE_MODE_WRITE, &viewer_refined));
3267538d537SJames Wright       PetscCall(VecView(Q_refined_loc, viewer_refined));
3277538d537SJames Wright       PetscCall(DMRestoreLocalVector(user->dm_viz, &Q_refined_loc));
3287538d537SJames Wright       PetscCall(DMRestoreGlobalVector(user->dm_viz, &Q_refined));
3297538d537SJames Wright       PetscCall(PetscViewerDestroy(&viewer_refined));
330a515125bSLeila Ghaffari     }
3317538d537SJames Wright     PetscCall(DMRestoreLocalVector(user->dm, &Q_loc));
332852e5969SJed Brown   }
333a515125bSLeila Ghaffari 
334a515125bSLeila Ghaffari   // Save data in a binary file for continuation of simulations
33591a36801SJames Wright   if (user->app_ctx->add_stepnum2bin) {
336852e5969SJed Brown     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", user->app_ctx->output_dir, step_no));
33791a36801SJames Wright   } else {
3382b916ea7SJeremy L Thompson     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", user->app_ctx->output_dir));
33991a36801SJames Wright   }
3402b916ea7SJeremy L Thompson   PetscCall(PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, &viewer));
3417538d537SJames Wright 
342e1233009SJames Wright   PetscInt32 token = PetscDefined(USE_64BIT_INDICES) ? FLUIDS_FILE_TOKEN_64 : FLUIDS_FILE_TOKEN_32;
343e1233009SJames Wright   PetscCall(PetscViewerBinaryWrite(viewer, &token, 1, PETSC_INT32));
3449293eaa1SJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &step_no, 1, PETSC_INT));
3459293eaa1SJed Brown   time /= user->units->second;  // Dimensionalize time back
3469293eaa1SJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL));
3477538d537SJames Wright   PetscCall(VecView(Q, viewer));
3487538d537SJames Wright   PetscCall(PetscViewerDestroy(&viewer));
349d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3507538d537SJames Wright }
3517538d537SJames Wright 
352c5e9980aSAdeleke O. Bankole // CSV Monitor
353c5e9980aSAdeleke O. Bankole PetscErrorCode TSMonitor_WallForce(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
354c5e9980aSAdeleke O. Bankole   User              user = ctx;
355c5e9980aSAdeleke O. Bankole   Vec               G_loc;
356c5e9980aSAdeleke O. Bankole   PetscInt          num_wall = user->app_ctx->wall_forces.num_wall, dim = 3;
357c5e9980aSAdeleke O. Bankole   const PetscInt   *walls  = user->app_ctx->wall_forces.walls;
358c5e9980aSAdeleke O. Bankole   PetscViewer       viewer = user->app_ctx->wall_forces.viewer;
359c5e9980aSAdeleke O. Bankole   PetscViewerFormat format = user->app_ctx->wall_forces.viewer_format;
360c5e9980aSAdeleke O. Bankole   PetscScalar      *reaction_force;
361c5e9980aSAdeleke O. Bankole   PetscBool         iascii;
362c5e9980aSAdeleke O. Bankole 
363c5e9980aSAdeleke O. Bankole   PetscFunctionBeginUser;
364d949ddfcSJames Wright   if (!viewer) PetscFunctionReturn(PETSC_SUCCESS);
365c5e9980aSAdeleke O. Bankole   PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
366c5e9980aSAdeleke O. Bankole   PetscCall(PetscMalloc1(num_wall * dim, &reaction_force));
367c5e9980aSAdeleke O. Bankole   PetscCall(Surface_Forces_NS(user->dm, G_loc, num_wall, walls, reaction_force));
368c5e9980aSAdeleke O. Bankole   PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
369c5e9980aSAdeleke O. Bankole 
370c5e9980aSAdeleke O. Bankole   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
371c5e9980aSAdeleke O. Bankole 
372c5e9980aSAdeleke O. Bankole   if (iascii) {
373c5e9980aSAdeleke O. Bankole     if (format == PETSC_VIEWER_ASCII_CSV && !user->app_ctx->wall_forces.header_written) {
374c5e9980aSAdeleke O. Bankole       PetscCall(PetscViewerASCIIPrintf(viewer, "Step,Time,Wall,ForceX,ForceY,ForceZ\n"));
375c5e9980aSAdeleke O. Bankole       user->app_ctx->wall_forces.header_written = PETSC_TRUE;
376c5e9980aSAdeleke O. Bankole     }
377c5e9980aSAdeleke O. Bankole     for (PetscInt w = 0; w < num_wall; w++) {
378c5e9980aSAdeleke O. Bankole       PetscInt wall = walls[w];
379c5e9980aSAdeleke O. Bankole       if (format == PETSC_VIEWER_ASCII_CSV) {
380c5e9980aSAdeleke O. Bankole         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT ",%g,%" PetscInt_FMT ",%g,%g,%g\n", step_no, time, wall,
381c5e9980aSAdeleke O. Bankole                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
382c5e9980aSAdeleke O. Bankole 
383c5e9980aSAdeleke O. Bankole       } else {
384c5e9980aSAdeleke O. Bankole         PetscCall(PetscViewerASCIIPrintf(viewer, "Wall %" PetscInt_FMT " Forces: Force_x = %12g, Force_y = %12g, Force_z = %12g\n", wall,
385c5e9980aSAdeleke O. Bankole                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
386c5e9980aSAdeleke O. Bankole       }
387c5e9980aSAdeleke O. Bankole     }
388c5e9980aSAdeleke O. Bankole   }
389c5e9980aSAdeleke O. Bankole   PetscCall(PetscFree(reaction_force));
390d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
391c5e9980aSAdeleke O. Bankole }
392c5e9980aSAdeleke O. Bankole 
3937538d537SJames Wright // User provided TS Monitor
3942b916ea7SJeremy L Thompson PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
3957538d537SJames Wright   User user = ctx;
3967538d537SJames Wright 
39706f41313SJames Wright   PetscFunctionBeginUser;
398852e5969SJed Brown   // Print every 'checkpoint_interval' steps
399c539088bSJames Wright   if (user->app_ctx->checkpoint_interval <= 0 || step_no % user->app_ctx->checkpoint_interval != 0 ||
400e419654dSJeremy L Thompson       (user->app_ctx->cont_steps == step_no && step_no != 0)) {
401d949ddfcSJames Wright     PetscFunctionReturn(PETSC_SUCCESS);
402e419654dSJeremy L Thompson   }
4037538d537SJames Wright 
4047538d537SJames Wright   PetscCall(WriteOutput(user, Q, step_no, time));
405d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
406a515125bSLeila Ghaffari }
407a515125bSLeila Ghaffari 
408a515125bSLeila Ghaffari // TS: Create, setup, and solve
4092b916ea7SJeremy L Thompson PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, Vec *Q, PetscScalar *f_time, TS *ts) {
410a515125bSLeila Ghaffari   MPI_Comm    comm = user->comm;
411a515125bSLeila Ghaffari   TSAdapt     adapt;
412a515125bSLeila Ghaffari   PetscScalar final_time;
413a515125bSLeila Ghaffari 
41406f41313SJames Wright   PetscFunctionBeginUser;
4152b916ea7SJeremy L Thompson   PetscCall(TSCreate(comm, ts));
4162b916ea7SJeremy L Thompson   PetscCall(TSSetDM(*ts, dm));
417632a41e1SJames Wright   PetscCall(TSSetApplicationContext(*ts, user));
418a515125bSLeila Ghaffari   if (phys->implicit) {
4192b916ea7SJeremy L Thompson     PetscCall(TSSetType(*ts, TSBDF));
420a515125bSLeila Ghaffari     if (user->op_ifunction) {
4212b916ea7SJeremy L Thompson       PetscCall(TSSetIFunction(*ts, NULL, IFunction_NS, &user));
422a515125bSLeila Ghaffari     } else {  // Implicit integrators can fall back to using an RHSFunction
4232b916ea7SJeremy L Thompson       PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
424a515125bSLeila Ghaffari     }
425f0b65372SJed Brown     if (user->op_ijacobian) {
4262b916ea7SJeremy L Thompson       PetscCall(DMTSSetIJacobian(dm, FormIJacobian_NS, &user));
427b107fddaSJed Brown       if (app_ctx->amat_type) {
428b107fddaSJed Brown         Mat Pmat, Amat;
4292b916ea7SJeremy L Thompson         PetscCall(DMCreateMatrix(dm, &Pmat));
4302b916ea7SJeremy L Thompson         PetscCall(DMSetMatType(dm, app_ctx->amat_type));
4312b916ea7SJeremy L Thompson         PetscCall(DMCreateMatrix(dm, &Amat));
4322b916ea7SJeremy L Thompson         PetscCall(TSSetIJacobian(*ts, Amat, Pmat, NULL, NULL));
4332b916ea7SJeremy L Thompson         PetscCall(MatDestroy(&Amat));
4342b916ea7SJeremy L Thompson         PetscCall(MatDestroy(&Pmat));
435b107fddaSJed Brown       }
436f0b65372SJed Brown     }
437a515125bSLeila Ghaffari   } else {
438da5fe0e4SJames Wright     PetscCheck(user->op_rhs_ctx, comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction");
4392b916ea7SJeremy L Thompson     PetscCall(TSSetType(*ts, TSRK));
4402b916ea7SJeremy L Thompson     PetscCall(TSRKSetType(*ts, TSRK5F));
4412b916ea7SJeremy L Thompson     PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
442a515125bSLeila Ghaffari   }
4432b916ea7SJeremy L Thompson   PetscCall(TSSetMaxTime(*ts, 500. * user->units->second));
4442b916ea7SJeremy L Thompson   PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER));
44522387d3aSJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(TSSetErrorIfStepFails(*ts, PETSC_FALSE));
4462b916ea7SJeremy L Thompson   PetscCall(TSSetTimeStep(*ts, 1.e-2 * user->units->second));
4472b916ea7SJeremy L Thompson   PetscCall(TSGetAdapt(*ts, &adapt));
4482b916ea7SJeremy L Thompson   PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * user->units->second, 1.e2 * user->units->second));
4492b916ea7SJeremy L Thompson   PetscCall(TSSetFromOptions(*ts));
45091f639d2SJames Wright   user->time_bc_set = -1.0;   // require all BCs be updated
451c26b555cSJames Wright   if (app_ctx->cont_steps) {  // continue from previous timestep data
452a515125bSLeila Ghaffari     PetscInt    count;
453a515125bSLeila Ghaffari     PetscViewer viewer;
4542b916ea7SJeremy L Thompson 
4559293eaa1SJed Brown     if (app_ctx->cont_time <= 0) {  // Legacy files did not include step number and time
4562b916ea7SJeremy L Thompson       PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_time_file, FILE_MODE_READ, &viewer));
4579293eaa1SJed Brown       PetscCall(PetscViewerBinaryRead(viewer, &app_ctx->cont_time, 1, &count, PETSC_REAL));
4582b916ea7SJeremy L Thompson       PetscCall(PetscViewerDestroy(&viewer));
4599293eaa1SJed Brown       PetscCheck(app_ctx->cont_steps != -1, comm, PETSC_ERR_ARG_INCOMP,
4609293eaa1SJed Brown                  "-continue step number not specified, but checkpoint file does not contain a step number (likely written by older code version)");
4619293eaa1SJed Brown     }
4629293eaa1SJed Brown     PetscCall(TSSetTime(*ts, app_ctx->cont_time * user->units->second));
46374a6f4ddSJed Brown     PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps));
464a515125bSLeila Ghaffari   }
4650e1e9333SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) {
4662b916ea7SJeremy L Thompson     PetscCall(TSMonitorSet(*ts, TSMonitor_NS, user, NULL));
4670e1e9333SJames Wright   }
468c5e9980aSAdeleke O. Bankole   if (app_ctx->wall_forces.viewer) {
469c5e9980aSAdeleke O. Bankole     PetscCall(TSMonitorSet(*ts, TSMonitor_WallForce, user, NULL));
470c5e9980aSAdeleke O. Bankole   }
471c931fa59SJames Wright   if (app_ctx->turb_spanstats_enable) {
47291933550SJames Wright     PetscCall(TSMonitorSet(*ts, TSMonitor_TurbulenceStatistics, user, NULL));
473b8daee98SJames Wright     CeedScalar previous_time = app_ctx->cont_time * user->units->second;
474b4c37c5cSJames Wright     PetscCallCeed(user->ceed,
475b4c37c5cSJames Wright                   CeedOperatorSetContextDouble(user->spanstats.op_stats_collect_ctx->op, user->spanstats.previous_time_label, &previous_time));
476b0488d1fSJames Wright   }
47788b07121SJames Wright   if (app_ctx->diff_filter_monitor) PetscCall(TSMonitorSet(*ts, TSMonitor_DifferentialFilter, user, NULL));
478a515125bSLeila Ghaffari 
4791c17f66aSJames Wright   if (app_ctx->sgs_train_enable) {
4801c17f66aSJames Wright     PetscCall(TSMonitorSet(*ts, TSMonitor_SGS_DD_Training, user, NULL));
4811c17f66aSJames Wright     PetscCall(TSSetPostStep(*ts, TSPostStep_SGS_DD_Training));
4821c17f66aSJames Wright   }
483a515125bSLeila Ghaffari   // Solve
48474a6f4ddSJed Brown   PetscReal start_time;
48574a6f4ddSJed Brown   PetscInt  start_step;
4862b916ea7SJeremy L Thompson   PetscCall(TSGetTime(*ts, &start_time));
48774a6f4ddSJed Brown   PetscCall(TSGetStepNumber(*ts, &start_step));
48891982731SJeremy L Thompson 
489df4304b5SJed Brown   PetscCall(PetscLogDefaultBegin());  // So we can use PetscLogStageGetPerfInfo without -log_view
49091982731SJeremy L Thompson   PetscPreLoadBegin(PETSC_FALSE, "Fluids Solve");
49191982731SJeremy L Thompson   PetscCall(TSSetTime(*ts, start_time));
49274a6f4ddSJed Brown   PetscCall(TSSetStepNumber(*ts, start_step));
49391982731SJeremy L Thompson   if (PetscPreLoadingOn) {
49491982731SJeremy L Thompson     // LCOV_EXCL_START
49591982731SJeremy L Thompson     SNES      snes;
49691982731SJeremy L Thompson     Vec       Q_preload;
49791982731SJeremy L Thompson     PetscReal rtol;
49891982731SJeremy L Thompson     PetscCall(VecDuplicate(*Q, &Q_preload));
49991982731SJeremy L Thompson     PetscCall(VecCopy(*Q, Q_preload));
50091982731SJeremy L Thompson     PetscCall(TSGetSNES(*ts, &snes));
50191982731SJeremy L Thompson     PetscCall(SNESGetTolerances(snes, NULL, &rtol, NULL, NULL, NULL));
5022b916ea7SJeremy L Thompson     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, .99, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
50322c1b34eSJames Wright     PetscCall(TSSetSolution(*ts, Q_preload));
50491982731SJeremy L Thompson     PetscCall(TSStep(*ts));
5052b916ea7SJeremy L Thompson     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, rtol, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
50691982731SJeremy L Thompson     PetscCall(VecDestroy(&Q_preload));
50791982731SJeremy L Thompson     // LCOV_EXCL_STOP
50891982731SJeremy L Thompson   } else {
5092b916ea7SJeremy L Thompson     PetscCall(PetscBarrier((PetscObject)*ts));
5102b916ea7SJeremy L Thompson     PetscCall(TSSolve(*ts, *Q));
51191982731SJeremy L Thompson   }
51291982731SJeremy L Thompson   PetscPreLoadEnd();
51391982731SJeremy L Thompson 
51491982731SJeremy L Thompson   PetscCall(TSGetSolveTime(*ts, &final_time));
515a515125bSLeila Ghaffari   *f_time = final_time;
51691982731SJeremy L Thompson 
5170e1e9333SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) {
5187538d537SJames Wright     PetscInt step_no;
5197538d537SJames Wright     PetscCall(TSGetStepNumber(*ts, &step_no));
520b0488d1fSJames Wright     if (user->app_ctx->checkpoint_interval > 0 || user->app_ctx->checkpoint_interval == -1) {
5217538d537SJames Wright       PetscCall(WriteOutput(user, *Q, step_no, final_time));
5227538d537SJames Wright     }
5237538d537SJames Wright 
5247eedc94cSJames Wright     PetscLogStage      stage_id;
525df4304b5SJed Brown     PetscEventPerfInfo stage_perf;
52691982731SJeremy L Thompson 
52791982731SJeremy L Thompson     PetscCall(PetscLogStageGetId("Fluids Solve", &stage_id));
528df4304b5SJed Brown     PetscCall(PetscLogStageGetPerfInfo(stage_id, &stage_perf));
529df4304b5SJed Brown     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time taken for solution (sec): %g\n", stage_perf.time));
530a515125bSLeila Ghaffari   }
531d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
532a515125bSLeila Ghaffari }
533