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