xref: /honee/src/setupts.c (revision 71efcdda15017563735c294bfc5f5b5db7281a45)
1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3 
4 /// @file
5 /// Time-stepping functions for HONEE
6 
7 #include <ceed.h>
8 #include <petscdmplex.h>
9 #include <petscts.h>
10 
11 #include <navierstokes.h>
12 #include "../qfunctions/newtonian_state.h"
13 
14 // @brief Insert Boundary values if it's a new time
15 PetscErrorCode UpdateBoundaryValues(Honee honee, Vec Q_loc, PetscReal t) {
16   PetscFunctionBeginUser;
17   if (honee->time_bc_set != t) {
18     PetscCall(DMPlexInsertBoundaryValues(honee->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL));
19     honee->time_bc_set = t;
20   }
21   PetscFunctionReturn(PETSC_SUCCESS);
22 }
23 
24 // RHS (Explicit time-stepper) function setup
25 //   This is the RHS of the ODE, given as u_t = G(t,u)
26 //   This function takes in a state vector Q and writes into G
27 PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) {
28   Honee        honee = *(Honee *)user_data;
29   Ceed         ceed  = honee->ceed;
30   PetscScalar  dt;
31   Vec          Q_loc = honee->Q_loc, R;
32   PetscMemType q_mem_type;
33 
34   PetscFunctionBeginUser;
35   // Update time dependent data
36   PetscCall(UpdateBoundaryValues(honee, Q_loc, t));
37   if (honee->phys->solution_time_label)
38     PetscCallCeed(ceed, CeedOperatorSetContextDouble(honee->op_rhs_ctx->op, honee->phys->solution_time_label, &t));
39   PetscCall(TSGetTimeStep(ts, &dt));
40   if (honee->phys->timestep_size_label)
41     PetscCallCeed(ceed, CeedOperatorSetContextDouble(honee->op_rhs_ctx->op, honee->phys->timestep_size_label, &dt));
42 
43   PetscCall(DMGetNamedGlobalVector(honee->dm, "RHS Residual", &R));
44   PetscCall(DMGlobalToLocal(honee->dm, Q, INSERT_VALUES, Q_loc));
45   if (honee->app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) PetscCall(DivDiffFluxProjectionApply(honee->diff_flux_proj, Q_loc));
46   PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, R, honee->op_rhs_ctx));
47 
48   // Inverse of the mass matrix
49   PetscCall(VecReadPetscToCeed(Q_loc, &q_mem_type, honee->q_ceed));
50   PetscCall(KSPSolve(honee->mass_ksp, R, G));
51   PetscCall(VecReadCeedToPetsc(honee->q_ceed, q_mem_type, Q_loc));
52 
53   PetscCall(DMRestoreNamedGlobalVector(honee->dm, "RHS Residual", &R));
54   PetscFunctionReturn(PETSC_SUCCESS);
55 }
56 
57 // Surface forces function setup
58 static PetscErrorCode Surface_Forces_NS(DM dm, Vec G_loc, PetscInt num_walls, const PetscInt walls[], PetscScalar *reaction_force) {
59   DMLabel            face_label;
60   const PetscScalar *g_array;
61   PetscInt           dim  = 3;
62   MPI_Comm           comm = PetscObjectComm((PetscObject)dm);
63   PetscSection       section;
64 
65   PetscFunctionBeginUser;
66   PetscCall(DMGetLabel(dm, "Face Sets", &face_label));
67   PetscCall(VecGetArrayRead(G_loc, &g_array));
68   for (PetscInt w = 0; w < num_walls; w++) {
69     const PetscInt wall = walls[w], *points;
70     IS             wall_is;
71     PetscInt       num_points, num_comp = 0;
72 
73     PetscCall(DMLabelGetStratumIS(face_label, wall, &wall_is));
74     if (!wall_is) continue;  // No wall points on this process, skip
75 
76     PetscCall(DMGetLocalSection(dm, &section));
77     PetscCall(PetscSectionGetFieldComponents(section, 0, &num_comp));
78     PetscCall(ISGetSize(wall_is, &num_points));
79     PetscCall(ISGetIndices(wall_is, &points));
80     for (PetscInt i = 0; i < num_points; i++) {
81       const PetscInt           p = points[i];
82       const StateConservative *r;
83       PetscInt                 dof;
84 
85       PetscCall(DMPlexPointLocalRead(dm, p, g_array, &r));
86       PetscCall(PetscSectionGetDof(section, p, &dof));
87       for (PetscInt node = 0; node < dof / num_comp; node++) {
88         for (PetscInt j = 0; j < dim; j++) {
89           reaction_force[w * dim + j] -= r[node].momentum[j];
90         }
91       }
92     }
93     PetscCall(ISRestoreIndices(wall_is, &points));
94     PetscCall(ISDestroy(&wall_is));
95   }
96   PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, reaction_force, dim * num_walls, MPIU_SCALAR, MPI_SUM, comm));
97   PetscCall(VecRestoreArrayRead(G_loc, &g_array));
98   PetscFunctionReturn(PETSC_SUCCESS);
99 }
100 
101 // Implicit time-stepper function setup
102 PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) {
103   Honee        honee = *(Honee *)user_data;
104   Ceed         ceed  = honee->ceed;
105   PetscScalar  dt;
106   Vec          Q_loc = honee->Q_loc, Q_dot_loc = honee->Q_dot_loc, G_loc;
107   PetscMemType q_mem_type, q_dot_mem_type, g_mem_type;
108 
109   PetscFunctionBeginUser;
110   PetscCall(DMGlobalToLocalBegin(honee->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
111   PetscCall(DMGetNamedLocalVector(honee->dm, "ResidualLocal", &G_loc));
112 
113   // Update time dependent data
114   PetscCall(UpdateBoundaryValues(honee, Q_loc, t));
115   if (honee->phys->solution_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(honee->op_ifunction, honee->phys->solution_time_label, &t));
116   PetscCall(TSGetTimeStep(ts, &dt));
117   if (honee->phys->timestep_size_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(honee->op_ifunction, honee->phys->timestep_size_label, &dt));
118 
119   // Global-to-local
120   PetscCall(DMGlobalToLocalBegin(honee->dm, Q, INSERT_VALUES, Q_loc));
121   PetscCall(DMGlobalToLocalEnd(honee->dm, Q, INSERT_VALUES, Q_loc));
122   if (honee->app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) PetscCall(DivDiffFluxProjectionApply(honee->diff_flux_proj, Q_loc));
123   PetscCall(DMGlobalToLocalEnd(honee->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
124 
125   // Place PETSc vectors in CEED vectors
126   PetscCall(VecReadPetscToCeed(Q_loc, &q_mem_type, honee->q_ceed));
127   PetscCall(VecReadPetscToCeed(Q_dot_loc, &q_dot_mem_type, honee->q_dot_ceed));
128   PetscCall(VecPetscToCeed(G_loc, &g_mem_type, honee->g_ceed));
129 
130   // Apply CEED operator
131   PetscCall(PetscLogEventBegin(HONEE_CeedOperatorApply, Q, G, 0, 0));
132   PetscCall(PetscLogGpuTimeBegin());
133   PetscCallCeed(honee->ceed, CeedOperatorApply(honee->op_ifunction, honee->q_ceed, honee->g_ceed, CEED_REQUEST_IMMEDIATE));
134   PetscCall(PetscLogGpuTimeEnd());
135   PetscCall(PetscLogEventEnd(HONEE_CeedOperatorApply, Q, G, 0, 0));
136 
137   // Restore vectors
138   PetscCall(VecReadCeedToPetsc(honee->q_ceed, q_mem_type, Q_loc));
139   PetscCall(VecReadCeedToPetsc(honee->q_dot_ceed, q_dot_mem_type, Q_dot_loc));
140   PetscCall(VecCeedToPetsc(honee->g_ceed, g_mem_type, G_loc));
141 
142   if (honee->app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) {
143     PetscCall(SgsDDApplyIFunction(honee, Q_loc, G_loc));
144   }
145 
146   // Local-to-Global
147   PetscCall(VecZeroEntries(G));
148   PetscCall(DMLocalToGlobal(honee->dm, G_loc, ADD_VALUES, G));
149 
150   // Restore vectors
151   PetscCall(DMRestoreNamedLocalVector(honee->dm, "ResidualLocal", &G_loc));
152   PetscFunctionReturn(PETSC_SUCCESS);
153 }
154 
155 PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, PetscReal shift, Mat J, Mat J_pre, void *user_data) {
156   Honee     honee = *(Honee *)user_data;
157   PetscBool J_is_matceed, J_is_mffd, J_pre_is_matceed, J_pre_is_mffd;
158 
159   PetscFunctionBeginUser;
160   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATMFFD, &J_is_mffd));
161   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATCEED, &J_is_matceed));
162   PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATMFFD, &J_pre_is_mffd));
163   PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATCEED, &J_pre_is_matceed));
164 
165   PetscCall(MatCeedSetContextReal(honee->mat_ijacobian, "ijacobian time shift", shift));
166 
167   if (J_is_matceed || J_is_mffd) {
168     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
169     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
170   } else PetscCall(MatCeedAssembleCOO(honee->mat_ijacobian, J));
171 
172   if (J_pre_is_matceed && J != J_pre) {
173     PetscCall(MatAssemblyBegin(J_pre, MAT_FINAL_ASSEMBLY));
174     PetscCall(MatAssemblyEnd(J_pre, MAT_FINAL_ASSEMBLY));
175   } else if (!J_pre_is_matceed && !J_pre_is_mffd && J != J_pre) {
176     PetscCall(MatCeedAssembleCOO(honee->mat_ijacobian, J_pre));
177   }
178   PetscFunctionReturn(PETSC_SUCCESS);
179 }
180 
181 PetscErrorCode WriteOutput(Honee honee, Vec Q, PetscInt step_no, PetscScalar time) {
182   Vec         Q_loc;
183   char        file_path[PETSC_MAX_PATH_LEN];
184   PetscViewer viewer;
185 
186   PetscFunctionBeginUser;
187   if (honee->app_ctx->checkpoint_vtk) {
188     // Set up output
189     PetscCall(DMGetLocalVector(honee->dm, &Q_loc));
190     PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec"));
191     PetscCall(VecZeroEntries(Q_loc));
192     PetscCall(DMGlobalToLocal(honee->dm, Q, INSERT_VALUES, Q_loc));
193 
194     // Output
195     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-%03" PetscInt_FMT ".vtu", honee->app_ctx->output_dir, step_no));
196 
197     PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path, FILE_MODE_WRITE, &viewer));
198     PetscCall(VecView(Q_loc, viewer));
199     PetscCall(PetscViewerDestroy(&viewer));
200     if (honee->dm_viz) {
201       Vec         Q_refined, Q_refined_loc;
202       char        file_path_refined[PETSC_MAX_PATH_LEN];
203       PetscViewer viewer_refined;
204 
205       PetscCall(DMGetGlobalVector(honee->dm_viz, &Q_refined));
206       PetscCall(DMGetLocalVector(honee->dm_viz, &Q_refined_loc));
207       PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined"));
208 
209       PetscCall(MatInterpolate(honee->interp_viz, Q, Q_refined));
210       PetscCall(VecZeroEntries(Q_refined_loc));
211       PetscCall(DMGlobalToLocal(honee->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc));
212 
213       PetscCall(
214           PetscSNPrintf(file_path_refined, sizeof file_path_refined, "%s/nsrefined-%03" PetscInt_FMT ".vtu", honee->app_ctx->output_dir, step_no));
215 
216       PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined), file_path_refined, FILE_MODE_WRITE, &viewer_refined));
217       PetscCall(VecView(Q_refined_loc, viewer_refined));
218       PetscCall(DMRestoreLocalVector(honee->dm_viz, &Q_refined_loc));
219       PetscCall(DMRestoreGlobalVector(honee->dm_viz, &Q_refined));
220       PetscCall(PetscViewerDestroy(&viewer_refined));
221     }
222     PetscCall(DMRestoreLocalVector(honee->dm, &Q_loc));
223   }
224 
225   // Save data in a binary file for continuation of simulations
226   if (honee->app_ctx->add_stepnum2bin) {
227     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", honee->app_ctx->output_dir, step_no));
228   } else {
229     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", honee->app_ctx->output_dir));
230   }
231   PetscCall(PetscViewerBinaryOpen(honee->comm, file_path, FILE_MODE_WRITE, &viewer));
232 
233   time /= honee->units->second;  // Dimensionalize time back
234   PetscCall(HoneeWriteBinaryVec(viewer, Q, time, step_no));
235   PetscCall(PetscViewerDestroy(&viewer));
236   PetscFunctionReturn(PETSC_SUCCESS);
237 }
238 
239 // CSV Monitor
240 PetscErrorCode TSMonitor_WallForce(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
241   Honee             honee = ctx;
242   Vec               G_loc;
243   PetscInt          num_wall = honee->app_ctx->wall_forces.num_wall, dim = 3;
244   const PetscInt   *walls  = honee->app_ctx->wall_forces.walls;
245   PetscViewer       viewer = honee->app_ctx->wall_forces.viewer;
246   PetscViewerFormat format = honee->app_ctx->wall_forces.viewer_format;
247   PetscScalar      *reaction_force;
248   PetscBool         is_ascii;
249 
250   PetscFunctionBeginUser;
251   if (!viewer) PetscFunctionReturn(PETSC_SUCCESS);
252   PetscCall(DMGetNamedLocalVector(honee->dm, "ResidualLocal", &G_loc));
253   PetscCall(PetscCalloc1(num_wall * dim, &reaction_force));
254   PetscCall(Surface_Forces_NS(honee->dm, G_loc, num_wall, walls, reaction_force));
255   PetscCall(DMRestoreNamedLocalVector(honee->dm, "ResidualLocal", &G_loc));
256 
257   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &is_ascii));
258 
259   if (is_ascii) {
260     if (format == PETSC_VIEWER_ASCII_CSV && !honee->app_ctx->wall_forces.header_written) {
261       PetscCall(PetscViewerASCIIPrintf(viewer, "Step,Time,Wall,ForceX,ForceY,ForceZ\n"));
262       honee->app_ctx->wall_forces.header_written = PETSC_TRUE;
263     }
264     for (PetscInt w = 0; w < num_wall; w++) {
265       PetscInt wall = walls[w];
266       if (format == PETSC_VIEWER_ASCII_CSV) {
267         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT ",%g,%" PetscInt_FMT ",%g,%g,%g\n", step_no, time, wall,
268                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
269 
270       } else {
271         PetscCall(PetscViewerASCIIPrintf(viewer, "Wall %" PetscInt_FMT " Forces: Force_x = %12g, Force_y = %12g, Force_z = %12g\n", wall,
272                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
273       }
274     }
275   }
276   PetscCall(PetscFree(reaction_force));
277   PetscFunctionReturn(PETSC_SUCCESS);
278 }
279 
280 // User provided TS Monitor
281 PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
282   Honee honee = ctx;
283 
284   PetscFunctionBeginUser;
285   // Print every 'checkpoint_interval' steps
286   if (honee->app_ctx->checkpoint_interval <= 0 || step_no % honee->app_ctx->checkpoint_interval != 0 ||
287       (honee->app_ctx->cont_steps == step_no && step_no != 0)) {
288     PetscFunctionReturn(PETSC_SUCCESS);
289   }
290 
291   PetscCall(WriteOutput(honee, Q, step_no, time));
292   PetscFunctionReturn(PETSC_SUCCESS);
293 }
294 
295 PetscErrorCode TSPostStep_CheckStep(TS ts) {
296   Honee     honee;
297   PetscReal norm;
298   PetscInt  step;
299   Vec       Q;
300 
301   PetscFunctionBeginUser;
302   PetscCall(TSGetApplicationContext(ts, &honee));
303   PetscCall(TSGetStepNumber(ts, &step));
304   PetscCall(TSGetSolution(ts, &Q));
305   if (step % honee->app_ctx->check_step_interval) PetscFunctionReturn(PETSC_SUCCESS);
306   PetscCall(VecNorm(Q, NORM_1, &norm));
307   if (PetscIsInfOrNanReal(norm)) {
308     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "Solution diverged: Nans found in solution\n"));
309     PetscCall(TSSetConvergedReason(ts, TS_DIVERGED_NONLINEAR_SOLVE));
310   }
311   PetscFunctionReturn(PETSC_SUCCESS);
312 }
313 
314 PetscErrorCode TSPostStep_MaxWallTime(TS ts) {
315   Honee       honee;
316   PetscInt    step;
317   PetscMPIInt rank;
318   MPI_Comm    comm;
319   PetscBool   is_wall_time_exceeded = PETSC_FALSE;
320 
321   PetscFunctionBeginUser;
322   PetscCall(TSGetApplicationContext(ts, &honee));
323   PetscCall(TSGetStepNumber(ts, &step));
324   if (step % honee->max_wall_time_interval) PetscFunctionReturn(PETSC_SUCCESS);
325   PetscCall(PetscObjectGetComm((PetscObject)ts, &comm));
326   PetscCallMPI(MPI_Comm_rank(comm, &rank));
327   if (rank == 0) is_wall_time_exceeded = time(NULL) > honee->max_wall_time ? PETSC_TRUE : PETSC_FALSE;
328   // Must broadcast to avoid race condition
329   PetscCallMPI(MPI_Bcast(&is_wall_time_exceeded, 1, MPIU_BOOL, 0, comm));
330   if (PetscUnlikely(is_wall_time_exceeded)) {
331     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "Stopping TSSolve: Set max wall time exceeded\n"));
332     PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_USER));
333   }
334   PetscFunctionReturn(PETSC_SUCCESS);
335 }
336 
337 /**
338   @brief TSPostStep for HONEE
339 
340   `TSSetPostStep()` only accepts a single function argument, so this function groups together all post-step
341   functionality needed for HONEE features
342 
343   @param[in] ts TS object
344 **/
345 PetscErrorCode TSPostStep_Honee(TS ts) {
346   Honee honee;
347 
348   PetscFunctionBeginUser;
349   PetscCall(TSGetApplicationContext(ts, &honee));
350   if (honee->max_wall_time > 0) PetscCall(TSPostStep_MaxWallTime(ts));
351   if (honee->app_ctx->sgs_train_enable) PetscCall(TSPostStep_SGS_DD_Training(ts));
352   if (honee->app_ctx->check_step_interval > 0) PetscCall(TSPostStep_CheckStep(ts));
353   PetscFunctionReturn(PETSC_SUCCESS);
354 }
355 
356 // TS: Create, setup, and solve
357 PetscErrorCode TSSolve_NS(DM dm, Honee honee, AppCtx app_ctx, Physics phys, ProblemData problem, Vec *Q, PetscScalar *f_time, TS *ts) {
358   MPI_Comm    comm = honee->comm;
359   TSAdapt     adapt;
360   PetscScalar final_time;
361 
362   PetscFunctionBeginUser;
363   PetscCall(TSCreate(comm, ts));
364   PetscCall(TSSetDM(*ts, dm));
365   PetscCall(TSSetApplicationContext(*ts, honee));
366   if (phys->implicit) {
367     PetscCall(TSSetType(*ts, TSBDF));
368     if (honee->op_ifunction) {
369       PetscCall(TSSetIFunction(*ts, NULL, IFunction_NS, &honee));
370     } else {  // Implicit integrators can fall back to using an RHSFunction
371       PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &honee));
372     }
373     if (honee->mat_ijacobian) {
374       PetscCall(DMTSSetIJacobian(dm, FormIJacobian_NS, &honee));
375     }
376   } else {
377     PetscCheck(honee->op_rhs_ctx, comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction");
378     PetscCall(TSSetType(*ts, TSRK));
379     PetscCall(TSRKSetType(*ts, TSRK5F));
380     PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &honee));
381   }
382   PetscCall(TSSetMaxTime(*ts, 500. * honee->units->second));
383   PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER));
384   if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(TSSetErrorIfStepFails(*ts, PETSC_FALSE));
385   PetscCall(TSSetTimeStep(*ts, 1.e-2 * honee->units->second));
386   PetscCall(TSGetAdapt(*ts, &adapt));
387   PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * honee->units->second, 1.e2 * honee->units->second));
388   PetscCall(TSSetFromOptions(*ts));
389   if (honee->mat_ijacobian) {
390     if (app_ctx->amat_type && !strcmp(app_ctx->amat_type, MATSHELL)) {
391       SNES snes;
392       KSP  ksp;
393       Mat  Pmat, Amat;
394 
395       PetscCall(TSGetSNES(*ts, &snes));
396       PetscCall(SNESGetKSP(snes, &ksp));
397       PetscCall(CreateSolveOperatorsFromMatCeed(ksp, honee->mat_ijacobian, PETSC_FALSE, &Amat, &Pmat));
398       PetscCall(TSSetIJacobian(*ts, honee->mat_ijacobian, Pmat, NULL, NULL));
399       PetscCall(MatDestroy(&Amat));
400       PetscCall(MatDestroy(&Pmat));
401     }
402   }
403   honee->time_bc_set = -1.0;  // require all BCs be updated
404   if (app_ctx->cont_steps) {  // continue from previous timestep data
405     PetscCall(TSSetTime(*ts, app_ctx->cont_time * honee->units->second));
406     PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps));
407   }
408   if (honee->set_poststep) PetscCall(TSSetPostStep(*ts, TSPostStep_Honee));
409   if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(TSMonitorSet(*ts, TSMonitor_NS, honee, NULL));
410   if (app_ctx->wall_forces.viewer) PetscCall(TSMonitorSet(*ts, TSMonitor_WallForce, honee, NULL));
411   if (app_ctx->turb_spanstats_enable) {
412     PetscCall(TSMonitorSet(*ts, TSMonitor_TurbulenceStatistics, honee, NULL));
413     CeedScalar previous_time = app_ctx->cont_time * honee->units->second;
414     PetscCallCeed(honee->ceed,
415                   CeedOperatorSetContextDouble(honee->spanstats.op_stats_collect_ctx->op, honee->spanstats.previous_time_label, &previous_time));
416   }
417   PetscCall(TSMonitorSetFromOptions(*ts, "-ts_monitor_total_kinetic_energy", "Monitor total kinetic energy balance terms in the domain", NULL,
418                                     TSMonitor_TotalKineticEnergy, SetupMontiorTotalKineticEnergy));
419   PetscCall(TSMonitorSetFromOptions(*ts, "-ts_monitor_cfl", "Monitor element CFL statistics", NULL, TSMonitor_Cfl, SetupMontiorCfl));
420   if (app_ctx->diff_filter_monitor) PetscCall(TSMonitorSet(*ts, TSMonitor_DifferentialFilter, honee, NULL));
421   if (app_ctx->sgs_train_enable) PetscCall(TSMonitorSet(*ts, TSMonitor_SGS_DD_Training, honee, NULL));
422 
423   if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(PrintRunInfo(honee, honee->phys, problem, *ts));
424   // Solve
425   PetscReal start_time;
426   PetscInt  start_step;
427   PetscCall(TSGetTime(*ts, &start_time));
428   PetscCall(TSGetStepNumber(*ts, &start_step));
429 
430   PetscCall(PetscLogDefaultBegin());  // So we can use PetscLogStageGetPerfInfo without -log_view
431   PetscPreLoadBegin(PETSC_FALSE, "HONEE Solve");
432   PetscCall(TSSetTime(*ts, start_time));
433   PetscCall(TSSetStepNumber(*ts, start_step));
434   if (PetscPreLoadingOn) {
435     // LCOV_EXCL_START
436     SNES      snes;
437     KSP       ksp;
438     Vec       Q_preload;
439     PetscReal rtol_snes, rtol_ksp;
440     PetscCall(VecDuplicate(*Q, &Q_preload));
441     PetscCall(VecCopy(*Q, Q_preload));
442     PetscCall(TSGetSNES(*ts, &snes));
443     PetscCall(SNESGetTolerances(snes, NULL, &rtol_snes, NULL, NULL, NULL));
444     PetscCall(SNESSetTolerances(snes, PETSC_CURRENT, .99, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT));
445     PetscCall(SNESGetKSP(snes, &ksp));
446     PetscCall(KSPGetTolerances(ksp, &rtol_ksp, NULL, NULL, NULL));
447     PetscCall(KSPSetTolerances(ksp, .99, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT));
448     PetscCall(TSSetSolution(*ts, Q_preload));
449     PetscCall(TSStep(*ts));
450     PetscCall(SNESSetTolerances(snes, PETSC_CURRENT, rtol_snes, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT));
451     PetscCall(KSPSetTolerances(ksp, rtol_ksp, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT));
452     PetscCall(VecDestroy(&Q_preload));
453     // LCOV_EXCL_STOP
454   } else {
455     PetscCall(PetscBarrier((PetscObject)*ts));
456     PetscCall(TSSolve(*ts, *Q));
457   }
458   PetscPreLoadEnd();
459 
460   PetscCall(TSGetSolveTime(*ts, &final_time));
461   *f_time = final_time;
462 
463   if (app_ctx->test_type == TESTTYPE_NONE) {
464     PetscInt step_no;
465     PetscCall(TSGetStepNumber(*ts, &step_no));
466     if (honee->app_ctx->checkpoint_interval > 0 || honee->app_ctx->checkpoint_interval == -1) {
467       PetscCall(WriteOutput(honee, *Q, step_no, final_time));
468     }
469 
470     PetscLogStage      stage_id;
471     PetscEventPerfInfo stage_perf;
472 
473     PetscCall(PetscLogStageGetId("HONEE Solve", &stage_id));
474     PetscCall(PetscLogStageGetPerfInfo(stage_id, &stage_perf));
475     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time taken for solution (sec): %g\n", stage_perf.time));
476   }
477   PetscFunctionReturn(PETSC_SUCCESS);
478 }
479