xref: /libCEED/examples/fluids/src/setupts.c (revision f17d818ddb5f0fdb5c67183d33cc31a0016906ee)
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 
26*f17d818dSJames Wright   PetscFunctionBeginUser;
2777841947SLeila Ghaffari   // CEED Restriction
28a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q));
29a424bcd0SJames Wright   PetscCallCeed(ceed, 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
35a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass));
36a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
37a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data));
38a424bcd0SJames Wright   PetscCallCeed(ceed, 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));
52a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
53a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
54ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
5577841947SLeila Ghaffari }
5677841947SLeila Ghaffari 
57a0b9a424SJames Wright // Insert Boundary values if it's a new time
58a0b9a424SJames Wright PetscErrorCode UpdateBoundaryValues(User user, Vec Q_loc, PetscReal t) {
59a0b9a424SJames Wright   PetscFunctionBeginUser;
60a0b9a424SJames Wright   if (user->time_bc_set != t) {
61a0b9a424SJames Wright     PetscCall(DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL));
62a0b9a424SJames Wright     user->time_bc_set = t;
63a0b9a424SJames Wright   }
64ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
65a0b9a424SJames Wright }
66a0b9a424SJames Wright 
6715637395SJames Wright // @brief Update the context label value to new value if necessary.
6815637395SJames Wright // @note This only supports labels with scalar label values (ie. not arrays)
69a424bcd0SJames Wright PetscErrorCode UpdateContextLabel(Ceed ceed, MPI_Comm comm, PetscScalar update_value, CeedOperator op, CeedContextFieldLabel label) {
7015637395SJames Wright   PetscScalar label_value;
71e42c1df6SJames Wright 
7215637395SJames Wright   PetscFunctionBeginUser;
7315637395SJames Wright   PetscCheck(label, comm, PETSC_ERR_ARG_BADPTR, "Label should be non-NULL");
7415637395SJames Wright 
7515637395SJames Wright   {
7615637395SJames Wright     size_t             num_elements;
7715637395SJames Wright     const PetscScalar *label_values;
78a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContextDoubleRead(op, label, &num_elements, &label_values));
7915637395SJames Wright     PetscCheck(num_elements == 1, comm, PETSC_ERR_SUP, "%s does not support labels with more than 1 value. Label has %zu values", __func__,
8015637395SJames Wright                num_elements);
8115637395SJames Wright     label_value = *label_values;
82a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorRestoreContextDoubleRead(op, label, &label_values));
83e42c1df6SJames Wright   }
8415637395SJames Wright 
8515637395SJames Wright   if (label_value != update_value) {
86a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetContextDouble(op, label, &update_value));
87e42c1df6SJames Wright   }
88ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
89e42c1df6SJames Wright }
90e42c1df6SJames Wright 
9177841947SLeila Ghaffari // RHS (Explicit time-stepper) function setup
9277841947SLeila Ghaffari //   This is the RHS of the ODE, given as u_t = G(t,u)
9377841947SLeila Ghaffari //   This function takes in a state vector Q and writes into G
9477841947SLeila Ghaffari PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) {
9577841947SLeila Ghaffari   User        user = *(User *)user_data;
9615637395SJames Wright   MPI_Comm    comm = PetscObjectComm((PetscObject)ts);
97c798d55aSJames Wright   PetscScalar dt;
989ad5e8e4SJames Wright   Vec         Q_loc = user->Q_loc;
9977841947SLeila Ghaffari 
100*f17d818dSJames Wright   PetscFunctionBeginUser;
1015e82a6e1SJeremy L Thompson   // Update time dependent data
102a0b9a424SJames Wright   PetscCall(UpdateBoundaryValues(user, Q_loc, t));
103a424bcd0SJames Wright   if (user->phys->solution_time_label) PetscCall(UpdateContextLabel(user->ceed, comm, t, user->op_rhs_ctx->op, user->phys->solution_time_label));
1042b730f8bSJeremy L Thompson   PetscCall(TSGetTimeStep(ts, &dt));
105a424bcd0SJames Wright   if (user->phys->timestep_size_label) PetscCall(UpdateContextLabel(user->ceed, comm, dt, user->op_rhs_ctx->op, user->phys->timestep_size_label));
10677841947SLeila Ghaffari 
1079ad5e8e4SJames Wright   PetscCall(ApplyCeedOperatorGlobalToGlobal(Q, G, user->op_rhs_ctx));
10877841947SLeila Ghaffari 
109d6e67e47SJames Wright   // Inverse of the lumped mass matrix
110186595e6SJames Wright   PetscCall(VecPointwiseMult(G, G, user->M_inv));
111ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
11277841947SLeila Ghaffari }
11377841947SLeila Ghaffari 
114ca69d878SAdeleke O. Bankole // Surface forces function setup
115ca69d878SAdeleke O. Bankole static PetscErrorCode Surface_Forces_NS(DM dm, Vec G_loc, PetscInt num_walls, const PetscInt walls[], PetscScalar *reaction_force) {
116ca69d878SAdeleke O. Bankole   DMLabel            face_label;
117ca69d878SAdeleke O. Bankole   const PetscScalar *g;
118d6734f85SAdeleke O. Bankole   PetscInt           dof, dim = 3;
119ca69d878SAdeleke O. Bankole   MPI_Comm           comm;
120d6734f85SAdeleke O. Bankole   PetscSection       s;
121ca69d878SAdeleke O. Bankole 
122ca69d878SAdeleke O. Bankole   PetscFunctionBeginUser;
123ca69d878SAdeleke O. Bankole   PetscCall(PetscArrayzero(reaction_force, num_walls * dim));
124ca69d878SAdeleke O. Bankole   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
125ca69d878SAdeleke O. Bankole   PetscCall(DMGetLabel(dm, "Face Sets", &face_label));
126ca69d878SAdeleke O. Bankole   PetscCall(VecGetArrayRead(G_loc, &g));
127ca69d878SAdeleke O. Bankole   for (PetscInt w = 0; w < num_walls; w++) {
128ca69d878SAdeleke O. Bankole     const PetscInt wall = walls[w];
129ca69d878SAdeleke O. Bankole     IS             wall_is;
130d6734f85SAdeleke O. Bankole     PetscCall(DMGetLocalSection(dm, &s));
131ca69d878SAdeleke O. Bankole     PetscCall(DMLabelGetStratumIS(face_label, wall, &wall_is));
132ca69d878SAdeleke O. Bankole     if (wall_is) {  // There exist such points on this process
133ca69d878SAdeleke O. Bankole       PetscInt        num_points;
134d6734f85SAdeleke O. Bankole       PetscInt        num_comp = 0;
135ca69d878SAdeleke O. Bankole       const PetscInt *points;
136d6734f85SAdeleke O. Bankole       PetscCall(PetscSectionGetFieldComponents(s, 0, &num_comp));
137ca69d878SAdeleke O. Bankole       PetscCall(ISGetSize(wall_is, &num_points));
138ca69d878SAdeleke O. Bankole       PetscCall(ISGetIndices(wall_is, &points));
139ca69d878SAdeleke O. Bankole       for (PetscInt i = 0; i < num_points; i++) {
140ca69d878SAdeleke O. Bankole         const PetscInt           p = points[i];
141ca69d878SAdeleke O. Bankole         const StateConservative *r;
142ca69d878SAdeleke O. Bankole         PetscCall(DMPlexPointLocalRead(dm, p, g, &r));
143d6734f85SAdeleke O. Bankole         PetscCall(PetscSectionGetDof(s, p, &dof));
144d6734f85SAdeleke O. Bankole         for (PetscInt node = 0; node < dof / num_comp; node++) {
145ca69d878SAdeleke O. Bankole           for (PetscInt j = 0; j < 3; j++) {
146d6734f85SAdeleke O. Bankole             reaction_force[w * dim + j] -= r[node].momentum[j];
147d6734f85SAdeleke O. Bankole           }
148ca69d878SAdeleke O. Bankole         }
149ca69d878SAdeleke O. Bankole       }
150ca69d878SAdeleke O. Bankole       PetscCall(ISRestoreIndices(wall_is, &points));
151ca69d878SAdeleke O. Bankole     }
152ca69d878SAdeleke O. Bankole     PetscCall(ISDestroy(&wall_is));
153ca69d878SAdeleke O. Bankole   }
154ca69d878SAdeleke O. Bankole   PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, reaction_force, dim * num_walls, MPIU_SCALAR, MPI_SUM, comm));
155ca69d878SAdeleke O. Bankole   //  Restore Vectors
156ca69d878SAdeleke O. Bankole   PetscCall(VecRestoreArrayRead(G_loc, &g));
157ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
158ca69d878SAdeleke O. Bankole }
159ca69d878SAdeleke O. Bankole 
16077841947SLeila Ghaffari // Implicit time-stepper function setup
1612b730f8bSJeremy L Thompson PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) {
16277841947SLeila Ghaffari   User         user = *(User *)user_data;
16315637395SJames Wright   MPI_Comm     comm = PetscObjectComm((PetscObject)ts);
164c798d55aSJames Wright   PetscScalar  dt;
1655e82a6e1SJeremy L Thompson   Vec          Q_loc = user->Q_loc, Q_dot_loc = user->Q_dot_loc, G_loc;
16677841947SLeila Ghaffari   PetscMemType q_mem_type, q_dot_mem_type, g_mem_type;
16777841947SLeila Ghaffari 
168*f17d818dSJames Wright   PetscFunctionBeginUser;
1695e82a6e1SJeremy L Thompson   // Get local vectors
170ca69d878SAdeleke O. Bankole   PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
1715e82a6e1SJeremy L Thompson 
1725e82a6e1SJeremy L Thompson   // Update time dependent data
173a0b9a424SJames Wright   PetscCall(UpdateBoundaryValues(user, Q_loc, t));
174a424bcd0SJames Wright   if (user->phys->solution_time_label) PetscCall(UpdateContextLabel(user->ceed, comm, t, user->op_ifunction, user->phys->solution_time_label));
1752b730f8bSJeremy L Thompson   PetscCall(TSGetTimeStep(ts, &dt));
176a424bcd0SJames Wright   if (user->phys->timestep_size_label) PetscCall(UpdateContextLabel(user->ceed, comm, dt, user->op_ifunction, user->phys->timestep_size_label));
17777841947SLeila Ghaffari 
17877841947SLeila Ghaffari   // Global-to-local
179878eb0e7SJames Wright   PetscCall(DMGlobalToLocalBegin(user->dm, Q, INSERT_VALUES, Q_loc));
180878eb0e7SJames Wright   PetscCall(DMGlobalToLocalBegin(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
181878eb0e7SJames Wright   PetscCall(DMGlobalToLocalEnd(user->dm, Q, INSERT_VALUES, Q_loc));
182878eb0e7SJames Wright   PetscCall(DMGlobalToLocalEnd(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
18377841947SLeila Ghaffari 
18477841947SLeila Ghaffari   // Place PETSc vectors in CEED vectors
185c798d55aSJames Wright   PetscCall(VecReadP2C(Q_loc, &q_mem_type, user->q_ceed));
186c798d55aSJames Wright   PetscCall(VecReadP2C(Q_dot_loc, &q_dot_mem_type, user->q_dot_ceed));
187c798d55aSJames Wright   PetscCall(VecP2C(G_loc, &g_mem_type, user->g_ceed));
18877841947SLeila Ghaffari 
18977841947SLeila Ghaffari   // Apply CEED operator
19075d1798cSJames Wright   PetscCall(PetscLogEventBegin(FLUIDS_CeedOperatorApply, Q, G, 0, 0));
19175d1798cSJames Wright   PetscCall(PetscLogGpuTimeBegin());
192a424bcd0SJames Wright   PetscCallCeed(user->ceed, CeedOperatorApply(user->op_ifunction, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE));
19375d1798cSJames Wright   PetscCall(PetscLogGpuTimeEnd());
19475d1798cSJames Wright   PetscCall(PetscLogEventEnd(FLUIDS_CeedOperatorApply, Q, G, 0, 0));
19577841947SLeila Ghaffari 
19677841947SLeila Ghaffari   // Restore vectors
197c798d55aSJames Wright   PetscCall(VecReadC2P(user->q_ceed, q_mem_type, Q_loc));
198c798d55aSJames Wright   PetscCall(VecReadC2P(user->q_dot_ceed, q_dot_mem_type, Q_dot_loc));
199c798d55aSJames Wright   PetscCall(VecC2P(user->g_ceed, g_mem_type, G_loc));
20077841947SLeila Ghaffari 
20119ffbc25SJames Wright   if (user->app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) {
202cbef7084SJames Wright     PetscCall(SgsDDModelApplyIFunction(user, Q_loc, G_loc));
20319ffbc25SJames Wright   }
204be34b3b7SJames Wright 
20577841947SLeila Ghaffari   // Local-to-Global
2062b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(G));
2072b730f8bSJeremy L Thompson   PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G));
20877841947SLeila Ghaffari 
20977841947SLeila Ghaffari   // Restore vectors
210ca69d878SAdeleke O. Bankole   PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
211ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
21277841947SLeila Ghaffari }
21377841947SLeila Ghaffari 
2142b730f8bSJeremy L Thompson static PetscErrorCode FormPreallocation(User user, PetscBool pbdiagonal, Mat J, CeedVector *coo_values) {
215544be873SJed Brown   PetscCount ncoo;
216457e73b2SJames Wright   PetscInt  *rows_petsc, *cols_petsc;
217544be873SJed Brown 
218544be873SJed Brown   PetscFunctionBeginUser;
219544be873SJed Brown   if (pbdiagonal) {
220544be873SJed Brown     CeedSize l_size;
221a424bcd0SJames Wright     PetscCallCeed(user->ceed, CeedOperatorGetActiveVectorLengths(user->op_ijacobian, &l_size, NULL));
222544be873SJed Brown     ncoo       = l_size * 5;
223457e73b2SJames Wright     rows_petsc = malloc(ncoo * sizeof(rows_petsc));
224457e73b2SJames Wright     cols_petsc = malloc(ncoo * sizeof(cols_petsc));
225544be873SJed Brown     for (PetscCount n = 0; n < l_size / 5; n++) {
226544be873SJed Brown       for (PetscInt i = 0; i < 5; i++) {
227544be873SJed Brown         for (PetscInt j = 0; j < 5; j++) {
228457e73b2SJames Wright           rows_petsc[(n * 5 + i) * 5 + j] = n * 5 + i;
229457e73b2SJames Wright           cols_petsc[(n * 5 + i) * 5 + j] = n * 5 + j;
230544be873SJed Brown         }
231544be873SJed Brown       }
232544be873SJed Brown     }
233544be873SJed Brown   } else {
234457e73b2SJames Wright     CeedInt *rows_ceed, *cols_ceed;
235a424bcd0SJames Wright     PetscCallCeed(user->ceed, CeedOperatorLinearAssembleSymbolic(user->op_ijacobian, &ncoo, &rows_ceed, &cols_ceed));
236457e73b2SJames Wright     PetscCall(IntArrayC2P(ncoo, &rows_ceed, &rows_petsc));
237457e73b2SJames Wright     PetscCall(IntArrayC2P(ncoo, &cols_ceed, &cols_petsc));
238544be873SJed Brown   }
239457e73b2SJames Wright   PetscCall(MatSetPreallocationCOOLocal(J, ncoo, rows_petsc, cols_petsc));
240457e73b2SJames Wright   free(rows_petsc);
241457e73b2SJames Wright   free(cols_petsc);
242a424bcd0SJames Wright   PetscCallCeed(user->ceed, CeedVectorCreate(user->ceed, ncoo, coo_values));
243ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
244544be873SJed Brown }
245544be873SJed Brown 
2462b730f8bSJeremy L Thompson static PetscErrorCode FormSetValues(User user, PetscBool pbdiagonal, Mat J, CeedVector coo_values) {
247544be873SJed Brown   CeedMemType        mem_type = CEED_MEM_HOST;
248544be873SJed Brown   const PetscScalar *values;
249544be873SJed Brown   MatType            mat_type;
250544be873SJed Brown 
251544be873SJed Brown   PetscFunctionBeginUser;
252544be873SJed Brown   PetscCall(MatGetType(J, &mat_type));
2532b730f8bSJeremy L Thompson   if (strstr(mat_type, "kokkos") || strstr(mat_type, "cusparse")) mem_type = CEED_MEM_DEVICE;
254fbb546ddSJames Wright   if (pbdiagonal) {
25575d1798cSJames Wright     PetscCall(PetscLogEventBegin(FLUIDS_CeedOperatorAssemblePointBlockDiagonal, J, 0, 0, 0));
25675d1798cSJames Wright     PetscCall(PetscLogGpuTimeBegin());
257a424bcd0SJames Wright     PetscCallCeed(user->ceed, CeedOperatorLinearAssemblePointBlockDiagonal(user->op_ijacobian, coo_values, CEED_REQUEST_IMMEDIATE));
25875d1798cSJames Wright     PetscCall(PetscLogGpuTimeEnd());
25975d1798cSJames Wright     PetscCall(PetscLogEventEnd(FLUIDS_CeedOperatorAssemblePointBlockDiagonal, J, 0, 0, 0));
260544be873SJed Brown   } else {
26175d1798cSJames Wright     PetscCall(PetscLogEventBegin(FLUIDS_CeedOperatorAssemble, J, 0, 0, 0));
26275d1798cSJames Wright     PetscCall(PetscLogGpuTimeBegin());
263a424bcd0SJames Wright     PetscCallCeed(user->ceed, CeedOperatorLinearAssemble(user->op_ijacobian, coo_values));
26475d1798cSJames Wright     PetscCall(PetscLogGpuTimeEnd());
26575d1798cSJames Wright     PetscCall(PetscLogEventEnd(FLUIDS_CeedOperatorAssemble, J, 0, 0, 0));
266544be873SJed Brown   }
267a424bcd0SJames Wright   PetscCallCeed(user->ceed, CeedVectorGetArrayRead(coo_values, mem_type, &values));
268544be873SJed Brown   PetscCall(MatSetValuesCOO(J, values, INSERT_VALUES));
269a424bcd0SJames Wright   PetscCallCeed(user->ceed, CeedVectorRestoreArrayRead(coo_values, &values));
270ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
271544be873SJed Brown }
272544be873SJed Brown 
2732b730f8bSJeremy L Thompson PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, PetscReal shift, Mat J, Mat J_pre, void *user_data) {
274e334ad8fSJed Brown   User      user = *(User *)user_data;
275a424bcd0SJames Wright   Ceed      ceed = user->ceed;
276d6bf345cSJed Brown   PetscBool J_is_shell, J_is_mffd, J_pre_is_shell;
277*f17d818dSJames Wright 
278e334ad8fSJed Brown   PetscFunctionBeginUser;
279a424bcd0SJames Wright   if (user->phys->ijacobian_time_shift_label)
280a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_ijacobian, user->phys->ijacobian_time_shift_label, &shift));
281d6bf345cSJed Brown   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATMFFD, &J_is_mffd));
282e334ad8fSJed Brown   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATSHELL, &J_is_shell));
2832b730f8bSJeremy L Thompson   PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATSHELL, &J_pre_is_shell));
284e334ad8fSJed Brown   if (!user->matrices_set_up) {
285e334ad8fSJed Brown     if (J_is_shell) {
286b07979f9SJames Wright       OperatorApplyContext op_ijacobian_ctx;
287b07979f9SJames Wright       OperatorApplyContextCreate(user->dm, user->dm, user->ceed, user->op_ijacobian, user->q_ceed, user->g_ceed, user->Q_dot_loc, NULL,
288b07979f9SJames Wright                                  &op_ijacobian_ctx);
289b07979f9SJames Wright       PetscCall(MatShellSetContext(J, op_ijacobian_ctx));
290b07979f9SJames Wright       PetscCall(MatShellSetContextDestroy(J, (PetscErrorCode(*)(void *))OperatorApplyContextDestroy));
291b07979f9SJames Wright       PetscCall(MatShellSetOperation(J, MATOP_MULT, (void (*)(void))MatMult_Ceed));
292b07979f9SJames Wright       PetscCall(MatShellSetOperation(J, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiag_Ceed));
293e334ad8fSJed Brown       PetscCall(MatSetUp(J));
294e334ad8fSJed Brown     }
295e334ad8fSJed Brown     if (!J_pre_is_shell) {
2962b730f8bSJeremy L Thompson       PetscCall(FormPreallocation(user, user->app_ctx->pmat_pbdiagonal, J_pre, &user->coo_values_pmat));
297544be873SJed Brown     }
298d6bf345cSJed Brown     if (J != J_pre && !J_is_shell && !J_is_mffd) {
299544be873SJed Brown       PetscCall(FormPreallocation(user, PETSC_FALSE, J, &user->coo_values_amat));
300544be873SJed Brown     }
301e334ad8fSJed Brown     user->matrices_set_up = true;
302e334ad8fSJed Brown   }
303e334ad8fSJed Brown   if (!J_pre_is_shell) {
3042b730f8bSJeremy L Thompson     PetscCall(FormSetValues(user, user->app_ctx->pmat_pbdiagonal, J_pre, user->coo_values_pmat));
305e334ad8fSJed Brown   }
306d6bf345cSJed Brown   if (user->coo_values_amat) {
307d6bf345cSJed Brown     PetscCall(FormSetValues(user, PETSC_FALSE, J, user->coo_values_amat));
308d6bf345cSJed Brown   } else if (J_is_mffd) {
309d6bf345cSJed Brown     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
310d6bf345cSJed Brown     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
311d6bf345cSJed Brown   }
312ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
313e334ad8fSJed Brown }
314e334ad8fSJed Brown 
3152b730f8bSJeremy L Thompson PetscErrorCode WriteOutput(User user, Vec Q, PetscInt step_no, PetscScalar time) {
31677841947SLeila Ghaffari   Vec         Q_loc;
31777841947SLeila Ghaffari   char        file_path[PETSC_MAX_PATH_LEN];
31877841947SLeila Ghaffari   PetscViewer viewer;
31977841947SLeila Ghaffari 
320*f17d818dSJames Wright   PetscFunctionBeginUser;
32137cbb16aSJed Brown   if (user->app_ctx->checkpoint_vtk) {
32277841947SLeila Ghaffari     // Set up output
323f14660a4SJames Wright     PetscCall(DMGetLocalVector(user->dm, &Q_loc));
324f14660a4SJames Wright     PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec"));
325f14660a4SJames Wright     PetscCall(VecZeroEntries(Q_loc));
326f14660a4SJames Wright     PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
32777841947SLeila Ghaffari 
32877841947SLeila Ghaffari     // Output
32937cbb16aSJed Brown     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
330f14660a4SJames Wright 
3312b730f8bSJeremy L Thompson     PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path, FILE_MODE_WRITE, &viewer));
332f14660a4SJames Wright     PetscCall(VecView(Q_loc, viewer));
333f14660a4SJames Wright     PetscCall(PetscViewerDestroy(&viewer));
33477841947SLeila Ghaffari     if (user->dm_viz) {
33577841947SLeila Ghaffari       Vec         Q_refined, Q_refined_loc;
33677841947SLeila Ghaffari       char        file_path_refined[PETSC_MAX_PATH_LEN];
33777841947SLeila Ghaffari       PetscViewer viewer_refined;
33877841947SLeila Ghaffari 
339f14660a4SJames Wright       PetscCall(DMGetGlobalVector(user->dm_viz, &Q_refined));
340f14660a4SJames Wright       PetscCall(DMGetLocalVector(user->dm_viz, &Q_refined_loc));
341f14660a4SJames Wright       PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined"));
342f14660a4SJames Wright 
343f14660a4SJames Wright       PetscCall(MatInterpolate(user->interp_viz, Q, Q_refined));
344f14660a4SJames Wright       PetscCall(VecZeroEntries(Q_refined_loc));
3452b730f8bSJeremy L Thompson       PetscCall(DMGlobalToLocal(user->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc));
346f14660a4SJames Wright 
34737cbb16aSJed Brown       PetscCall(
34837cbb16aSJed Brown           PetscSNPrintf(file_path_refined, sizeof file_path_refined, "%s/nsrefined-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
349f14660a4SJames Wright 
3502b730f8bSJeremy L Thompson       PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined), file_path_refined, FILE_MODE_WRITE, &viewer_refined));
351f14660a4SJames Wright       PetscCall(VecView(Q_refined_loc, viewer_refined));
352f14660a4SJames Wright       PetscCall(DMRestoreLocalVector(user->dm_viz, &Q_refined_loc));
353f14660a4SJames Wright       PetscCall(DMRestoreGlobalVector(user->dm_viz, &Q_refined));
354f14660a4SJames Wright       PetscCall(PetscViewerDestroy(&viewer_refined));
35577841947SLeila Ghaffari     }
356f14660a4SJames Wright     PetscCall(DMRestoreLocalVector(user->dm, &Q_loc));
35737cbb16aSJed Brown   }
35877841947SLeila Ghaffari 
35977841947SLeila Ghaffari   // Save data in a binary file for continuation of simulations
36069293791SJames Wright   if (user->app_ctx->add_stepnum2bin) {
36137cbb16aSJed Brown     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", user->app_ctx->output_dir, step_no));
36269293791SJames Wright   } else {
3632b730f8bSJeremy L Thompson     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", user->app_ctx->output_dir));
36469293791SJames Wright   }
3652b730f8bSJeremy L Thompson   PetscCall(PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, &viewer));
366f14660a4SJames Wright 
3670de6a49fSJames Wright   PetscInt32 token = PetscDefined(USE_64BIT_INDICES) ? FLUIDS_FILE_TOKEN_64 : FLUIDS_FILE_TOKEN_32;
3680de6a49fSJames Wright   PetscCall(PetscViewerBinaryWrite(viewer, &token, 1, PETSC_INT32));
3694de8550aSJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &step_no, 1, PETSC_INT));
3704de8550aSJed Brown   time /= user->units->second;  // Dimensionalize time back
3714de8550aSJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL));
372f14660a4SJames Wright   PetscCall(VecView(Q, viewer));
373f14660a4SJames Wright   PetscCall(PetscViewerDestroy(&viewer));
374ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
375f14660a4SJames Wright }
376f14660a4SJames Wright 
377ca69d878SAdeleke O. Bankole // CSV Monitor
378ca69d878SAdeleke O. Bankole PetscErrorCode TSMonitor_WallForce(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
379ca69d878SAdeleke O. Bankole   User              user = ctx;
380ca69d878SAdeleke O. Bankole   Vec               G_loc;
381ca69d878SAdeleke O. Bankole   PetscInt          num_wall = user->app_ctx->wall_forces.num_wall, dim = 3;
382ca69d878SAdeleke O. Bankole   const PetscInt   *walls  = user->app_ctx->wall_forces.walls;
383ca69d878SAdeleke O. Bankole   PetscViewer       viewer = user->app_ctx->wall_forces.viewer;
384ca69d878SAdeleke O. Bankole   PetscViewerFormat format = user->app_ctx->wall_forces.viewer_format;
385ca69d878SAdeleke O. Bankole   PetscScalar      *reaction_force;
386ca69d878SAdeleke O. Bankole   PetscBool         iascii;
387ca69d878SAdeleke O. Bankole 
388ca69d878SAdeleke O. Bankole   PetscFunctionBeginUser;
389ee4ca9cbSJames Wright   if (!viewer) PetscFunctionReturn(PETSC_SUCCESS);
390ca69d878SAdeleke O. Bankole   PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
391ca69d878SAdeleke O. Bankole   PetscCall(PetscMalloc1(num_wall * dim, &reaction_force));
392ca69d878SAdeleke O. Bankole   PetscCall(Surface_Forces_NS(user->dm, G_loc, num_wall, walls, reaction_force));
393ca69d878SAdeleke O. Bankole   PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
394ca69d878SAdeleke O. Bankole 
395ca69d878SAdeleke O. Bankole   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
396ca69d878SAdeleke O. Bankole 
397ca69d878SAdeleke O. Bankole   if (iascii) {
398ca69d878SAdeleke O. Bankole     if (format == PETSC_VIEWER_ASCII_CSV && !user->app_ctx->wall_forces.header_written) {
399ca69d878SAdeleke O. Bankole       PetscCall(PetscViewerASCIIPrintf(viewer, "Step,Time,Wall,ForceX,ForceY,ForceZ\n"));
400ca69d878SAdeleke O. Bankole       user->app_ctx->wall_forces.header_written = PETSC_TRUE;
401ca69d878SAdeleke O. Bankole     }
402ca69d878SAdeleke O. Bankole     for (PetscInt w = 0; w < num_wall; w++) {
403ca69d878SAdeleke O. Bankole       PetscInt wall = walls[w];
404ca69d878SAdeleke O. Bankole       if (format == PETSC_VIEWER_ASCII_CSV) {
405ca69d878SAdeleke O. Bankole         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT ",%g,%" PetscInt_FMT ",%g,%g,%g\n", step_no, time, wall,
406ca69d878SAdeleke O. Bankole                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
407ca69d878SAdeleke O. Bankole 
408ca69d878SAdeleke O. Bankole       } else {
409ca69d878SAdeleke O. Bankole         PetscCall(PetscViewerASCIIPrintf(viewer, "Wall %" PetscInt_FMT " Forces: Force_x = %12g, Force_y = %12g, Force_z = %12g\n", wall,
410ca69d878SAdeleke O. Bankole                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
411ca69d878SAdeleke O. Bankole       }
412ca69d878SAdeleke O. Bankole     }
413ca69d878SAdeleke O. Bankole   }
414ca69d878SAdeleke O. Bankole   PetscCall(PetscFree(reaction_force));
415ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
416ca69d878SAdeleke O. Bankole }
417ca69d878SAdeleke O. Bankole 
418f14660a4SJames Wright // User provided TS Monitor
4192b730f8bSJeremy L Thompson PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
420f14660a4SJames Wright   User user = ctx;
421f14660a4SJames Wright 
422*f17d818dSJames Wright   PetscFunctionBeginUser;
42337cbb16aSJed Brown   // Print every 'checkpoint_interval' steps
424894de27cSJames Wright   if (user->app_ctx->checkpoint_interval <= 0 || step_no % user->app_ctx->checkpoint_interval != 0 ||
42549aac155SJeremy L Thompson       (user->app_ctx->cont_steps == step_no && step_no != 0)) {
426ee4ca9cbSJames Wright     PetscFunctionReturn(PETSC_SUCCESS);
42749aac155SJeremy L Thompson   }
428f14660a4SJames Wright 
429f14660a4SJames Wright   PetscCall(WriteOutput(user, Q, step_no, time));
430ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
43177841947SLeila Ghaffari }
43277841947SLeila Ghaffari 
43377841947SLeila Ghaffari // TS: Create, setup, and solve
4342b730f8bSJeremy L Thompson PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, Vec *Q, PetscScalar *f_time, TS *ts) {
43577841947SLeila Ghaffari   MPI_Comm    comm = user->comm;
43677841947SLeila Ghaffari   TSAdapt     adapt;
43777841947SLeila Ghaffari   PetscScalar final_time;
43877841947SLeila Ghaffari 
439*f17d818dSJames Wright   PetscFunctionBeginUser;
4402b730f8bSJeremy L Thompson   PetscCall(TSCreate(comm, ts));
4412b730f8bSJeremy L Thompson   PetscCall(TSSetDM(*ts, dm));
44277841947SLeila Ghaffari   if (phys->implicit) {
4432b730f8bSJeremy L Thompson     PetscCall(TSSetType(*ts, TSBDF));
44477841947SLeila Ghaffari     if (user->op_ifunction) {
4452b730f8bSJeremy L Thompson       PetscCall(TSSetIFunction(*ts, NULL, IFunction_NS, &user));
44677841947SLeila Ghaffari     } else {  // Implicit integrators can fall back to using an RHSFunction
4472b730f8bSJeremy L Thompson       PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
44877841947SLeila Ghaffari     }
449e334ad8fSJed Brown     if (user->op_ijacobian) {
4502b730f8bSJeremy L Thompson       PetscCall(DMTSSetIJacobian(dm, FormIJacobian_NS, &user));
451544be873SJed Brown       if (app_ctx->amat_type) {
452544be873SJed Brown         Mat Pmat, Amat;
4532b730f8bSJeremy L Thompson         PetscCall(DMCreateMatrix(dm, &Pmat));
4542b730f8bSJeremy L Thompson         PetscCall(DMSetMatType(dm, app_ctx->amat_type));
4552b730f8bSJeremy L Thompson         PetscCall(DMCreateMatrix(dm, &Amat));
4562b730f8bSJeremy L Thompson         PetscCall(TSSetIJacobian(*ts, Amat, Pmat, NULL, NULL));
4572b730f8bSJeremy L Thompson         PetscCall(MatDestroy(&Amat));
4582b730f8bSJeremy L Thompson         PetscCall(MatDestroy(&Pmat));
459544be873SJed Brown       }
460e334ad8fSJed Brown     }
46177841947SLeila Ghaffari   } else {
4629ad5e8e4SJames Wright     PetscCheck(user->op_rhs_ctx, comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction");
4632b730f8bSJeremy L Thompson     PetscCall(TSSetType(*ts, TSRK));
4642b730f8bSJeremy L Thompson     PetscCall(TSRKSetType(*ts, TSRK5F));
4652b730f8bSJeremy L Thompson     PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
46677841947SLeila Ghaffari   }
4672b730f8bSJeremy L Thompson   PetscCall(TSSetMaxTime(*ts, 500. * user->units->second));
4682b730f8bSJeremy L Thompson   PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER));
46951b00d91SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(TSSetErrorIfStepFails(*ts, PETSC_FALSE));
4702b730f8bSJeremy L Thompson   PetscCall(TSSetTimeStep(*ts, 1.e-2 * user->units->second));
4718fb33541SJames Wright   if (app_ctx->test_type != TESTTYPE_NONE) {
4722b730f8bSJeremy L Thompson     PetscCall(TSSetMaxSteps(*ts, 10));
4732b730f8bSJeremy L Thompson   }
4742b730f8bSJeremy L Thompson   PetscCall(TSGetAdapt(*ts, &adapt));
4752b730f8bSJeremy L Thompson   PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * user->units->second, 1.e2 * user->units->second));
4762b730f8bSJeremy L Thompson   PetscCall(TSSetFromOptions(*ts));
47715637395SJames Wright   user->time_bc_set = -1.0;   // require all BCs be updated
478d55646a4SJames Wright   if (app_ctx->cont_steps) {  // continue from previous timestep data
47977841947SLeila Ghaffari     PetscInt    count;
48077841947SLeila Ghaffari     PetscViewer viewer;
4812b730f8bSJeremy L Thompson 
4824de8550aSJed Brown     if (app_ctx->cont_time <= 0) {  // Legacy files did not include step number and time
4832b730f8bSJeremy L Thompson       PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_time_file, FILE_MODE_READ, &viewer));
4844de8550aSJed Brown       PetscCall(PetscViewerBinaryRead(viewer, &app_ctx->cont_time, 1, &count, PETSC_REAL));
4852b730f8bSJeremy L Thompson       PetscCall(PetscViewerDestroy(&viewer));
4864de8550aSJed Brown       PetscCheck(app_ctx->cont_steps != -1, comm, PETSC_ERR_ARG_INCOMP,
4874de8550aSJed Brown                  "-continue step number not specified, but checkpoint file does not contain a step number (likely written by older code version)");
4884de8550aSJed Brown     }
4894de8550aSJed Brown     PetscCall(TSSetTime(*ts, app_ctx->cont_time * user->units->second));
490d8e0aecdSJed Brown     PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps));
49177841947SLeila Ghaffari   }
4928fb33541SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) {
4932b730f8bSJeremy L Thompson     PetscCall(TSMonitorSet(*ts, TSMonitor_NS, user, NULL));
4948fb33541SJames Wright   }
495ca69d878SAdeleke O. Bankole   if (app_ctx->wall_forces.viewer) {
496ca69d878SAdeleke O. Bankole     PetscCall(TSMonitorSet(*ts, TSMonitor_WallForce, user, NULL));
497ca69d878SAdeleke O. Bankole   }
498b7d66439SJames Wright   if (app_ctx->turb_spanstats_enable) {
499f5452247SJames Wright     PetscCall(TSMonitorSet(*ts, TSMonitor_TurbulenceStatistics, user, NULL));
500495b9769SJames Wright     CeedScalar previous_time = app_ctx->cont_time * user->units->second;
501a424bcd0SJames Wright     PetscCallCeed(user->ceed,
502a424bcd0SJames Wright                   CeedOperatorSetContextDouble(user->spanstats.op_stats_collect_ctx->op, user->spanstats.previous_time_label, &previous_time));
503a175e481SJames Wright   }
5044e9802d1SJames Wright   if (app_ctx->diff_filter_monitor) PetscCall(TSMonitorSet(*ts, TSMonitor_DifferentialFilter, user, NULL));
50577841947SLeila Ghaffari 
50677841947SLeila Ghaffari   // Solve
507d8e0aecdSJed Brown   PetscReal start_time;
508d8e0aecdSJed Brown   PetscInt  start_step;
5092b730f8bSJeremy L Thompson   PetscCall(TSGetTime(*ts, &start_time));
510d8e0aecdSJed Brown   PetscCall(TSGetStepNumber(*ts, &start_step));
5117b39487dSJeremy L Thompson 
512711a423aSJed Brown   PetscCall(PetscLogDefaultBegin());  // So we can use PetscLogStageGetPerfInfo without -log_view
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));
526fc4c3d69SJames Wright     PetscCall(TSSetSolution(*ts, Q_preload));
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;
548711a423aSJed Brown     PetscEventPerfInfo stage_perf;
5497b39487dSJeremy L Thompson 
5507b39487dSJeremy L Thompson     PetscCall(PetscLogStageGetId("Fluids Solve", &stage_id));
551711a423aSJed Brown     PetscCall(PetscLogStageGetPerfInfo(stage_id, &stage_perf));
552711a423aSJed Brown     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time taken for solution (sec): %g\n", stage_perf.time));
55377841947SLeila Ghaffari   }
554ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
55577841947SLeila Ghaffari }
556