xref: /honee/src/setupts.c (revision 8d78d7c8b0f3930b8658b7ec1d601fd78b6e01eb)
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 
11*8d78d7c8SJames Wright #include <differential_filter.h>
12149fb536SJames Wright #include <navierstokes.h>
13c5e9980aSAdeleke O. Bankole #include "../qfunctions/newtonian_state.h"
14a515125bSLeila Ghaffari 
154cb09fddSJames Wright // @brief Insert Boundary values if it's a new time
160c373b74SJames Wright PetscErrorCode UpdateBoundaryValues(Honee honee, Vec Q_loc, PetscReal t) {
17c996854bSJames Wright   PetscFunctionBeginUser;
180c373b74SJames Wright   if (honee->time_bc_set != t) {
190c373b74SJames Wright     PetscCall(DMPlexInsertBoundaryValues(honee->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL));
200c373b74SJames Wright     honee->time_bc_set = t;
21c996854bSJames Wright   }
22d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
23c996854bSJames Wright }
24c996854bSJames Wright 
25a515125bSLeila Ghaffari // RHS (Explicit time-stepper) function setup
26a515125bSLeila Ghaffari //   This is the RHS of the ODE, given as u_t = G(t,u)
27a515125bSLeila Ghaffari //   This function takes in a state vector Q and writes into G
28a515125bSLeila Ghaffari PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) {
290c373b74SJames Wright   Honee        honee = *(Honee *)user_data;
300c373b74SJames Wright   Ceed         ceed  = honee->ceed;
31fd969b44SJames Wright   PetscScalar  dt;
32dc3c760aSJames Wright   Vec          Q_loc = honee->Q_loc, R;
33a78efa86SJames Wright   PetscMemType q_mem_type;
34a515125bSLeila Ghaffari 
3506f41313SJames Wright   PetscFunctionBeginUser;
36e2f84137SJeremy L Thompson   // Update time dependent data
370c373b74SJames Wright   PetscCall(UpdateBoundaryValues(honee, Q_loc, t));
380c373b74SJames Wright   if (honee->phys->solution_time_label)
390c373b74SJames Wright     PetscCallCeed(ceed, CeedOperatorSetContextDouble(honee->op_rhs_ctx->op, honee->phys->solution_time_label, &t));
402b916ea7SJeremy L Thompson   PetscCall(TSGetTimeStep(ts, &dt));
410c373b74SJames Wright   if (honee->phys->timestep_size_label)
420c373b74SJames Wright     PetscCallCeed(ceed, CeedOperatorSetContextDouble(honee->op_rhs_ctx->op, honee->phys->timestep_size_label, &dt));
43a515125bSLeila Ghaffari 
44dc3c760aSJames Wright   PetscCall(DMGetNamedGlobalVector(honee->dm, "RHS Residual", &R));
450c373b74SJames Wright   PetscCall(DMGlobalToLocal(honee->dm, Q, INSERT_VALUES, Q_loc));
460c373b74SJames Wright   if (honee->app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) PetscCall(DivDiffFluxProjectionApply(honee->diff_flux_proj, Q_loc));
47dc3c760aSJames Wright   PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, R, honee->op_rhs_ctx));
48a78efa86SJames Wright 
49a78efa86SJames Wright   // Inverse of the mass matrix
50dc3c760aSJames Wright   PetscCall(VecReadPetscToCeed(Q_loc, &q_mem_type, honee->q_ceed));
51bc5573aaSJames Wright   {
52bc5573aaSJames Wright     // Run PCApply manually if using ksp_type preonly -pc_type jacobi
53bc5573aaSJames Wright     // This is to avoid an AllReduce call in KSPSolve_Preonly, which causes significant slowdowns for lumped mass matrix solves.
54bc5573aaSJames Wright     // See https://gitlab.com/petsc/petsc/-/merge_requests/8048 for more details and a possible fix
55bc5573aaSJames Wright     PC        pc;
56bc5573aaSJames Wright     PetscBool ispreonly, isjacobi;
57bc5573aaSJames Wright     PetscCall(KSPGetPC(honee->mass_ksp, &pc));
58bc5573aaSJames Wright     PetscCall(PetscObjectTypeCompare((PetscObject)honee->mass_ksp, KSPPREONLY, &ispreonly));
59bc5573aaSJames Wright     PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &isjacobi));
60bc5573aaSJames Wright     if (ispreonly && isjacobi) PetscCall(PCApply(pc, R, G));
61bc5573aaSJames Wright     else PetscCall(KSPSolve(honee->mass_ksp, R, G));
62bc5573aaSJames Wright   }
630c373b74SJames Wright   PetscCall(VecReadCeedToPetsc(honee->q_ceed, q_mem_type, Q_loc));
64dc3c760aSJames Wright 
65dc3c760aSJames Wright   PetscCall(DMRestoreNamedGlobalVector(honee->dm, "RHS Residual", &R));
66d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
67a515125bSLeila Ghaffari }
68a515125bSLeila Ghaffari 
69c5e9980aSAdeleke O. Bankole // Surface forces function setup
70c5e9980aSAdeleke O. Bankole static PetscErrorCode Surface_Forces_NS(DM dm, Vec G_loc, PetscInt num_walls, const PetscInt walls[], PetscScalar *reaction_force) {
71c5e9980aSAdeleke O. Bankole   DMLabel            face_label;
7262ab3c0bSJames Wright   const PetscScalar *g_array;
7362ab3c0bSJames Wright   PetscInt           dim  = 3;
7462ab3c0bSJames Wright   MPI_Comm           comm = PetscObjectComm((PetscObject)dm);
7562ab3c0bSJames Wright   PetscSection       section;
76c5e9980aSAdeleke O. Bankole 
77c5e9980aSAdeleke O. Bankole   PetscFunctionBeginUser;
78c5e9980aSAdeleke O. Bankole   PetscCall(DMGetLabel(dm, "Face Sets", &face_label));
7962ab3c0bSJames Wright   PetscCall(VecGetArrayRead(G_loc, &g_array));
80c5e9980aSAdeleke O. Bankole   for (PetscInt w = 0; w < num_walls; w++) {
8162ab3c0bSJames Wright     const PetscInt wall = walls[w], *points;
82c5e9980aSAdeleke O. Bankole     IS             wall_is;
8362ab3c0bSJames Wright     PetscInt       num_points, num_comp = 0;
8462ab3c0bSJames Wright 
85c5e9980aSAdeleke O. Bankole     PetscCall(DMLabelGetStratumIS(face_label, wall, &wall_is));
8662ab3c0bSJames Wright     if (!wall_is) continue;  // No wall points on this process, skip
8762ab3c0bSJames Wright 
8862ab3c0bSJames Wright     PetscCall(DMGetLocalSection(dm, &section));
8962ab3c0bSJames Wright     PetscCall(PetscSectionGetFieldComponents(section, 0, &num_comp));
90c5e9980aSAdeleke O. Bankole     PetscCall(ISGetSize(wall_is, &num_points));
91c5e9980aSAdeleke O. Bankole     PetscCall(ISGetIndices(wall_is, &points));
92c5e9980aSAdeleke O. Bankole     for (PetscInt i = 0; i < num_points; i++) {
93c5e9980aSAdeleke O. Bankole       const PetscInt           p = points[i];
94c5e9980aSAdeleke O. Bankole       const StateConservative *r;
9562ab3c0bSJames Wright       PetscInt                 dof;
9662ab3c0bSJames Wright 
9762ab3c0bSJames Wright       PetscCall(DMPlexPointLocalRead(dm, p, g_array, &r));
9862ab3c0bSJames Wright       PetscCall(PetscSectionGetDof(section, p, &dof));
992004e3acSAdeleke O. Bankole       for (PetscInt node = 0; node < dof / num_comp; node++) {
10062ab3c0bSJames Wright         for (PetscInt j = 0; j < dim; j++) {
1012004e3acSAdeleke O. Bankole           reaction_force[w * dim + j] -= r[node].momentum[j];
1022004e3acSAdeleke O. Bankole         }
103c5e9980aSAdeleke O. Bankole       }
104c5e9980aSAdeleke O. Bankole     }
105c5e9980aSAdeleke O. Bankole     PetscCall(ISRestoreIndices(wall_is, &points));
106c5e9980aSAdeleke O. Bankole     PetscCall(ISDestroy(&wall_is));
107c5e9980aSAdeleke O. Bankole   }
108c5e9980aSAdeleke O. Bankole   PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, reaction_force, dim * num_walls, MPIU_SCALAR, MPI_SUM, comm));
10962ab3c0bSJames Wright   PetscCall(VecRestoreArrayRead(G_loc, &g_array));
110d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
111c5e9980aSAdeleke O. Bankole }
112c5e9980aSAdeleke O. Bankole 
113a515125bSLeila Ghaffari // Implicit time-stepper function setup
1142b916ea7SJeremy L Thompson PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) {
1150c373b74SJames Wright   Honee        honee = *(Honee *)user_data;
1160c373b74SJames Wright   Ceed         ceed  = honee->ceed;
117fd969b44SJames Wright   PetscScalar  dt;
1180c373b74SJames Wright   Vec          Q_loc = honee->Q_loc, Q_dot_loc = honee->Q_dot_loc, G_loc;
119a515125bSLeila Ghaffari   PetscMemType q_mem_type, q_dot_mem_type, g_mem_type;
120a515125bSLeila Ghaffari 
12106f41313SJames Wright   PetscFunctionBeginUser;
122dabd2275SJames Wright   PetscCall(DMGlobalToLocalBegin(honee->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
1230c373b74SJames Wright   PetscCall(DMGetNamedLocalVector(honee->dm, "ResidualLocal", &G_loc));
124e2f84137SJeremy L Thompson 
125e2f84137SJeremy L Thompson   // Update time dependent data
1260c373b74SJames Wright   PetscCall(UpdateBoundaryValues(honee, Q_loc, t));
1270c373b74SJames Wright   if (honee->phys->solution_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(honee->op_ifunction, honee->phys->solution_time_label, &t));
1282b916ea7SJeremy L Thompson   PetscCall(TSGetTimeStep(ts, &dt));
1290c373b74SJames Wright   if (honee->phys->timestep_size_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(honee->op_ifunction, honee->phys->timestep_size_label, &dt));
130a515125bSLeila Ghaffari 
131a515125bSLeila Ghaffari   // Global-to-local
1320c373b74SJames Wright   PetscCall(DMGlobalToLocalBegin(honee->dm, Q, INSERT_VALUES, Q_loc));
1330c373b74SJames Wright   PetscCall(DMGlobalToLocalEnd(honee->dm, Q, INSERT_VALUES, Q_loc));
1340c373b74SJames Wright   if (honee->app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) PetscCall(DivDiffFluxProjectionApply(honee->diff_flux_proj, Q_loc));
1350c373b74SJames Wright   PetscCall(DMGlobalToLocalEnd(honee->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
136a515125bSLeila Ghaffari 
137a515125bSLeila Ghaffari   // Place PETSc vectors in CEED vectors
1380c373b74SJames Wright   PetscCall(VecReadPetscToCeed(Q_loc, &q_mem_type, honee->q_ceed));
1390c373b74SJames Wright   PetscCall(VecReadPetscToCeed(Q_dot_loc, &q_dot_mem_type, honee->q_dot_ceed));
1400c373b74SJames Wright   PetscCall(VecPetscToCeed(G_loc, &g_mem_type, honee->g_ceed));
141a515125bSLeila Ghaffari 
142a515125bSLeila Ghaffari   // Apply CEED operator
143ea615d4cSJames Wright   PetscCall(PetscLogEventBegin(HONEE_CeedOperatorApply, Q, G, 0, 0));
1447eedc94cSJames Wright   PetscCall(PetscLogGpuTimeBegin());
1450c373b74SJames Wright   PetscCallCeed(honee->ceed, CeedOperatorApply(honee->op_ifunction, honee->q_ceed, honee->g_ceed, CEED_REQUEST_IMMEDIATE));
1467eedc94cSJames Wright   PetscCall(PetscLogGpuTimeEnd());
147ea615d4cSJames Wright   PetscCall(PetscLogEventEnd(HONEE_CeedOperatorApply, Q, G, 0, 0));
148a515125bSLeila Ghaffari 
149a515125bSLeila Ghaffari   // Restore vectors
1500c373b74SJames Wright   PetscCall(VecReadCeedToPetsc(honee->q_ceed, q_mem_type, Q_loc));
1510c373b74SJames Wright   PetscCall(VecReadCeedToPetsc(honee->q_dot_ceed, q_dot_mem_type, Q_dot_loc));
1520c373b74SJames Wright   PetscCall(VecCeedToPetsc(honee->g_ceed, g_mem_type, G_loc));
153a515125bSLeila Ghaffari 
1540c373b74SJames Wright   if (honee->app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) {
1550c373b74SJames Wright     PetscCall(SgsDDApplyIFunction(honee, Q_loc, G_loc));
15601ab89c1SJames Wright   }
1579c678832SJames Wright 
158a515125bSLeila Ghaffari   // Local-to-Global
1592b916ea7SJeremy L Thompson   PetscCall(VecZeroEntries(G));
1600c373b74SJames Wright   PetscCall(DMLocalToGlobal(honee->dm, G_loc, ADD_VALUES, G));
161a515125bSLeila Ghaffari 
162a515125bSLeila Ghaffari   // Restore vectors
1630c373b74SJames Wright   PetscCall(DMRestoreNamedLocalVector(honee->dm, "ResidualLocal", &G_loc));
164d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
165a515125bSLeila Ghaffari }
166a515125bSLeila Ghaffari 
1672b916ea7SJeremy L Thompson PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, PetscReal shift, Mat J, Mat J_pre, void *user_data) {
1680c373b74SJames Wright   Honee     honee = *(Honee *)user_data;
169ebfabadfSJames Wright   PetscBool J_is_matceed, J_is_mffd, J_pre_is_matceed, J_pre_is_mffd;
17006f41313SJames Wright 
171f0b65372SJed Brown   PetscFunctionBeginUser;
17204855949SJed Brown   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATMFFD, &J_is_mffd));
173ebfabadfSJames Wright   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATCEED, &J_is_matceed));
174ebfabadfSJames Wright   PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATMFFD, &J_pre_is_mffd));
175ebfabadfSJames Wright   PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATCEED, &J_pre_is_matceed));
176ebfabadfSJames Wright 
1770c373b74SJames Wright   PetscCall(MatCeedSetContextReal(honee->mat_ijacobian, "ijacobian time shift", shift));
178ebfabadfSJames Wright 
179ebfabadfSJames Wright   if (J_is_matceed || J_is_mffd) {
18004855949SJed Brown     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
18104855949SJed Brown     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
1820c373b74SJames Wright   } else PetscCall(MatCeedAssembleCOO(honee->mat_ijacobian, J));
183ebfabadfSJames Wright 
184ebfabadfSJames Wright   if (J_pre_is_matceed && J != J_pre) {
185ebfabadfSJames Wright     PetscCall(MatAssemblyBegin(J_pre, MAT_FINAL_ASSEMBLY));
186ebfabadfSJames Wright     PetscCall(MatAssemblyEnd(J_pre, MAT_FINAL_ASSEMBLY));
187ebfabadfSJames Wright   } else if (!J_pre_is_matceed && !J_pre_is_mffd && J != J_pre) {
1880c373b74SJames Wright     PetscCall(MatCeedAssembleCOO(honee->mat_ijacobian, J_pre));
18904855949SJed Brown   }
190d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
191f0b65372SJed Brown }
192f0b65372SJed Brown 
1930c373b74SJames Wright PetscErrorCode WriteOutput(Honee honee, Vec Q, PetscInt step_no, PetscScalar time) {
194a515125bSLeila Ghaffari   Vec         Q_loc;
195a515125bSLeila Ghaffari   char        file_path[PETSC_MAX_PATH_LEN];
196a515125bSLeila Ghaffari   PetscViewer viewer;
197a515125bSLeila Ghaffari 
19806f41313SJames Wright   PetscFunctionBeginUser;
1990c373b74SJames Wright   if (honee->app_ctx->checkpoint_vtk) {
200a515125bSLeila Ghaffari     // Set up output
2010c373b74SJames Wright     PetscCall(DMGetLocalVector(honee->dm, &Q_loc));
2027538d537SJames Wright     PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec"));
2037538d537SJames Wright     PetscCall(VecZeroEntries(Q_loc));
2040c373b74SJames Wright     PetscCall(DMGlobalToLocal(honee->dm, Q, INSERT_VALUES, Q_loc));
205a515125bSLeila Ghaffari 
206a515125bSLeila Ghaffari     // Output
2070c373b74SJames Wright     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-%03" PetscInt_FMT ".vtu", honee->app_ctx->output_dir, step_no));
2087538d537SJames Wright 
2092b916ea7SJeremy L Thompson     PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path, FILE_MODE_WRITE, &viewer));
2107538d537SJames Wright     PetscCall(VecView(Q_loc, viewer));
2117538d537SJames Wright     PetscCall(PetscViewerDestroy(&viewer));
2120c373b74SJames Wright     if (honee->dm_viz) {
213a515125bSLeila Ghaffari       Vec         Q_refined, Q_refined_loc;
214a515125bSLeila Ghaffari       char        file_path_refined[PETSC_MAX_PATH_LEN];
215a515125bSLeila Ghaffari       PetscViewer viewer_refined;
216a515125bSLeila Ghaffari 
2170c373b74SJames Wright       PetscCall(DMGetGlobalVector(honee->dm_viz, &Q_refined));
2180c373b74SJames Wright       PetscCall(DMGetLocalVector(honee->dm_viz, &Q_refined_loc));
2197538d537SJames Wright       PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined"));
2207538d537SJames Wright 
2210c373b74SJames Wright       PetscCall(MatInterpolate(honee->interp_viz, Q, Q_refined));
2227538d537SJames Wright       PetscCall(VecZeroEntries(Q_refined_loc));
2230c373b74SJames Wright       PetscCall(DMGlobalToLocal(honee->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc));
2247538d537SJames Wright 
225852e5969SJed Brown       PetscCall(
2260c373b74SJames Wright           PetscSNPrintf(file_path_refined, sizeof file_path_refined, "%s/nsrefined-%03" PetscInt_FMT ".vtu", honee->app_ctx->output_dir, step_no));
2277538d537SJames Wright 
2282b916ea7SJeremy L Thompson       PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined), file_path_refined, FILE_MODE_WRITE, &viewer_refined));
2297538d537SJames Wright       PetscCall(VecView(Q_refined_loc, viewer_refined));
2300c373b74SJames Wright       PetscCall(DMRestoreLocalVector(honee->dm_viz, &Q_refined_loc));
2310c373b74SJames Wright       PetscCall(DMRestoreGlobalVector(honee->dm_viz, &Q_refined));
2327538d537SJames Wright       PetscCall(PetscViewerDestroy(&viewer_refined));
233a515125bSLeila Ghaffari     }
2340c373b74SJames Wright     PetscCall(DMRestoreLocalVector(honee->dm, &Q_loc));
235852e5969SJed Brown   }
236a515125bSLeila Ghaffari 
237a515125bSLeila Ghaffari   // Save data in a binary file for continuation of simulations
2380c373b74SJames Wright   if (honee->app_ctx->add_stepnum2bin) {
2390c373b74SJames Wright     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", honee->app_ctx->output_dir, step_no));
24091a36801SJames Wright   } else {
2410c373b74SJames Wright     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", honee->app_ctx->output_dir));
24291a36801SJames Wright   }
2430c373b74SJames Wright   PetscCall(PetscViewerBinaryOpen(honee->comm, file_path, FILE_MODE_WRITE, &viewer));
2447538d537SJames Wright 
2450c373b74SJames Wright   time /= honee->units->second;  // Dimensionalize time back
246b237916aSJames Wright   PetscCall(HoneeWriteBinaryVec(viewer, Q, time, step_no));
2477538d537SJames Wright   PetscCall(PetscViewerDestroy(&viewer));
248d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2497538d537SJames Wright }
2507538d537SJames Wright 
251c5e9980aSAdeleke O. Bankole // CSV Monitor
252c5e9980aSAdeleke O. Bankole PetscErrorCode TSMonitor_WallForce(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
2530c373b74SJames Wright   Honee             honee = ctx;
254c5e9980aSAdeleke O. Bankole   Vec               G_loc;
2550c373b74SJames Wright   PetscInt          num_wall = honee->app_ctx->wall_forces.num_wall, dim = 3;
2560c373b74SJames Wright   const PetscInt   *walls  = honee->app_ctx->wall_forces.walls;
2570c373b74SJames Wright   PetscViewer       viewer = honee->app_ctx->wall_forces.viewer;
2580c373b74SJames Wright   PetscViewerFormat format = honee->app_ctx->wall_forces.viewer_format;
259c5e9980aSAdeleke O. Bankole   PetscScalar      *reaction_force;
26062ab3c0bSJames Wright   PetscBool         is_ascii;
261c5e9980aSAdeleke O. Bankole 
262c5e9980aSAdeleke O. Bankole   PetscFunctionBeginUser;
263d949ddfcSJames Wright   if (!viewer) PetscFunctionReturn(PETSC_SUCCESS);
2640c373b74SJames Wright   PetscCall(DMGetNamedLocalVector(honee->dm, "ResidualLocal", &G_loc));
26562ab3c0bSJames Wright   PetscCall(PetscCalloc1(num_wall * dim, &reaction_force));
2660c373b74SJames Wright   PetscCall(Surface_Forces_NS(honee->dm, G_loc, num_wall, walls, reaction_force));
2670c373b74SJames Wright   PetscCall(DMRestoreNamedLocalVector(honee->dm, "ResidualLocal", &G_loc));
268c5e9980aSAdeleke O. Bankole 
26962ab3c0bSJames Wright   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &is_ascii));
270c5e9980aSAdeleke O. Bankole 
27162ab3c0bSJames Wright   if (is_ascii) {
2720c373b74SJames Wright     if (format == PETSC_VIEWER_ASCII_CSV && !honee->app_ctx->wall_forces.header_written) {
273c5e9980aSAdeleke O. Bankole       PetscCall(PetscViewerASCIIPrintf(viewer, "Step,Time,Wall,ForceX,ForceY,ForceZ\n"));
2740c373b74SJames Wright       honee->app_ctx->wall_forces.header_written = PETSC_TRUE;
275c5e9980aSAdeleke O. Bankole     }
276c5e9980aSAdeleke O. Bankole     for (PetscInt w = 0; w < num_wall; w++) {
277c5e9980aSAdeleke O. Bankole       PetscInt wall = walls[w];
278c5e9980aSAdeleke O. Bankole       if (format == PETSC_VIEWER_ASCII_CSV) {
279c5e9980aSAdeleke O. Bankole         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT ",%g,%" PetscInt_FMT ",%g,%g,%g\n", step_no, time, wall,
280c5e9980aSAdeleke O. Bankole                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
281c5e9980aSAdeleke O. Bankole 
282c5e9980aSAdeleke O. Bankole       } else {
283c5e9980aSAdeleke O. Bankole         PetscCall(PetscViewerASCIIPrintf(viewer, "Wall %" PetscInt_FMT " Forces: Force_x = %12g, Force_y = %12g, Force_z = %12g\n", wall,
284c5e9980aSAdeleke O. Bankole                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
285c5e9980aSAdeleke O. Bankole       }
286c5e9980aSAdeleke O. Bankole     }
287c5e9980aSAdeleke O. Bankole   }
288c5e9980aSAdeleke O. Bankole   PetscCall(PetscFree(reaction_force));
289d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
290c5e9980aSAdeleke O. Bankole }
291c5e9980aSAdeleke O. Bankole 
2927538d537SJames Wright // User provided TS Monitor
2932b916ea7SJeremy L Thompson PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
2940c373b74SJames Wright   Honee honee = ctx;
2957538d537SJames Wright 
29606f41313SJames Wright   PetscFunctionBeginUser;
297852e5969SJed Brown   // Print every 'checkpoint_interval' steps
2980c373b74SJames Wright   if (honee->app_ctx->checkpoint_interval <= 0 || step_no % honee->app_ctx->checkpoint_interval != 0 ||
2990c373b74SJames Wright       (honee->app_ctx->cont_steps == step_no && step_no != 0)) {
300d949ddfcSJames Wright     PetscFunctionReturn(PETSC_SUCCESS);
301e419654dSJeremy L Thompson   }
3027538d537SJames Wright 
3030c373b74SJames Wright   PetscCall(WriteOutput(honee, Q, step_no, time));
304d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
305a515125bSLeila Ghaffari }
306a515125bSLeila Ghaffari 
3078b774af8SJames Wright PetscErrorCode TSPostStep_CheckStep(TS ts) {
3088b774af8SJames Wright   Honee     honee;
3098b774af8SJames Wright   PetscReal norm;
3108b774af8SJames Wright   PetscInt  step;
3118b774af8SJames Wright   Vec       Q;
3128b774af8SJames Wright 
3138b774af8SJames Wright   PetscFunctionBeginUser;
3148b774af8SJames Wright   PetscCall(TSGetApplicationContext(ts, &honee));
3158b774af8SJames Wright   PetscCall(TSGetStepNumber(ts, &step));
3168b774af8SJames Wright   PetscCall(TSGetSolution(ts, &Q));
3178b774af8SJames Wright   if (step % honee->app_ctx->check_step_interval) PetscFunctionReturn(PETSC_SUCCESS);
3188b774af8SJames Wright   PetscCall(VecNorm(Q, NORM_1, &norm));
319827052a7SJames Wright   if (PetscIsInfOrNanReal(norm)) {
3208b774af8SJames Wright     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "Solution diverged: Nans found in solution\n"));
3218b774af8SJames Wright     PetscCall(TSSetConvergedReason(ts, TS_DIVERGED_NONLINEAR_SOLVE));
3228b774af8SJames Wright   }
3238b774af8SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3248b774af8SJames Wright }
3258b774af8SJames Wright 
326354560d1SJames Wright PetscErrorCode TSPostStep_MaxWallTime(TS ts) {
327354560d1SJames Wright   Honee       honee;
328354560d1SJames Wright   PetscInt    step;
329354560d1SJames Wright   PetscMPIInt rank;
330354560d1SJames Wright   MPI_Comm    comm;
331354560d1SJames Wright   PetscBool   is_wall_time_exceeded = PETSC_FALSE;
332354560d1SJames Wright 
333354560d1SJames Wright   PetscFunctionBeginUser;
334354560d1SJames Wright   PetscCall(TSGetApplicationContext(ts, &honee));
335354560d1SJames Wright   PetscCall(TSGetStepNumber(ts, &step));
336354560d1SJames Wright   if (step % honee->max_wall_time_interval) PetscFunctionReturn(PETSC_SUCCESS);
337354560d1SJames Wright   PetscCall(PetscObjectGetComm((PetscObject)ts, &comm));
338354560d1SJames Wright   PetscCallMPI(MPI_Comm_rank(comm, &rank));
339354560d1SJames Wright   if (rank == 0) is_wall_time_exceeded = time(NULL) > honee->max_wall_time ? PETSC_TRUE : PETSC_FALSE;
340354560d1SJames Wright   // Must broadcast to avoid race condition
341354560d1SJames Wright   PetscCallMPI(MPI_Bcast(&is_wall_time_exceeded, 1, MPIU_BOOL, 0, comm));
342354560d1SJames Wright   if (PetscUnlikely(is_wall_time_exceeded)) {
343354560d1SJames Wright     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "Stopping TSSolve: Set max wall time exceeded\n"));
344354560d1SJames Wright     PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_USER));
345354560d1SJames Wright   }
346354560d1SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
347354560d1SJames Wright }
348354560d1SJames Wright 
3494c6ae86eSJames Wright /**
3504c6ae86eSJames Wright   @brief TSPostStep for HONEE
3514c6ae86eSJames Wright 
3524c6ae86eSJames Wright   `TSSetPostStep()` only accepts a single function argument, so this function groups together all post-step
3534c6ae86eSJames Wright   functionality needed for HONEE features
3544c6ae86eSJames Wright 
3554c6ae86eSJames Wright   @param[in] ts TS object
3564c6ae86eSJames Wright **/
3574c6ae86eSJames Wright PetscErrorCode TSPostStep_Honee(TS ts) {
3584c6ae86eSJames Wright   Honee honee;
3594c6ae86eSJames Wright 
3604c6ae86eSJames Wright   PetscFunctionBeginUser;
3614c6ae86eSJames Wright   PetscCall(TSGetApplicationContext(ts, &honee));
3624c6ae86eSJames Wright   if (honee->max_wall_time > 0) PetscCall(TSPostStep_MaxWallTime(ts));
3634c6ae86eSJames Wright   if (honee->app_ctx->sgs_train_enable) PetscCall(TSPostStep_SGS_DD_Training(ts));
3644c6ae86eSJames Wright   if (honee->app_ctx->check_step_interval > 0) PetscCall(TSPostStep_CheckStep(ts));
3654c6ae86eSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3664c6ae86eSJames Wright }
3674c6ae86eSJames Wright 
3684e1dcb99SJames Wright /**
3694e1dcb99SJames Wright   @brief Save `TSEvaluationSolutions()` to file defined by `-ts_eval_solutions_view`
3704e1dcb99SJames Wright 
3714e1dcb99SJames Wright   Intended to be used with `TSSetEvaluationTimes()`/`-ts_eval_times`
3724e1dcb99SJames Wright 
3734e1dcb99SJames Wright   @param[in] ts TS object that has the evaluation solutions
3744e1dcb99SJames Wright **/
3754e1dcb99SJames Wright static PetscErrorCode HoneeTSEvaluationSolutions(TS ts) {
3764e1dcb99SJames Wright   MPI_Comm          comm = PetscObjectComm((PetscObject)ts);
3774e1dcb99SJames Wright   const PetscReal  *sol_times;
3784e1dcb99SJames Wright   PetscReal         orig_time;
3794e1dcb99SJames Wright   PetscInt          num_sols, orig_step;
3804e1dcb99SJames Wright   Vec              *sol_vecs;
3814e1dcb99SJames Wright   DM                dm;
3824e1dcb99SJames Wright   PetscBool         is_viewer_set;
3834e1dcb99SJames Wright   PetscViewer       viewer;
3844e1dcb99SJames Wright   PetscViewerFormat format;
3854e1dcb99SJames Wright 
3864e1dcb99SJames Wright   PetscFunctionBeginUser;
3874e1dcb99SJames Wright   PetscCall(TSGetEvaluationSolutions(ts, &num_sols, &sol_times, &sol_vecs));
3884e1dcb99SJames Wright   if (num_sols == 0) PetscFunctionReturn(PETSC_SUCCESS);
3894e1dcb99SJames Wright   PetscCall(TSGetDM(ts, &dm));
3904e1dcb99SJames Wright   PetscCall(DMGetOutputSequenceNumber(dm, &orig_step, &orig_time));
3914e1dcb99SJames Wright 
3924e1dcb99SJames Wright   const char *option = "-ts_eval_solutions_view";
3934e1dcb99SJames Wright   PetscCall(PetscOptionsCreateViewer(comm, NULL, NULL, option, &viewer, &format, &is_viewer_set));
3944e1dcb99SJames Wright   if (!is_viewer_set) {
3954e1dcb99SJames Wright     PetscCall(PetscPrintf(comm, "\n\nWARNING: Viewer not set for TSEvaluationSolutions. Set %s to save this data. Now throwing it away!\n", option));
3964e1dcb99SJames Wright     PetscFunctionReturn(PETSC_SUCCESS);
3974e1dcb99SJames Wright   }
3984e1dcb99SJames Wright   for (PetscInt i = 0; i < num_sols; i++) {
3994e1dcb99SJames Wright     PetscCall(DMSetOutputSequenceNumber(dm, i, sol_times[i]));
4004e1dcb99SJames Wright     PetscCall(VecView(sol_vecs[i], viewer));
4014e1dcb99SJames Wright   }
4024e1dcb99SJames Wright   PetscCall(PetscViewerPopFormat(viewer));
4034e1dcb99SJames Wright   PetscCall(PetscViewerDestroy(&viewer));
4044e1dcb99SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
4054e1dcb99SJames Wright }
4064e1dcb99SJames Wright 
407a515125bSLeila Ghaffari // TS: Create, setup, and solve
4082a9a4b51SJames Wright PetscErrorCode TSSolve_NS(DM dm, Honee honee, AppCtx app_ctx, Physics phys, ProblemData problem, Vec Q, PetscScalar *f_time, TS *ts) {
4090c373b74SJames Wright   MPI_Comm    comm = honee->comm;
410a515125bSLeila Ghaffari   TSAdapt     adapt;
411a515125bSLeila Ghaffari   PetscScalar final_time;
412a515125bSLeila Ghaffari 
41306f41313SJames Wright   PetscFunctionBeginUser;
4142b916ea7SJeremy L Thompson   PetscCall(TSCreate(comm, ts));
4152b916ea7SJeremy L Thompson   PetscCall(TSSetDM(*ts, dm));
4160c373b74SJames Wright   PetscCall(TSSetApplicationContext(*ts, honee));
417a515125bSLeila Ghaffari   if (phys->implicit) {
4182b916ea7SJeremy L Thompson     PetscCall(TSSetType(*ts, TSBDF));
4190c373b74SJames Wright     if (honee->op_ifunction) {
4200c373b74SJames Wright       PetscCall(TSSetIFunction(*ts, NULL, IFunction_NS, &honee));
421a515125bSLeila Ghaffari     } else {  // Implicit integrators can fall back to using an RHSFunction
4220c373b74SJames Wright       PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &honee));
423a515125bSLeila Ghaffari     }
4240c373b74SJames Wright     if (honee->mat_ijacobian) {
4250c373b74SJames Wright       PetscCall(DMTSSetIJacobian(dm, FormIJacobian_NS, &honee));
426f0b65372SJed Brown     }
427a515125bSLeila Ghaffari   } else {
4280c373b74SJames Wright     PetscCheck(honee->op_rhs_ctx, comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction");
4292b916ea7SJeremy L Thompson     PetscCall(TSSetType(*ts, TSRK));
4302b916ea7SJeremy L Thompson     PetscCall(TSRKSetType(*ts, TSRK5F));
4310c373b74SJames Wright     PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &honee));
432a515125bSLeila Ghaffari   }
4330c373b74SJames Wright   PetscCall(TSSetMaxTime(*ts, 500. * honee->units->second));
4342b916ea7SJeremy L Thompson   PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER));
43522387d3aSJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(TSSetErrorIfStepFails(*ts, PETSC_FALSE));
4360c373b74SJames Wright   PetscCall(TSSetTimeStep(*ts, 1.e-2 * honee->units->second));
4372b916ea7SJeremy L Thompson   PetscCall(TSGetAdapt(*ts, &adapt));
4380c373b74SJames Wright   PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * honee->units->second, 1.e2 * honee->units->second));
4392b916ea7SJeremy L Thompson   PetscCall(TSSetFromOptions(*ts));
4400c373b74SJames Wright   if (honee->mat_ijacobian) {
441ebfabadfSJames Wright     if (app_ctx->amat_type && !strcmp(app_ctx->amat_type, MATSHELL)) {
442ebfabadfSJames Wright       SNES snes;
443ebfabadfSJames Wright       KSP  ksp;
4443170c09fSJames Wright       Mat  Pmat, Amat;
445ebfabadfSJames Wright 
446ebfabadfSJames Wright       PetscCall(TSGetSNES(*ts, &snes));
447ebfabadfSJames Wright       PetscCall(SNESGetKSP(snes, &ksp));
4480c373b74SJames Wright       PetscCall(CreateSolveOperatorsFromMatCeed(ksp, honee->mat_ijacobian, PETSC_FALSE, &Amat, &Pmat));
4490c373b74SJames Wright       PetscCall(TSSetIJacobian(*ts, honee->mat_ijacobian, Pmat, NULL, NULL));
4503170c09fSJames Wright       PetscCall(MatDestroy(&Amat));
4513170c09fSJames Wright       PetscCall(MatDestroy(&Pmat));
452ebfabadfSJames Wright     }
453ebfabadfSJames Wright   }
4540c373b74SJames Wright   honee->time_bc_set = -1.0;  // require all BCs be updated
455c26b555cSJames Wright   if (app_ctx->cont_steps) {  // continue from previous timestep data
4560c373b74SJames Wright     PetscCall(TSSetTime(*ts, app_ctx->cont_time * honee->units->second));
45774a6f4ddSJed Brown     PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps));
458a515125bSLeila Ghaffari   }
45916cb6b6bSJames Wright 
46016cb6b6bSJames Wright   PetscBool add_ksp_postsolve_residual = PETSC_FALSE;
46116cb6b6bSJames Wright   PetscCall(PetscOptionsGetBool(NULL, NULL, "-ksp_post_solve_residual", &add_ksp_postsolve_residual, NULL));
46216cb6b6bSJames Wright   if (add_ksp_postsolve_residual) {
46316cb6b6bSJames Wright     SNES snes;
46416cb6b6bSJames Wright     KSP  ksp;
46516cb6b6bSJames Wright 
46616cb6b6bSJames Wright     PetscCall(TSGetSNES(*ts, &snes));
46716cb6b6bSJames Wright     PetscCall(SNESGetKSP(snes, &ksp));
46816cb6b6bSJames Wright     PetscCall(KSPSetResidualHistory(ksp, NULL, PETSC_DECIDE, PETSC_TRUE));
46916cb6b6bSJames Wright     PetscCall(KSPSetPostSolve(ksp, KSPPostSolve_Honee, NULL));
47016cb6b6bSJames Wright   }
4714c6ae86eSJames Wright   if (honee->set_poststep) PetscCall(TSSetPostStep(*ts, TSPostStep_Honee));
4724c6ae86eSJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(TSMonitorSet(*ts, TSMonitor_NS, honee, NULL));
4734c6ae86eSJames Wright   if (app_ctx->wall_forces.viewer) PetscCall(TSMonitorSet(*ts, TSMonitor_WallForce, honee, NULL));
4748fbf1d3fSJames Wright 
47578c5b8e5SJames Wright   PetscOptionsBegin(comm, NULL, "HONEE TS Monitor Options", NULL);
47625125139SJames Wright   PetscCall(TSMonitorSetFromOptions(*ts, "-ts_monitor_total_kinetic_energy", "Monitor total kinetic energy balance terms in the domain", NULL,
47725125139SJames Wright                                     TSMonitor_TotalKineticEnergy, SetupMontiorTotalKineticEnergy));
47887fd7f33SJames Wright   PetscCall(TSMonitorSetFromOptions(*ts, "-ts_monitor_cfl", "Monitor element CFL statistics", NULL, TSMonitor_Cfl, SetupMontiorCfl));
47978c5b8e5SJames Wright 
48078c5b8e5SJames Wright   PetscCall(
48178c5b8e5SJames Wright       PetscOptionsDeprecated("-ts_monitor_turbulence_spanstats_viewer_interval", "-ts_monitor_spanstats_turbulence_interval", "HONEE 0.0", NULL));
48278c5b8e5SJames Wright   PetscCall(PetscOptionsDeprecated("-ts_monitor_turbulence_spanstats_viewer", "-ts_monitor_spanstats_turbulence", "HONEE 0.0", NULL));
48378c5b8e5SJames Wright   PetscCall(TSMonitorSetFromOptions(*ts, "-ts_monitor_spanstats_turbulence", "Setup spanwise statistics collection", NULL,
48478c5b8e5SJames Wright                                     TSMonitor_TurbulenceSpanwiseStatistics, SpanwiseStatisticsSetup_Turbulence));
485cb8a476cSJames Wright   PetscCall(TSMonitorSetFromOptions(*ts, "-diff_filter_monitor", "Run differential filtering for every timestep (used for testing)", NULL,
486cb8a476cSJames Wright                                     TSMonitor_DifferentialFilter, TSMonitor_DifferentialFilterSetup));
48778c5b8e5SJames Wright   PetscOptionsEnd();
48878c5b8e5SJames Wright 
4894c6ae86eSJames Wright   if (app_ctx->sgs_train_enable) PetscCall(TSMonitorSet(*ts, TSMonitor_SGS_DD_Training, honee, NULL));
490e539bbbaSJames Wright 
4910c373b74SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(PrintRunInfo(honee, honee->phys, problem, *ts));
49271f2ed29SJames Wright   if (honee->mass_ksp) PetscCall(KSPViewFromOptions(honee->mass_ksp, NULL, "-ksp_view_pre_ts_solve"));
49371f2ed29SJames Wright 
494a515125bSLeila Ghaffari   // Solve
49574a6f4ddSJed Brown   PetscReal start_time;
49674a6f4ddSJed Brown   PetscInt  start_step;
4972b916ea7SJeremy L Thompson   PetscCall(TSGetTime(*ts, &start_time));
49874a6f4ddSJed Brown   PetscCall(TSGetStepNumber(*ts, &start_step));
49991982731SJeremy L Thompson 
500df4304b5SJed Brown   PetscCall(PetscLogDefaultBegin());  // So we can use PetscLogStageGetPerfInfo without -log_view
501ea615d4cSJames Wright   PetscPreLoadBegin(PETSC_FALSE, "HONEE Solve");
50291982731SJeremy L Thompson   PetscCall(TSSetTime(*ts, start_time));
50374a6f4ddSJed Brown   PetscCall(TSSetStepNumber(*ts, start_step));
50491982731SJeremy L Thompson   if (PetscPreLoadingOn) {
50591982731SJeremy L Thompson     // LCOV_EXCL_START
50691982731SJeremy L Thompson     SNES      snes;
5074f25acd9SJames Wright     KSP       ksp;
50891982731SJeremy L Thompson     Vec       Q_preload;
5094f25acd9SJames Wright     PetscReal rtol_snes, rtol_ksp;
5102a9a4b51SJames Wright     PetscCall(VecDuplicate(Q, &Q_preload));
5112a9a4b51SJames Wright     PetscCall(VecCopy(Q, Q_preload));
51291982731SJeremy L Thompson     PetscCall(TSGetSNES(*ts, &snes));
5134f25acd9SJames Wright     PetscCall(SNESGetTolerances(snes, NULL, &rtol_snes, NULL, NULL, NULL));
5144f25acd9SJames Wright     PetscCall(SNESSetTolerances(snes, PETSC_CURRENT, .99, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT));
5154f25acd9SJames Wright     PetscCall(SNESGetKSP(snes, &ksp));
5164f25acd9SJames Wright     PetscCall(KSPGetTolerances(ksp, &rtol_ksp, NULL, NULL, NULL));
5174f25acd9SJames Wright     PetscCall(KSPSetTolerances(ksp, .99, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT));
51822c1b34eSJames Wright     PetscCall(TSSetSolution(*ts, Q_preload));
51991982731SJeremy L Thompson     PetscCall(TSStep(*ts));
5204f25acd9SJames Wright     PetscCall(SNESSetTolerances(snes, PETSC_CURRENT, rtol_snes, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT));
5214f25acd9SJames Wright     PetscCall(KSPSetTolerances(ksp, rtol_ksp, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT));
52291982731SJeremy L Thompson     PetscCall(VecDestroy(&Q_preload));
52391982731SJeremy L Thompson     // LCOV_EXCL_STOP
52491982731SJeremy L Thompson   } else {
5252b916ea7SJeremy L Thompson     PetscCall(PetscBarrier((PetscObject)*ts));
5262a9a4b51SJames Wright     PetscCall(TSSolve(*ts, Q));
52791982731SJeremy L Thompson   }
52891982731SJeremy L Thompson   PetscPreLoadEnd();
52991982731SJeremy L Thompson 
53091982731SJeremy L Thompson   PetscCall(TSGetSolveTime(*ts, &final_time));
531a515125bSLeila Ghaffari   *f_time = final_time;
53291982731SJeremy L Thompson 
5330e1e9333SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) {
5347538d537SJames Wright     PetscInt           step_no;
5354e1dcb99SJames Wright     PetscLogStage      stage_id;
5364e1dcb99SJames Wright     PetscEventPerfInfo stage_perf;
5374e1dcb99SJames Wright 
5387538d537SJames Wright     PetscCall(TSGetStepNumber(*ts, &step_no));
5390c373b74SJames Wright     if (honee->app_ctx->checkpoint_interval > 0 || honee->app_ctx->checkpoint_interval == -1) {
5402a9a4b51SJames Wright       PetscCall(WriteOutput(honee, Q, step_no, final_time));
5417538d537SJames Wright     }
5424e1dcb99SJames Wright     PetscCall(HoneeTSEvaluationSolutions(*ts));
54391982731SJeremy L Thompson 
544ea615d4cSJames Wright     PetscCall(PetscLogStageGetId("HONEE Solve", &stage_id));
545df4304b5SJed Brown     PetscCall(PetscLogStageGetPerfInfo(stage_id, &stage_perf));
546df4304b5SJed Brown     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time taken for solution (sec): %g\n", stage_perf.time));
547a515125bSLeila Ghaffari   }
548d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
549a515125bSLeila Ghaffari }
550