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