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