xref: /honee/src/setupts.c (revision 4e1dcb991ec2adc889afc518394faf42beafdaca)
1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3a515125bSLeila Ghaffari 
4a515125bSLeila Ghaffari /// @file
5ea615d4cSJames Wright /// Time-stepping functions for HONEE
6a515125bSLeila Ghaffari 
7e419654dSJeremy L Thompson #include <ceed.h>
8e419654dSJeremy L Thompson #include <petscdmplex.h>
9e419654dSJeremy L Thompson #include <petscts.h>
10e419654dSJeremy L Thompson 
11149fb536SJames Wright #include <navierstokes.h>
12c5e9980aSAdeleke O. Bankole #include "../qfunctions/newtonian_state.h"
13a515125bSLeila Ghaffari 
144cb09fddSJames Wright // @brief Insert Boundary values if it's a new time
150c373b74SJames Wright PetscErrorCode UpdateBoundaryValues(Honee honee, Vec Q_loc, PetscReal t) {
16c996854bSJames Wright   PetscFunctionBeginUser;
170c373b74SJames Wright   if (honee->time_bc_set != t) {
180c373b74SJames Wright     PetscCall(DMPlexInsertBoundaryValues(honee->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL));
190c373b74SJames Wright     honee->time_bc_set = t;
20c996854bSJames Wright   }
21d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
22c996854bSJames Wright }
23c996854bSJames Wright 
24a515125bSLeila Ghaffari // RHS (Explicit time-stepper) function setup
25a515125bSLeila Ghaffari //   This is the RHS of the ODE, given as u_t = G(t,u)
26a515125bSLeila Ghaffari //   This function takes in a state vector Q and writes into G
27a515125bSLeila Ghaffari PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) {
280c373b74SJames Wright   Honee        honee = *(Honee *)user_data;
290c373b74SJames Wright   Ceed         ceed  = honee->ceed;
30fd969b44SJames Wright   PetscScalar  dt;
31dc3c760aSJames Wright   Vec          Q_loc = honee->Q_loc, R;
32a78efa86SJames Wright   PetscMemType q_mem_type;
33a515125bSLeila Ghaffari 
3406f41313SJames Wright   PetscFunctionBeginUser;
35e2f84137SJeremy L Thompson   // Update time dependent data
360c373b74SJames Wright   PetscCall(UpdateBoundaryValues(honee, Q_loc, t));
370c373b74SJames Wright   if (honee->phys->solution_time_label)
380c373b74SJames Wright     PetscCallCeed(ceed, CeedOperatorSetContextDouble(honee->op_rhs_ctx->op, honee->phys->solution_time_label, &t));
392b916ea7SJeremy L Thompson   PetscCall(TSGetTimeStep(ts, &dt));
400c373b74SJames Wright   if (honee->phys->timestep_size_label)
410c373b74SJames Wright     PetscCallCeed(ceed, CeedOperatorSetContextDouble(honee->op_rhs_ctx->op, honee->phys->timestep_size_label, &dt));
42a515125bSLeila Ghaffari 
43dc3c760aSJames Wright   PetscCall(DMGetNamedGlobalVector(honee->dm, "RHS Residual", &R));
440c373b74SJames Wright   PetscCall(DMGlobalToLocal(honee->dm, Q, INSERT_VALUES, Q_loc));
450c373b74SJames Wright   if (honee->app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) PetscCall(DivDiffFluxProjectionApply(honee->diff_flux_proj, Q_loc));
46dc3c760aSJames Wright   PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, R, honee->op_rhs_ctx));
47a78efa86SJames Wright 
48a78efa86SJames Wright   // Inverse of the mass matrix
49dc3c760aSJames Wright   PetscCall(VecReadPetscToCeed(Q_loc, &q_mem_type, honee->q_ceed));
50bc5573aaSJames Wright   {
51bc5573aaSJames Wright     // Run PCApply manually if using ksp_type preonly -pc_type jacobi
52bc5573aaSJames Wright     // This is to avoid an AllReduce call in KSPSolve_Preonly, which causes significant slowdowns for lumped mass matrix solves.
53bc5573aaSJames Wright     // See https://gitlab.com/petsc/petsc/-/merge_requests/8048 for more details and a possible fix
54bc5573aaSJames Wright     PC        pc;
55bc5573aaSJames Wright     PetscBool ispreonly, isjacobi;
56bc5573aaSJames Wright     PetscCall(KSPGetPC(honee->mass_ksp, &pc));
57bc5573aaSJames Wright     PetscCall(PetscObjectTypeCompare((PetscObject)honee->mass_ksp, KSPPREONLY, &ispreonly));
58bc5573aaSJames Wright     PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &isjacobi));
59bc5573aaSJames Wright     if (ispreonly && isjacobi) PetscCall(PCApply(pc, R, G));
60bc5573aaSJames Wright     else PetscCall(KSPSolve(honee->mass_ksp, R, G));
61bc5573aaSJames Wright   }
620c373b74SJames Wright   PetscCall(VecReadCeedToPetsc(honee->q_ceed, q_mem_type, Q_loc));
63dc3c760aSJames Wright 
64dc3c760aSJames Wright   PetscCall(DMRestoreNamedGlobalVector(honee->dm, "RHS Residual", &R));
65d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
66a515125bSLeila Ghaffari }
67a515125bSLeila Ghaffari 
68c5e9980aSAdeleke O. Bankole // Surface forces function setup
69c5e9980aSAdeleke O. Bankole static PetscErrorCode Surface_Forces_NS(DM dm, Vec G_loc, PetscInt num_walls, const PetscInt walls[], PetscScalar *reaction_force) {
70c5e9980aSAdeleke O. Bankole   DMLabel            face_label;
7162ab3c0bSJames Wright   const PetscScalar *g_array;
7262ab3c0bSJames Wright   PetscInt           dim  = 3;
7362ab3c0bSJames Wright   MPI_Comm           comm = PetscObjectComm((PetscObject)dm);
7462ab3c0bSJames Wright   PetscSection       section;
75c5e9980aSAdeleke O. Bankole 
76c5e9980aSAdeleke O. Bankole   PetscFunctionBeginUser;
77c5e9980aSAdeleke O. Bankole   PetscCall(DMGetLabel(dm, "Face Sets", &face_label));
7862ab3c0bSJames Wright   PetscCall(VecGetArrayRead(G_loc, &g_array));
79c5e9980aSAdeleke O. Bankole   for (PetscInt w = 0; w < num_walls; w++) {
8062ab3c0bSJames Wright     const PetscInt wall = walls[w], *points;
81c5e9980aSAdeleke O. Bankole     IS             wall_is;
8262ab3c0bSJames Wright     PetscInt       num_points, num_comp = 0;
8362ab3c0bSJames Wright 
84c5e9980aSAdeleke O. Bankole     PetscCall(DMLabelGetStratumIS(face_label, wall, &wall_is));
8562ab3c0bSJames Wright     if (!wall_is) continue;  // No wall points on this process, skip
8662ab3c0bSJames Wright 
8762ab3c0bSJames Wright     PetscCall(DMGetLocalSection(dm, &section));
8862ab3c0bSJames Wright     PetscCall(PetscSectionGetFieldComponents(section, 0, &num_comp));
89c5e9980aSAdeleke O. Bankole     PetscCall(ISGetSize(wall_is, &num_points));
90c5e9980aSAdeleke O. Bankole     PetscCall(ISGetIndices(wall_is, &points));
91c5e9980aSAdeleke O. Bankole     for (PetscInt i = 0; i < num_points; i++) {
92c5e9980aSAdeleke O. Bankole       const PetscInt           p = points[i];
93c5e9980aSAdeleke O. Bankole       const StateConservative *r;
9462ab3c0bSJames Wright       PetscInt                 dof;
9562ab3c0bSJames Wright 
9662ab3c0bSJames Wright       PetscCall(DMPlexPointLocalRead(dm, p, g_array, &r));
9762ab3c0bSJames Wright       PetscCall(PetscSectionGetDof(section, p, &dof));
982004e3acSAdeleke O. Bankole       for (PetscInt node = 0; node < dof / num_comp; node++) {
9962ab3c0bSJames Wright         for (PetscInt j = 0; j < dim; j++) {
1002004e3acSAdeleke O. Bankole           reaction_force[w * dim + j] -= r[node].momentum[j];
1012004e3acSAdeleke O. Bankole         }
102c5e9980aSAdeleke O. Bankole       }
103c5e9980aSAdeleke O. Bankole     }
104c5e9980aSAdeleke O. Bankole     PetscCall(ISRestoreIndices(wall_is, &points));
105c5e9980aSAdeleke O. Bankole     PetscCall(ISDestroy(&wall_is));
106c5e9980aSAdeleke O. Bankole   }
107c5e9980aSAdeleke O. Bankole   PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, reaction_force, dim * num_walls, MPIU_SCALAR, MPI_SUM, comm));
10862ab3c0bSJames Wright   PetscCall(VecRestoreArrayRead(G_loc, &g_array));
109d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
110c5e9980aSAdeleke O. Bankole }
111c5e9980aSAdeleke O. Bankole 
112a515125bSLeila Ghaffari // Implicit time-stepper function setup
1132b916ea7SJeremy L Thompson PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) {
1140c373b74SJames Wright   Honee        honee = *(Honee *)user_data;
1150c373b74SJames Wright   Ceed         ceed  = honee->ceed;
116fd969b44SJames Wright   PetscScalar  dt;
1170c373b74SJames Wright   Vec          Q_loc = honee->Q_loc, Q_dot_loc = honee->Q_dot_loc, G_loc;
118a515125bSLeila Ghaffari   PetscMemType q_mem_type, q_dot_mem_type, g_mem_type;
119a515125bSLeila Ghaffari 
12006f41313SJames Wright   PetscFunctionBeginUser;
121dabd2275SJames Wright   PetscCall(DMGlobalToLocalBegin(honee->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
1220c373b74SJames Wright   PetscCall(DMGetNamedLocalVector(honee->dm, "ResidualLocal", &G_loc));
123e2f84137SJeremy L Thompson 
124e2f84137SJeremy L Thompson   // Update time dependent data
1250c373b74SJames Wright   PetscCall(UpdateBoundaryValues(honee, Q_loc, t));
1260c373b74SJames Wright   if (honee->phys->solution_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(honee->op_ifunction, honee->phys->solution_time_label, &t));
1272b916ea7SJeremy L Thompson   PetscCall(TSGetTimeStep(ts, &dt));
1280c373b74SJames Wright   if (honee->phys->timestep_size_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(honee->op_ifunction, honee->phys->timestep_size_label, &dt));
129a515125bSLeila Ghaffari 
130a515125bSLeila Ghaffari   // Global-to-local
1310c373b74SJames Wright   PetscCall(DMGlobalToLocalBegin(honee->dm, Q, INSERT_VALUES, Q_loc));
1320c373b74SJames Wright   PetscCall(DMGlobalToLocalEnd(honee->dm, Q, INSERT_VALUES, Q_loc));
1330c373b74SJames Wright   if (honee->app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) PetscCall(DivDiffFluxProjectionApply(honee->diff_flux_proj, Q_loc));
1340c373b74SJames Wright   PetscCall(DMGlobalToLocalEnd(honee->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
135a515125bSLeila Ghaffari 
136a515125bSLeila Ghaffari   // Place PETSc vectors in CEED vectors
1370c373b74SJames Wright   PetscCall(VecReadPetscToCeed(Q_loc, &q_mem_type, honee->q_ceed));
1380c373b74SJames Wright   PetscCall(VecReadPetscToCeed(Q_dot_loc, &q_dot_mem_type, honee->q_dot_ceed));
1390c373b74SJames Wright   PetscCall(VecPetscToCeed(G_loc, &g_mem_type, honee->g_ceed));
140a515125bSLeila Ghaffari 
141a515125bSLeila Ghaffari   // Apply CEED operator
142ea615d4cSJames Wright   PetscCall(PetscLogEventBegin(HONEE_CeedOperatorApply, Q, G, 0, 0));
1437eedc94cSJames Wright   PetscCall(PetscLogGpuTimeBegin());
1440c373b74SJames Wright   PetscCallCeed(honee->ceed, CeedOperatorApply(honee->op_ifunction, honee->q_ceed, honee->g_ceed, CEED_REQUEST_IMMEDIATE));
1457eedc94cSJames Wright   PetscCall(PetscLogGpuTimeEnd());
146ea615d4cSJames Wright   PetscCall(PetscLogEventEnd(HONEE_CeedOperatorApply, Q, G, 0, 0));
147a515125bSLeila Ghaffari 
148a515125bSLeila Ghaffari   // Restore vectors
1490c373b74SJames Wright   PetscCall(VecReadCeedToPetsc(honee->q_ceed, q_mem_type, Q_loc));
1500c373b74SJames Wright   PetscCall(VecReadCeedToPetsc(honee->q_dot_ceed, q_dot_mem_type, Q_dot_loc));
1510c373b74SJames Wright   PetscCall(VecCeedToPetsc(honee->g_ceed, g_mem_type, G_loc));
152a515125bSLeila Ghaffari 
1530c373b74SJames Wright   if (honee->app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) {
1540c373b74SJames Wright     PetscCall(SgsDDApplyIFunction(honee, Q_loc, G_loc));
15501ab89c1SJames Wright   }
1569c678832SJames Wright 
157a515125bSLeila Ghaffari   // Local-to-Global
1582b916ea7SJeremy L Thompson   PetscCall(VecZeroEntries(G));
1590c373b74SJames Wright   PetscCall(DMLocalToGlobal(honee->dm, G_loc, ADD_VALUES, G));
160a515125bSLeila Ghaffari 
161a515125bSLeila Ghaffari   // Restore vectors
1620c373b74SJames Wright   PetscCall(DMRestoreNamedLocalVector(honee->dm, "ResidualLocal", &G_loc));
163d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
164a515125bSLeila Ghaffari }
165a515125bSLeila Ghaffari 
1662b916ea7SJeremy L Thompson PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, PetscReal shift, Mat J, Mat J_pre, void *user_data) {
1670c373b74SJames Wright   Honee     honee = *(Honee *)user_data;
168ebfabadfSJames Wright   PetscBool J_is_matceed, J_is_mffd, J_pre_is_matceed, J_pre_is_mffd;
16906f41313SJames Wright 
170f0b65372SJed Brown   PetscFunctionBeginUser;
17104855949SJed Brown   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATMFFD, &J_is_mffd));
172ebfabadfSJames Wright   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATCEED, &J_is_matceed));
173ebfabadfSJames Wright   PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATMFFD, &J_pre_is_mffd));
174ebfabadfSJames Wright   PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATCEED, &J_pre_is_matceed));
175ebfabadfSJames Wright 
1760c373b74SJames Wright   PetscCall(MatCeedSetContextReal(honee->mat_ijacobian, "ijacobian time shift", shift));
177ebfabadfSJames Wright 
178ebfabadfSJames Wright   if (J_is_matceed || J_is_mffd) {
17904855949SJed Brown     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
18004855949SJed Brown     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
1810c373b74SJames Wright   } else PetscCall(MatCeedAssembleCOO(honee->mat_ijacobian, J));
182ebfabadfSJames Wright 
183ebfabadfSJames Wright   if (J_pre_is_matceed && J != J_pre) {
184ebfabadfSJames Wright     PetscCall(MatAssemblyBegin(J_pre, MAT_FINAL_ASSEMBLY));
185ebfabadfSJames Wright     PetscCall(MatAssemblyEnd(J_pre, MAT_FINAL_ASSEMBLY));
186ebfabadfSJames Wright   } else if (!J_pre_is_matceed && !J_pre_is_mffd && J != J_pre) {
1870c373b74SJames Wright     PetscCall(MatCeedAssembleCOO(honee->mat_ijacobian, J_pre));
18804855949SJed Brown   }
189d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
190f0b65372SJed Brown }
191f0b65372SJed Brown 
1920c373b74SJames Wright PetscErrorCode WriteOutput(Honee honee, Vec Q, PetscInt step_no, PetscScalar time) {
193a515125bSLeila Ghaffari   Vec         Q_loc;
194a515125bSLeila Ghaffari   char        file_path[PETSC_MAX_PATH_LEN];
195a515125bSLeila Ghaffari   PetscViewer viewer;
196a515125bSLeila Ghaffari 
19706f41313SJames Wright   PetscFunctionBeginUser;
1980c373b74SJames Wright   if (honee->app_ctx->checkpoint_vtk) {
199a515125bSLeila Ghaffari     // Set up output
2000c373b74SJames Wright     PetscCall(DMGetLocalVector(honee->dm, &Q_loc));
2017538d537SJames Wright     PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec"));
2027538d537SJames Wright     PetscCall(VecZeroEntries(Q_loc));
2030c373b74SJames Wright     PetscCall(DMGlobalToLocal(honee->dm, Q, INSERT_VALUES, Q_loc));
204a515125bSLeila Ghaffari 
205a515125bSLeila Ghaffari     // Output
2060c373b74SJames Wright     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-%03" PetscInt_FMT ".vtu", honee->app_ctx->output_dir, step_no));
2077538d537SJames Wright 
2082b916ea7SJeremy L Thompson     PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path, FILE_MODE_WRITE, &viewer));
2097538d537SJames Wright     PetscCall(VecView(Q_loc, viewer));
2107538d537SJames Wright     PetscCall(PetscViewerDestroy(&viewer));
2110c373b74SJames Wright     if (honee->dm_viz) {
212a515125bSLeila Ghaffari       Vec         Q_refined, Q_refined_loc;
213a515125bSLeila Ghaffari       char        file_path_refined[PETSC_MAX_PATH_LEN];
214a515125bSLeila Ghaffari       PetscViewer viewer_refined;
215a515125bSLeila Ghaffari 
2160c373b74SJames Wright       PetscCall(DMGetGlobalVector(honee->dm_viz, &Q_refined));
2170c373b74SJames Wright       PetscCall(DMGetLocalVector(honee->dm_viz, &Q_refined_loc));
2187538d537SJames Wright       PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined"));
2197538d537SJames Wright 
2200c373b74SJames Wright       PetscCall(MatInterpolate(honee->interp_viz, Q, Q_refined));
2217538d537SJames Wright       PetscCall(VecZeroEntries(Q_refined_loc));
2220c373b74SJames Wright       PetscCall(DMGlobalToLocal(honee->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc));
2237538d537SJames Wright 
224852e5969SJed Brown       PetscCall(
2250c373b74SJames Wright           PetscSNPrintf(file_path_refined, sizeof file_path_refined, "%s/nsrefined-%03" PetscInt_FMT ".vtu", honee->app_ctx->output_dir, step_no));
2267538d537SJames Wright 
2272b916ea7SJeremy L Thompson       PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined), file_path_refined, FILE_MODE_WRITE, &viewer_refined));
2287538d537SJames Wright       PetscCall(VecView(Q_refined_loc, viewer_refined));
2290c373b74SJames Wright       PetscCall(DMRestoreLocalVector(honee->dm_viz, &Q_refined_loc));
2300c373b74SJames Wright       PetscCall(DMRestoreGlobalVector(honee->dm_viz, &Q_refined));
2317538d537SJames Wright       PetscCall(PetscViewerDestroy(&viewer_refined));
232a515125bSLeila Ghaffari     }
2330c373b74SJames Wright     PetscCall(DMRestoreLocalVector(honee->dm, &Q_loc));
234852e5969SJed Brown   }
235a515125bSLeila Ghaffari 
236a515125bSLeila Ghaffari   // Save data in a binary file for continuation of simulations
2370c373b74SJames Wright   if (honee->app_ctx->add_stepnum2bin) {
2380c373b74SJames Wright     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", honee->app_ctx->output_dir, step_no));
23991a36801SJames Wright   } else {
2400c373b74SJames Wright     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", honee->app_ctx->output_dir));
24191a36801SJames Wright   }
2420c373b74SJames Wright   PetscCall(PetscViewerBinaryOpen(honee->comm, file_path, FILE_MODE_WRITE, &viewer));
2437538d537SJames Wright 
2440c373b74SJames Wright   time /= honee->units->second;  // Dimensionalize time back
245b237916aSJames Wright   PetscCall(HoneeWriteBinaryVec(viewer, Q, time, step_no));
2467538d537SJames Wright   PetscCall(PetscViewerDestroy(&viewer));
247d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2487538d537SJames Wright }
2497538d537SJames Wright 
250c5e9980aSAdeleke O. Bankole // CSV Monitor
251c5e9980aSAdeleke O. Bankole PetscErrorCode TSMonitor_WallForce(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
2520c373b74SJames Wright   Honee             honee = ctx;
253c5e9980aSAdeleke O. Bankole   Vec               G_loc;
2540c373b74SJames Wright   PetscInt          num_wall = honee->app_ctx->wall_forces.num_wall, dim = 3;
2550c373b74SJames Wright   const PetscInt   *walls  = honee->app_ctx->wall_forces.walls;
2560c373b74SJames Wright   PetscViewer       viewer = honee->app_ctx->wall_forces.viewer;
2570c373b74SJames Wright   PetscViewerFormat format = honee->app_ctx->wall_forces.viewer_format;
258c5e9980aSAdeleke O. Bankole   PetscScalar      *reaction_force;
25962ab3c0bSJames Wright   PetscBool         is_ascii;
260c5e9980aSAdeleke O. Bankole 
261c5e9980aSAdeleke O. Bankole   PetscFunctionBeginUser;
262d949ddfcSJames Wright   if (!viewer) PetscFunctionReturn(PETSC_SUCCESS);
2630c373b74SJames Wright   PetscCall(DMGetNamedLocalVector(honee->dm, "ResidualLocal", &G_loc));
26462ab3c0bSJames Wright   PetscCall(PetscCalloc1(num_wall * dim, &reaction_force));
2650c373b74SJames Wright   PetscCall(Surface_Forces_NS(honee->dm, G_loc, num_wall, walls, reaction_force));
2660c373b74SJames Wright   PetscCall(DMRestoreNamedLocalVector(honee->dm, "ResidualLocal", &G_loc));
267c5e9980aSAdeleke O. Bankole 
26862ab3c0bSJames Wright   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &is_ascii));
269c5e9980aSAdeleke O. Bankole 
27062ab3c0bSJames Wright   if (is_ascii) {
2710c373b74SJames Wright     if (format == PETSC_VIEWER_ASCII_CSV && !honee->app_ctx->wall_forces.header_written) {
272c5e9980aSAdeleke O. Bankole       PetscCall(PetscViewerASCIIPrintf(viewer, "Step,Time,Wall,ForceX,ForceY,ForceZ\n"));
2730c373b74SJames Wright       honee->app_ctx->wall_forces.header_written = PETSC_TRUE;
274c5e9980aSAdeleke O. Bankole     }
275c5e9980aSAdeleke O. Bankole     for (PetscInt w = 0; w < num_wall; w++) {
276c5e9980aSAdeleke O. Bankole       PetscInt wall = walls[w];
277c5e9980aSAdeleke O. Bankole       if (format == PETSC_VIEWER_ASCII_CSV) {
278c5e9980aSAdeleke O. Bankole         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT ",%g,%" PetscInt_FMT ",%g,%g,%g\n", step_no, time, wall,
279c5e9980aSAdeleke O. Bankole                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
280c5e9980aSAdeleke O. Bankole 
281c5e9980aSAdeleke O. Bankole       } else {
282c5e9980aSAdeleke O. Bankole         PetscCall(PetscViewerASCIIPrintf(viewer, "Wall %" PetscInt_FMT " Forces: Force_x = %12g, Force_y = %12g, Force_z = %12g\n", wall,
283c5e9980aSAdeleke O. Bankole                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
284c5e9980aSAdeleke O. Bankole       }
285c5e9980aSAdeleke O. Bankole     }
286c5e9980aSAdeleke O. Bankole   }
287c5e9980aSAdeleke O. Bankole   PetscCall(PetscFree(reaction_force));
288d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
289c5e9980aSAdeleke O. Bankole }
290c5e9980aSAdeleke O. Bankole 
2917538d537SJames Wright // User provided TS Monitor
2922b916ea7SJeremy L Thompson PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
2930c373b74SJames Wright   Honee honee = ctx;
2947538d537SJames Wright 
29506f41313SJames Wright   PetscFunctionBeginUser;
296852e5969SJed Brown   // Print every 'checkpoint_interval' steps
2970c373b74SJames Wright   if (honee->app_ctx->checkpoint_interval <= 0 || step_no % honee->app_ctx->checkpoint_interval != 0 ||
2980c373b74SJames Wright       (honee->app_ctx->cont_steps == step_no && step_no != 0)) {
299d949ddfcSJames Wright     PetscFunctionReturn(PETSC_SUCCESS);
300e419654dSJeremy L Thompson   }
3017538d537SJames Wright 
3020c373b74SJames Wright   PetscCall(WriteOutput(honee, Q, step_no, time));
303d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
304a515125bSLeila Ghaffari }
305a515125bSLeila Ghaffari 
3068b774af8SJames Wright PetscErrorCode TSPostStep_CheckStep(TS ts) {
3078b774af8SJames Wright   Honee     honee;
3088b774af8SJames Wright   PetscReal norm;
3098b774af8SJames Wright   PetscInt  step;
3108b774af8SJames Wright   Vec       Q;
3118b774af8SJames Wright 
3128b774af8SJames Wright   PetscFunctionBeginUser;
3138b774af8SJames Wright   PetscCall(TSGetApplicationContext(ts, &honee));
3148b774af8SJames Wright   PetscCall(TSGetStepNumber(ts, &step));
3158b774af8SJames Wright   PetscCall(TSGetSolution(ts, &Q));
3168b774af8SJames Wright   if (step % honee->app_ctx->check_step_interval) PetscFunctionReturn(PETSC_SUCCESS);
3178b774af8SJames Wright   PetscCall(VecNorm(Q, NORM_1, &norm));
318827052a7SJames Wright   if (PetscIsInfOrNanReal(norm)) {
3198b774af8SJames Wright     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "Solution diverged: Nans found in solution\n"));
3208b774af8SJames Wright     PetscCall(TSSetConvergedReason(ts, TS_DIVERGED_NONLINEAR_SOLVE));
3218b774af8SJames Wright   }
3228b774af8SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3238b774af8SJames Wright }
3248b774af8SJames Wright 
325354560d1SJames Wright PetscErrorCode TSPostStep_MaxWallTime(TS ts) {
326354560d1SJames Wright   Honee       honee;
327354560d1SJames Wright   PetscInt    step;
328354560d1SJames Wright   PetscMPIInt rank;
329354560d1SJames Wright   MPI_Comm    comm;
330354560d1SJames Wright   PetscBool   is_wall_time_exceeded = PETSC_FALSE;
331354560d1SJames Wright 
332354560d1SJames Wright   PetscFunctionBeginUser;
333354560d1SJames Wright   PetscCall(TSGetApplicationContext(ts, &honee));
334354560d1SJames Wright   PetscCall(TSGetStepNumber(ts, &step));
335354560d1SJames Wright   if (step % honee->max_wall_time_interval) PetscFunctionReturn(PETSC_SUCCESS);
336354560d1SJames Wright   PetscCall(PetscObjectGetComm((PetscObject)ts, &comm));
337354560d1SJames Wright   PetscCallMPI(MPI_Comm_rank(comm, &rank));
338354560d1SJames Wright   if (rank == 0) is_wall_time_exceeded = time(NULL) > honee->max_wall_time ? PETSC_TRUE : PETSC_FALSE;
339354560d1SJames Wright   // Must broadcast to avoid race condition
340354560d1SJames Wright   PetscCallMPI(MPI_Bcast(&is_wall_time_exceeded, 1, MPIU_BOOL, 0, comm));
341354560d1SJames Wright   if (PetscUnlikely(is_wall_time_exceeded)) {
342354560d1SJames Wright     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "Stopping TSSolve: Set max wall time exceeded\n"));
343354560d1SJames Wright     PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_USER));
344354560d1SJames Wright   }
345354560d1SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
346354560d1SJames Wright }
347354560d1SJames Wright 
3484c6ae86eSJames Wright /**
3494c6ae86eSJames Wright   @brief TSPostStep for HONEE
3504c6ae86eSJames Wright 
3514c6ae86eSJames Wright   `TSSetPostStep()` only accepts a single function argument, so this function groups together all post-step
3524c6ae86eSJames Wright   functionality needed for HONEE features
3534c6ae86eSJames Wright 
3544c6ae86eSJames Wright   @param[in] ts TS object
3554c6ae86eSJames Wright **/
3564c6ae86eSJames Wright PetscErrorCode TSPostStep_Honee(TS ts) {
3574c6ae86eSJames Wright   Honee honee;
3584c6ae86eSJames Wright 
3594c6ae86eSJames Wright   PetscFunctionBeginUser;
3604c6ae86eSJames Wright   PetscCall(TSGetApplicationContext(ts, &honee));
3614c6ae86eSJames Wright   if (honee->max_wall_time > 0) PetscCall(TSPostStep_MaxWallTime(ts));
3624c6ae86eSJames Wright   if (honee->app_ctx->sgs_train_enable) PetscCall(TSPostStep_SGS_DD_Training(ts));
3634c6ae86eSJames Wright   if (honee->app_ctx->check_step_interval > 0) PetscCall(TSPostStep_CheckStep(ts));
3644c6ae86eSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3654c6ae86eSJames Wright }
3664c6ae86eSJames Wright 
367*4e1dcb99SJames Wright /**
368*4e1dcb99SJames Wright   @brief Save `TSEvaluationSolutions()` to file defined by `-ts_eval_solutions_view`
369*4e1dcb99SJames Wright 
370*4e1dcb99SJames Wright   Intended to be used with `TSSetEvaluationTimes()`/`-ts_eval_times`
371*4e1dcb99SJames Wright 
372*4e1dcb99SJames Wright   @param[in] ts TS object that has the evaluation solutions
373*4e1dcb99SJames Wright **/
374*4e1dcb99SJames Wright static PetscErrorCode HoneeTSEvaluationSolutions(TS ts) {
375*4e1dcb99SJames Wright   MPI_Comm          comm = PetscObjectComm((PetscObject)ts);
376*4e1dcb99SJames Wright   const PetscReal  *sol_times;
377*4e1dcb99SJames Wright   PetscReal         orig_time;
378*4e1dcb99SJames Wright   PetscInt          num_sols, orig_step;
379*4e1dcb99SJames Wright   Vec              *sol_vecs;
380*4e1dcb99SJames Wright   DM                dm;
381*4e1dcb99SJames Wright   PetscBool         is_viewer_set;
382*4e1dcb99SJames Wright   PetscViewer       viewer;
383*4e1dcb99SJames Wright   PetscViewerFormat format;
384*4e1dcb99SJames Wright 
385*4e1dcb99SJames Wright   PetscFunctionBeginUser;
386*4e1dcb99SJames Wright   PetscCall(TSGetEvaluationSolutions(ts, &num_sols, &sol_times, &sol_vecs));
387*4e1dcb99SJames Wright   if (num_sols == 0) PetscFunctionReturn(PETSC_SUCCESS);
388*4e1dcb99SJames Wright   PetscCall(TSGetDM(ts, &dm));
389*4e1dcb99SJames Wright   PetscCall(DMGetOutputSequenceNumber(dm, &orig_step, &orig_time));
390*4e1dcb99SJames Wright 
391*4e1dcb99SJames Wright   const char *option = "-ts_eval_solutions_view";
392*4e1dcb99SJames Wright   PetscCall(PetscOptionsCreateViewer(comm, NULL, NULL, option, &viewer, &format, &is_viewer_set));
393*4e1dcb99SJames Wright   if (!is_viewer_set) {
394*4e1dcb99SJames Wright     PetscCall(PetscPrintf(comm, "\n\nWARNING: Viewer not set for TSEvaluationSolutions. Set %s to save this data. Now throwing it away!\n", option));
395*4e1dcb99SJames Wright     PetscFunctionReturn(PETSC_SUCCESS);
396*4e1dcb99SJames Wright   }
397*4e1dcb99SJames Wright   for (PetscInt i = 0; i < num_sols; i++) {
398*4e1dcb99SJames Wright     PetscCall(DMSetOutputSequenceNumber(dm, i, sol_times[i]));
399*4e1dcb99SJames Wright     PetscCall(VecView(sol_vecs[i], viewer));
400*4e1dcb99SJames Wright   }
401*4e1dcb99SJames Wright   PetscCall(PetscViewerPopFormat(viewer));
402*4e1dcb99SJames Wright   PetscCall(PetscViewerDestroy(&viewer));
403*4e1dcb99SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
404*4e1dcb99SJames Wright }
405*4e1dcb99SJames Wright 
406a515125bSLeila Ghaffari // TS: Create, setup, and solve
4072a9a4b51SJames Wright PetscErrorCode TSSolve_NS(DM dm, Honee honee, AppCtx app_ctx, Physics phys, ProblemData problem, Vec Q, PetscScalar *f_time, TS *ts) {
4080c373b74SJames Wright   MPI_Comm    comm = honee->comm;
409a515125bSLeila Ghaffari   TSAdapt     adapt;
410a515125bSLeila Ghaffari   PetscScalar final_time;
411a515125bSLeila Ghaffari 
41206f41313SJames Wright   PetscFunctionBeginUser;
4132b916ea7SJeremy L Thompson   PetscCall(TSCreate(comm, ts));
4142b916ea7SJeremy L Thompson   PetscCall(TSSetDM(*ts, dm));
4150c373b74SJames Wright   PetscCall(TSSetApplicationContext(*ts, honee));
416a515125bSLeila Ghaffari   if (phys->implicit) {
4172b916ea7SJeremy L Thompson     PetscCall(TSSetType(*ts, TSBDF));
4180c373b74SJames Wright     if (honee->op_ifunction) {
4190c373b74SJames Wright       PetscCall(TSSetIFunction(*ts, NULL, IFunction_NS, &honee));
420a515125bSLeila Ghaffari     } else {  // Implicit integrators can fall back to using an RHSFunction
4210c373b74SJames Wright       PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &honee));
422a515125bSLeila Ghaffari     }
4230c373b74SJames Wright     if (honee->mat_ijacobian) {
4240c373b74SJames Wright       PetscCall(DMTSSetIJacobian(dm, FormIJacobian_NS, &honee));
425f0b65372SJed Brown     }
426a515125bSLeila Ghaffari   } else {
4270c373b74SJames Wright     PetscCheck(honee->op_rhs_ctx, comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction");
4282b916ea7SJeremy L Thompson     PetscCall(TSSetType(*ts, TSRK));
4292b916ea7SJeremy L Thompson     PetscCall(TSRKSetType(*ts, TSRK5F));
4300c373b74SJames Wright     PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &honee));
431a515125bSLeila Ghaffari   }
4320c373b74SJames Wright   PetscCall(TSSetMaxTime(*ts, 500. * honee->units->second));
4332b916ea7SJeremy L Thompson   PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER));
43422387d3aSJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(TSSetErrorIfStepFails(*ts, PETSC_FALSE));
4350c373b74SJames Wright   PetscCall(TSSetTimeStep(*ts, 1.e-2 * honee->units->second));
4362b916ea7SJeremy L Thompson   PetscCall(TSGetAdapt(*ts, &adapt));
4370c373b74SJames Wright   PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * honee->units->second, 1.e2 * honee->units->second));
4382b916ea7SJeremy L Thompson   PetscCall(TSSetFromOptions(*ts));
4390c373b74SJames Wright   if (honee->mat_ijacobian) {
440ebfabadfSJames Wright     if (app_ctx->amat_type && !strcmp(app_ctx->amat_type, MATSHELL)) {
441ebfabadfSJames Wright       SNES snes;
442ebfabadfSJames Wright       KSP  ksp;
4433170c09fSJames Wright       Mat  Pmat, Amat;
444ebfabadfSJames Wright 
445ebfabadfSJames Wright       PetscCall(TSGetSNES(*ts, &snes));
446ebfabadfSJames Wright       PetscCall(SNESGetKSP(snes, &ksp));
4470c373b74SJames Wright       PetscCall(CreateSolveOperatorsFromMatCeed(ksp, honee->mat_ijacobian, PETSC_FALSE, &Amat, &Pmat));
4480c373b74SJames Wright       PetscCall(TSSetIJacobian(*ts, honee->mat_ijacobian, Pmat, NULL, NULL));
4493170c09fSJames Wright       PetscCall(MatDestroy(&Amat));
4503170c09fSJames Wright       PetscCall(MatDestroy(&Pmat));
451ebfabadfSJames Wright     }
452ebfabadfSJames Wright   }
4530c373b74SJames Wright   honee->time_bc_set = -1.0;  // require all BCs be updated
454c26b555cSJames Wright   if (app_ctx->cont_steps) {  // continue from previous timestep data
4550c373b74SJames Wright     PetscCall(TSSetTime(*ts, app_ctx->cont_time * honee->units->second));
45674a6f4ddSJed Brown     PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps));
457a515125bSLeila Ghaffari   }
45816cb6b6bSJames Wright 
45916cb6b6bSJames Wright   PetscBool add_ksp_postsolve_residual = PETSC_FALSE;
46016cb6b6bSJames Wright   PetscCall(PetscOptionsGetBool(NULL, NULL, "-ksp_post_solve_residual", &add_ksp_postsolve_residual, NULL));
46116cb6b6bSJames Wright   if (add_ksp_postsolve_residual) {
46216cb6b6bSJames Wright     SNES snes;
46316cb6b6bSJames Wright     KSP  ksp;
46416cb6b6bSJames Wright 
46516cb6b6bSJames Wright     PetscCall(TSGetSNES(*ts, &snes));
46616cb6b6bSJames Wright     PetscCall(SNESGetKSP(snes, &ksp));
46716cb6b6bSJames Wright     PetscCall(KSPSetResidualHistory(ksp, NULL, PETSC_DECIDE, PETSC_TRUE));
46816cb6b6bSJames Wright     PetscCall(KSPSetPostSolve(ksp, KSPPostSolve_Honee, NULL));
46916cb6b6bSJames Wright   }
4704c6ae86eSJames Wright   if (honee->set_poststep) PetscCall(TSSetPostStep(*ts, TSPostStep_Honee));
4714c6ae86eSJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(TSMonitorSet(*ts, TSMonitor_NS, honee, NULL));
4724c6ae86eSJames Wright   if (app_ctx->wall_forces.viewer) PetscCall(TSMonitorSet(*ts, TSMonitor_WallForce, honee, NULL));
473c931fa59SJames Wright   if (app_ctx->turb_spanstats_enable) {
4748fbf1d3fSJames Wright     PetscReal start_time;
4758fbf1d3fSJames Wright     PetscInt  start_step;
4768fbf1d3fSJames Wright 
4770c373b74SJames Wright     PetscCall(TSMonitorSet(*ts, TSMonitor_TurbulenceStatistics, honee, NULL));
4788fbf1d3fSJames Wright     PetscCall(TSGetTime(*ts, &start_time));
4798fbf1d3fSJames Wright     PetscCall(TSGetStepNumber(*ts, &start_step));
4808fbf1d3fSJames Wright     CeedScalar initial_solution_time = start_time;  // done for type conversion
481e987420bSJames Wright     PetscCallCeed(honee->ceed, CeedOperatorSetContextDouble(honee->spanstats->op_stats_collect_ctx->op, honee->spanstats->previous_time_label,
4828fbf1d3fSJames Wright                                                             &initial_solution_time));
483e987420bSJames Wright     honee->spanstats->initial_solution_time = start_time;
484e987420bSJames Wright     honee->spanstats->initial_solution_step = start_step;
485b0488d1fSJames Wright   }
48625125139SJames Wright   PetscCall(TSMonitorSetFromOptions(*ts, "-ts_monitor_total_kinetic_energy", "Monitor total kinetic energy balance terms in the domain", NULL,
48725125139SJames Wright                                     TSMonitor_TotalKineticEnergy, SetupMontiorTotalKineticEnergy));
48887fd7f33SJames Wright   PetscCall(TSMonitorSetFromOptions(*ts, "-ts_monitor_cfl", "Monitor element CFL statistics", NULL, TSMonitor_Cfl, SetupMontiorCfl));
4890c373b74SJames Wright   if (app_ctx->diff_filter_monitor) PetscCall(TSMonitorSet(*ts, TSMonitor_DifferentialFilter, honee, NULL));
4904c6ae86eSJames Wright   if (app_ctx->sgs_train_enable) PetscCall(TSMonitorSet(*ts, TSMonitor_SGS_DD_Training, honee, NULL));
491e539bbbaSJames Wright 
4920c373b74SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(PrintRunInfo(honee, honee->phys, problem, *ts));
49371f2ed29SJames Wright   if (honee->mass_ksp) PetscCall(KSPViewFromOptions(honee->mass_ksp, NULL, "-ksp_view_pre_ts_solve"));
49471f2ed29SJames Wright 
495a515125bSLeila Ghaffari   // Solve
49674a6f4ddSJed Brown   PetscReal start_time;
49774a6f4ddSJed Brown   PetscInt  start_step;
4982b916ea7SJeremy L Thompson   PetscCall(TSGetTime(*ts, &start_time));
49974a6f4ddSJed Brown   PetscCall(TSGetStepNumber(*ts, &start_step));
50091982731SJeremy L Thompson 
501df4304b5SJed Brown   PetscCall(PetscLogDefaultBegin());  // So we can use PetscLogStageGetPerfInfo without -log_view
502ea615d4cSJames Wright   PetscPreLoadBegin(PETSC_FALSE, "HONEE Solve");
50391982731SJeremy L Thompson   PetscCall(TSSetTime(*ts, start_time));
50474a6f4ddSJed Brown   PetscCall(TSSetStepNumber(*ts, start_step));
50591982731SJeremy L Thompson   if (PetscPreLoadingOn) {
50691982731SJeremy L Thompson     // LCOV_EXCL_START
50791982731SJeremy L Thompson     SNES      snes;
5084f25acd9SJames Wright     KSP       ksp;
50991982731SJeremy L Thompson     Vec       Q_preload;
5104f25acd9SJames Wright     PetscReal rtol_snes, rtol_ksp;
5112a9a4b51SJames Wright     PetscCall(VecDuplicate(Q, &Q_preload));
5122a9a4b51SJames Wright     PetscCall(VecCopy(Q, Q_preload));
51391982731SJeremy L Thompson     PetscCall(TSGetSNES(*ts, &snes));
5144f25acd9SJames Wright     PetscCall(SNESGetTolerances(snes, NULL, &rtol_snes, NULL, NULL, NULL));
5154f25acd9SJames Wright     PetscCall(SNESSetTolerances(snes, PETSC_CURRENT, .99, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT));
5164f25acd9SJames Wright     PetscCall(SNESGetKSP(snes, &ksp));
5174f25acd9SJames Wright     PetscCall(KSPGetTolerances(ksp, &rtol_ksp, NULL, NULL, NULL));
5184f25acd9SJames Wright     PetscCall(KSPSetTolerances(ksp, .99, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT));
51922c1b34eSJames Wright     PetscCall(TSSetSolution(*ts, Q_preload));
52091982731SJeremy L Thompson     PetscCall(TSStep(*ts));
5214f25acd9SJames Wright     PetscCall(SNESSetTolerances(snes, PETSC_CURRENT, rtol_snes, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT));
5224f25acd9SJames Wright     PetscCall(KSPSetTolerances(ksp, rtol_ksp, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT));
52391982731SJeremy L Thompson     PetscCall(VecDestroy(&Q_preload));
52491982731SJeremy L Thompson     // LCOV_EXCL_STOP
52591982731SJeremy L Thompson   } else {
5262b916ea7SJeremy L Thompson     PetscCall(PetscBarrier((PetscObject)*ts));
5272a9a4b51SJames Wright     PetscCall(TSSolve(*ts, Q));
52891982731SJeremy L Thompson   }
52991982731SJeremy L Thompson   PetscPreLoadEnd();
53091982731SJeremy L Thompson 
53191982731SJeremy L Thompson   PetscCall(TSGetSolveTime(*ts, &final_time));
532a515125bSLeila Ghaffari   *f_time = final_time;
53391982731SJeremy L Thompson 
5340e1e9333SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) {
5357538d537SJames Wright     PetscInt           step_no;
536*4e1dcb99SJames Wright     PetscLogStage      stage_id;
537*4e1dcb99SJames Wright     PetscEventPerfInfo stage_perf;
538*4e1dcb99SJames Wright 
5397538d537SJames Wright     PetscCall(TSGetStepNumber(*ts, &step_no));
5400c373b74SJames Wright     if (honee->app_ctx->checkpoint_interval > 0 || honee->app_ctx->checkpoint_interval == -1) {
5412a9a4b51SJames Wright       PetscCall(WriteOutput(honee, Q, step_no, final_time));
5427538d537SJames Wright     }
543*4e1dcb99SJames Wright     PetscCall(HoneeTSEvaluationSolutions(*ts));
54491982731SJeremy L Thompson 
545ea615d4cSJames Wright     PetscCall(PetscLogStageGetId("HONEE Solve", &stage_id));
546df4304b5SJed Brown     PetscCall(PetscLogStageGetPerfInfo(stage_id, &stage_perf));
547df4304b5SJed Brown     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time taken for solution (sec): %g\n", stage_perf.time));
548a515125bSLeila Ghaffari   }
549d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
550a515125bSLeila Ghaffari }
551