xref: /libCEED/examples/fluids/src/setupts.c (revision 8b85812c06baf6831ed3ba6d22d2276f4fca36d7)
1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 /// @file
9 /// Time-stepping functions for Navier-Stokes example using PETSc
10 
11 #include "../navierstokes.h"
12 #include "../qfunctions/mass.h"
13 
14 // Compute mass matrix for explicit scheme
15 PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm, CeedData ceed_data,
16                                        Vec M) {
17   Vec            M_loc;
18   CeedQFunction  qf_mass;
19   CeedOperator   op_mass;
20   CeedVector     m_ceed, ones_vec;
21   CeedInt        num_comp_q, q_data_size;
22   PetscErrorCode ierr;
23   PetscFunctionBeginUser;
24 
25   // CEED Restriction
26   CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q);
27   CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &q_data_size);
28   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &m_ceed, NULL);
29   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &ones_vec, NULL);
30   CeedVectorSetValue(ones_vec, 1.0);
31 
32   // CEED QFunction
33   CeedQFunctionCreateInterior(ceed, 1, Mass, Mass_loc, &qf_mass);
34   CeedQFunctionAddInput(qf_mass, "q", num_comp_q, CEED_EVAL_INTERP);
35   CeedQFunctionAddInput(qf_mass, "qdata", q_data_size, CEED_EVAL_NONE);
36   CeedQFunctionAddOutput(qf_mass, "v", num_comp_q, CEED_EVAL_INTERP);
37 
38   // CEED Operator
39   CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass);
40   CeedOperatorSetField(op_mass, "q", ceed_data->elem_restr_q, ceed_data->basis_q,
41                        CEED_VECTOR_ACTIVE);
42   CeedOperatorSetField(op_mass, "qdata", ceed_data->elem_restr_qd_i,
43                        CEED_BASIS_COLLOCATED, ceed_data->q_data);
44   CeedOperatorSetField(op_mass, "v", ceed_data->elem_restr_q, ceed_data->basis_q,
45                        CEED_VECTOR_ACTIVE);
46 
47   // Place PETSc vector in CEED vector
48   CeedScalar *m;
49   PetscMemType m_mem_type;
50   ierr = DMGetLocalVector(dm, &M_loc); CHKERRQ(ierr);
51   ierr = VecGetArrayAndMemType(M_loc, (PetscScalar **)&m, &m_mem_type);
52   CHKERRQ(ierr);
53   CeedVectorSetArray(m_ceed, MemTypeP2C(m_mem_type), CEED_USE_POINTER, m);
54 
55   // Apply CEED Operator
56   CeedOperatorApply(op_mass, ones_vec, m_ceed, CEED_REQUEST_IMMEDIATE);
57 
58   // Restore vectors
59   CeedVectorTakeArray(m_ceed, MemTypeP2C(m_mem_type), NULL);
60   ierr = VecRestoreArrayReadAndMemType(M_loc, (const PetscScalar **)&m);
61   CHKERRQ(ierr);
62 
63   // Local-to-Global
64   ierr = VecZeroEntries(M); CHKERRQ(ierr);
65   ierr = DMLocalToGlobal(dm, M_loc, ADD_VALUES, M); CHKERRQ(ierr);
66   ierr = DMRestoreLocalVector(dm, &M_loc); CHKERRQ(ierr);
67 
68   // Invert diagonally lumped mass vector for RHS function
69   ierr = VecReciprocal(M); CHKERRQ(ierr);
70 
71   // Cleanup
72   CeedVectorDestroy(&ones_vec);
73   CeedVectorDestroy(&m_ceed);
74   CeedQFunctionDestroy(&qf_mass);
75   CeedOperatorDestroy(&op_mass);
76 
77   PetscFunctionReturn(0);
78 }
79 
80 // RHS (Explicit time-stepper) function setup
81 //   This is the RHS of the ODE, given as u_t = G(t,u)
82 //   This function takes in a state vector Q and writes into G
83 PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) {
84   User           user = *(User *)user_data;
85   PetscScalar    *q, *g;
86   Vec            Q_loc = user->Q_loc, G_loc;
87   PetscMemType   q_mem_type, g_mem_type;
88   PetscErrorCode ierr;
89   PetscFunctionBeginUser;
90 
91   // Get local vector
92   ierr = DMGetLocalVector(user->dm, &G_loc); CHKERRQ(ierr);
93 
94   // Update time dependent data
95   if (user->time != t) {
96     ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t,
97                                       NULL, NULL, NULL); CHKERRQ(ierr);
98     if (user->phys->solution_time_label) {
99       CeedOperatorContextSetDouble(user->op_rhs, user->phys->solution_time_label, &t);
100     }
101     user->time = t;
102   }
103   if (user->phys->timestep_size_label) {
104     PetscScalar dt;
105     ierr = TSGetTimeStep(ts, &dt); CHKERRQ(ierr);
106     if (user->dt != dt) {
107       CeedOperatorContextSetDouble(user->op_rhs, user->phys->timestep_size_label,
108                                    &dt);
109       user->dt = dt;
110     }
111   }
112 
113   // Global-to-local
114   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc); CHKERRQ(ierr);
115 
116   // Place PETSc vectors in CEED vectors
117   ierr = VecGetArrayReadAndMemType(Q_loc, (const PetscScalar **)&q, &q_mem_type);
118   CHKERRQ(ierr);
119   ierr = VecGetArrayAndMemType(G_loc, &g, &g_mem_type); CHKERRQ(ierr);
120   CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER, q);
121   CeedVectorSetArray(user->g_ceed, MemTypeP2C(g_mem_type), CEED_USE_POINTER, g);
122 
123   // Apply CEED operator
124   CeedOperatorApply(user->op_rhs, user->q_ceed, user->g_ceed,
125                     CEED_REQUEST_IMMEDIATE);
126 
127   // Restore vectors
128   CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL);
129   CeedVectorTakeArray(user->g_ceed, MemTypeP2C(g_mem_type), NULL);
130   ierr = VecRestoreArrayReadAndMemType(Q_loc, (const PetscScalar **)&q);
131   CHKERRQ(ierr);
132   ierr = VecRestoreArrayAndMemType(G_loc, &g); CHKERRQ(ierr);
133 
134   // Local-to-Global
135   ierr = VecZeroEntries(G); CHKERRQ(ierr);
136   ierr = DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G); CHKERRQ(ierr);
137 
138   // Inverse of the lumped mass matrix (M is Minv)
139   ierr = VecPointwiseMult(G, G, user->M); CHKERRQ(ierr);
140 
141   // Restore vectors
142   ierr = DMRestoreLocalVector(user->dm, &G_loc); CHKERRQ(ierr);
143 
144   PetscFunctionReturn(0);
145 }
146 
147 // Implicit time-stepper function setup
148 PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G,
149                             void *user_data) {
150   User              user = *(User *)user_data;
151   const PetscScalar *q, *q_dot;
152   PetscScalar       *g;
153   Vec               Q_loc = user->Q_loc, Q_dot_loc = user->Q_dot_loc, G_loc;
154   PetscMemType      q_mem_type, q_dot_mem_type, g_mem_type;
155   PetscErrorCode    ierr;
156   PetscFunctionBeginUser;
157 
158   // Get local vectors
159   ierr = DMGetLocalVector(user->dm, &G_loc); CHKERRQ(ierr);
160 
161   // Update time dependent data
162   if (user->time != t) {
163     ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t,
164                                       NULL, NULL, NULL); CHKERRQ(ierr);
165     if (user->phys->solution_time_label) {
166       CeedOperatorContextSetDouble(user->op_ifunction,
167                                    user->phys->solution_time_label, &t);
168     }
169     user->time = t;
170   }
171   if (user->phys->timestep_size_label) {
172     PetscScalar dt;
173     ierr = TSGetTimeStep(ts, &dt); CHKERRQ(ierr);
174     if (user->dt != dt) {
175       CeedOperatorContextSetDouble(user->op_ifunction,
176                                    user->phys->timestep_size_label, &dt);
177       user->dt = dt;
178     }
179   }
180 
181   // Global-to-local
182   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc); CHKERRQ(ierr);
183   ierr = DMGlobalToLocal(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc);
184   CHKERRQ(ierr);
185 
186   // Place PETSc vectors in CEED vectors
187   ierr = VecGetArrayReadAndMemType(Q_loc, &q, &q_mem_type); CHKERRQ(ierr);
188   ierr = VecGetArrayReadAndMemType(Q_dot_loc, &q_dot, &q_dot_mem_type);
189   CHKERRQ(ierr);
190   ierr = VecGetArrayAndMemType(G_loc, &g, &g_mem_type); CHKERRQ(ierr);
191   CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER,
192                      (PetscScalar *)q);
193   CeedVectorSetArray(user->q_dot_ceed, MemTypeP2C(q_dot_mem_type),
194                      CEED_USE_POINTER, (PetscScalar *)q_dot);
195   CeedVectorSetArray(user->g_ceed, MemTypeP2C(g_mem_type), CEED_USE_POINTER, g);
196 
197   // Apply CEED operator
198   CeedOperatorApply(user->op_ifunction, user->q_ceed, user->g_ceed,
199                     CEED_REQUEST_IMMEDIATE);
200 
201   // Restore vectors
202   CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL);
203   CeedVectorTakeArray(user->q_dot_ceed, MemTypeP2C(q_dot_mem_type), NULL);
204   CeedVectorTakeArray(user->g_ceed, MemTypeP2C(g_mem_type), NULL);
205   ierr = VecRestoreArrayReadAndMemType(Q_loc, &q); CHKERRQ(ierr);
206   ierr = VecRestoreArrayReadAndMemType(Q_dot_loc, &q_dot); CHKERRQ(ierr);
207   ierr = VecRestoreArrayAndMemType(G_loc, &g); CHKERRQ(ierr);
208 
209   // Local-to-Global
210   ierr = VecZeroEntries(G); CHKERRQ(ierr);
211   ierr = DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G); CHKERRQ(ierr);
212 
213   // Restore vectors
214   ierr = DMRestoreLocalVector(user->dm, &G_loc); CHKERRQ(ierr);
215 
216   PetscFunctionReturn(0);
217 }
218 
219 static PetscErrorCode MatMult_NS_IJacobian(Mat J, Vec Q, Vec G) {
220   User user;
221   const PetscScalar *q;
222   PetscScalar       *g;
223   PetscMemType      q_mem_type, g_mem_type;
224   PetscErrorCode    ierr;
225   PetscFunctionBeginUser;
226   ierr = MatShellGetContext(J, &user); CHKERRQ(ierr);
227   Vec               Q_loc = user->Q_dot_loc, // Note - Q_dot_loc has zero BCs
228                     G_loc;
229 
230   // Get local vectors
231   ierr = DMGetLocalVector(user->dm, &G_loc); CHKERRQ(ierr);
232 
233   // Global-to-local
234   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc); CHKERRQ(ierr);
235 
236   // Place PETSc vectors in CEED vectors
237   ierr = VecGetArrayReadAndMemType(Q_loc, &q, &q_mem_type); CHKERRQ(ierr);
238   ierr = VecGetArrayAndMemType(G_loc, &g, &g_mem_type); CHKERRQ(ierr);
239   CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER,
240                      (PetscScalar *)q);
241   CeedVectorSetArray(user->g_ceed, MemTypeP2C(g_mem_type), CEED_USE_POINTER, g);
242 
243   // Apply CEED operator
244   CeedOperatorApply(user->op_ijacobian, user->q_ceed, user->g_ceed,
245                     CEED_REQUEST_IMMEDIATE);
246 
247   // Restore vectors
248   CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL);
249   CeedVectorTakeArray(user->g_ceed, MemTypeP2C(g_mem_type), NULL);
250   ierr = VecRestoreArrayReadAndMemType(Q_loc, &q); CHKERRQ(ierr);
251   ierr = VecRestoreArrayAndMemType(G_loc, &g); CHKERRQ(ierr);
252 
253   // Local-to-Global
254   ierr = VecZeroEntries(G); CHKERRQ(ierr);
255   ierr = DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G); CHKERRQ(ierr);
256 
257   // Restore vectors
258   ierr = DMRestoreLocalVector(user->dm, &G_loc); CHKERRQ(ierr);
259   PetscFunctionReturn(0);
260 }
261 
262 PetscErrorCode MatGetDiagonal_NS_IJacobian(Mat A, Vec D) {
263   User user;
264   Vec D_loc;
265   PetscScalar *d;
266   PetscMemType mem_type;
267 
268   PetscFunctionBeginUser;
269   MatShellGetContext(A, &user);
270   PetscCall(DMGetLocalVector(user->dm, &D_loc));
271   PetscCall(VecGetArrayAndMemType(D_loc, &d, &mem_type));
272   CeedVectorSetArray(user->g_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, d);
273   CeedOperatorLinearAssembleDiagonal(user->op_ijacobian, user->g_ceed,
274                                      CEED_REQUEST_IMMEDIATE);
275   CeedVectorTakeArray(user->g_ceed, MemTypeP2C(mem_type), NULL);
276   PetscCall(VecRestoreArrayAndMemType(D_loc, &d));
277   PetscCall(VecZeroEntries(D));
278   PetscCall(DMLocalToGlobal(user->dm, D_loc, ADD_VALUES, D));
279   PetscCall(DMRestoreLocalVector(user->dm, &D_loc));
280   VecViewFromOptions(D, NULL, "-diag_vec_view");
281   PetscFunctionReturn(0);
282 }
283 
284 PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot,
285                                 PetscReal shift, Mat J, Mat J_pre,
286                                 void *user_data) {
287   User user = *(User *)user_data;
288   PetscBool J_is_shell, J_pre_is_shell;
289   PetscFunctionBeginUser;
290   if (user->phys->ijacobian_time_shift_label)
291     CeedOperatorContextSetDouble(user->op_ijacobian,
292                                  user->phys->ijacobian_time_shift_label, &shift);
293   PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
294   PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
295   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATSHELL, &J_is_shell));
296   PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATSHELL,
297                                    &J_pre_is_shell));
298   if (!user->matrices_set_up) {
299     if (J_is_shell) {
300       PetscCall(MatShellSetContext(J, user));
301       PetscCall(MatShellSetOperation(J, MATOP_MULT,
302                                      (void (*)(void))MatMult_NS_IJacobian));
303       PetscCall(MatShellSetOperation(J, MATOP_GET_DIAGONAL,
304                                      (void (*)(void))MatGetDiagonal_NS_IJacobian));
305       PetscCall(MatSetUp(J));
306     }
307     if (!J_pre_is_shell) {
308       PetscCount ncoo;
309       PetscInt *rows, *cols;
310       PetscCall(CeedOperatorLinearAssembleSymbolic(user->op_ijacobian, &ncoo, &rows,
311                 &cols));
312       PetscCall(MatSetPreallocationCOOLocal(J_pre, ncoo, rows, cols));
313       free(rows);
314       free(cols);
315       CeedVectorCreate(user->ceed, ncoo, &user->coo_values);
316       user->matrices_set_up = true;
317     }
318   }
319   if (!J_pre_is_shell) {
320     CeedMemType mem_type = CEED_MEM_HOST;
321     const PetscScalar *values;
322     MatType mat_type;
323     PetscCall(MatGetType(J_pre, &mat_type));
324     if (strstr(mat_type, "kokkos") || strstr(mat_type, "cusparse"))
325       mem_type = CEED_MEM_DEVICE;
326     CeedOperatorLinearAssemble(user->op_ijacobian, user->coo_values);
327     CeedVectorGetArrayRead(user->coo_values, mem_type, &values);
328     PetscCall(MatSetValuesCOO(J_pre, values, INSERT_VALUES));
329     CeedVectorRestoreArrayRead(user->coo_values, &values);
330   }
331   PetscFunctionReturn(0);
332 }
333 
334 // User provided TS Monitor
335 PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time,
336                             Vec Q, void *ctx) {
337   User           user = ctx;
338   Vec            Q_loc;
339   char           file_path[PETSC_MAX_PATH_LEN];
340   PetscViewer    viewer;
341   PetscErrorCode ierr;
342   PetscFunctionBeginUser;
343 
344   // Print every 'output_freq' steps
345   if (step_no % user->app_ctx->output_freq != 0)
346     PetscFunctionReturn(0);
347 
348   // Set up output
349   ierr = DMGetLocalVector(user->dm, &Q_loc); CHKERRQ(ierr);
350   ierr = PetscObjectSetName((PetscObject)Q_loc, "StateVec"); CHKERRQ(ierr);
351   ierr = VecZeroEntries(Q_loc); CHKERRQ(ierr);
352   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc); CHKERRQ(ierr);
353 
354   // Output
355   ierr = PetscSNPrintf(file_path, sizeof file_path,
356                        "%s/ns-%03" PetscInt_FMT ".vtu",
357                        user->app_ctx->output_dir, step_no + user->app_ctx->cont_steps);
358   CHKERRQ(ierr);
359   ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path,
360                             FILE_MODE_WRITE, &viewer); CHKERRQ(ierr);
361   ierr = VecView(Q_loc, viewer); CHKERRQ(ierr);
362   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
363   if (user->dm_viz) {
364     Vec         Q_refined, Q_refined_loc;
365     char        file_path_refined[PETSC_MAX_PATH_LEN];
366     PetscViewer viewer_refined;
367 
368     ierr = DMGetGlobalVector(user->dm_viz, &Q_refined); CHKERRQ(ierr);
369     ierr = DMGetLocalVector(user->dm_viz, &Q_refined_loc); CHKERRQ(ierr);
370     ierr = PetscObjectSetName((PetscObject)Q_refined_loc, "Refined");
371     CHKERRQ(ierr);
372     ierr = MatInterpolate(user->interp_viz, Q, Q_refined); CHKERRQ(ierr);
373     ierr = VecZeroEntries(Q_refined_loc); CHKERRQ(ierr);
374     ierr = DMGlobalToLocal(user->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc);
375     CHKERRQ(ierr);
376     ierr = PetscSNPrintf(file_path_refined, sizeof file_path_refined,
377                          "%s/nsrefined-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir,
378                          step_no + user->app_ctx->cont_steps);
379     CHKERRQ(ierr);
380     ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined),
381                               file_path_refined, FILE_MODE_WRITE, &viewer_refined); CHKERRQ(ierr);
382     ierr = VecView(Q_refined_loc, viewer_refined); CHKERRQ(ierr);
383     ierr = DMRestoreLocalVector(user->dm_viz, &Q_refined_loc); CHKERRQ(ierr);
384     ierr = DMRestoreGlobalVector(user->dm_viz, &Q_refined); CHKERRQ(ierr);
385     ierr = PetscViewerDestroy(&viewer_refined); CHKERRQ(ierr);
386   }
387   ierr = DMRestoreLocalVector(user->dm, &Q_loc); CHKERRQ(ierr);
388 
389   // Save data in a binary file for continuation of simulations
390   ierr = PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin",
391                        user->app_ctx->output_dir); CHKERRQ(ierr);
392   ierr = PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, &viewer);
393   CHKERRQ(ierr);
394   ierr = VecView(Q, viewer); CHKERRQ(ierr);
395   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
396 
397   // Save time stamp
398   // Dimensionalize time back
399   time /= user->units->second;
400   ierr = PetscSNPrintf(file_path, sizeof file_path, "%s/ns-time.bin",
401                        user->app_ctx->output_dir); CHKERRQ(ierr);
402   ierr = PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, &viewer);
403   CHKERRQ(ierr);
404   #if PETSC_VERSION_GE(3,13,0)
405   ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL);
406   #else
407   ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL, true);
408   #endif
409   CHKERRQ(ierr);
410   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
411 
412   PetscFunctionReturn(0);
413 }
414 
415 // TS: Create, setup, and solve
416 PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys,
417                           Vec *Q, PetscScalar *f_time, TS *ts) {
418   MPI_Comm       comm = user->comm;
419   TSAdapt        adapt;
420   PetscScalar    final_time;
421   PetscErrorCode ierr;
422   PetscFunctionBeginUser;
423 
424   ierr = TSCreate(comm, ts); CHKERRQ(ierr);
425   ierr = TSSetDM(*ts, dm); CHKERRQ(ierr);
426   if (phys->implicit) {
427     ierr = TSSetType(*ts, TSBDF); CHKERRQ(ierr);
428     if (user->op_ifunction) {
429       ierr = TSSetIFunction(*ts, NULL, IFunction_NS, &user); CHKERRQ(ierr);
430     } else { // Implicit integrators can fall back to using an RHSFunction
431       ierr = TSSetRHSFunction(*ts, NULL, RHS_NS, &user); CHKERRQ(ierr);
432     }
433     if (user->op_ijacobian) {
434       ierr = DMTSSetIJacobian(dm, FormIJacobian_NS, &user); CHKERRQ(ierr);
435     }
436   } else {
437     if (!user->op_rhs) SETERRQ(comm, PETSC_ERR_ARG_NULL,
438                                  "Problem does not provide RHSFunction");
439     ierr = TSSetType(*ts, TSRK); CHKERRQ(ierr);
440     ierr = TSRKSetType(*ts, TSRK5F); CHKERRQ(ierr);
441     ierr = TSSetRHSFunction(*ts, NULL, RHS_NS, &user); CHKERRQ(ierr);
442   }
443   ierr = TSSetMaxTime(*ts, 500. * user->units->second); CHKERRQ(ierr);
444   ierr = TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER); CHKERRQ(ierr);
445   ierr = TSSetTimeStep(*ts, 1.e-2 * user->units->second); CHKERRQ(ierr);
446   if (app_ctx->test_mode) {ierr = TSSetMaxSteps(*ts, 10); CHKERRQ(ierr);}
447   ierr = TSGetAdapt(*ts, &adapt); CHKERRQ(ierr);
448   ierr = TSAdaptSetStepLimits(adapt,
449                               1.e-12 * user->units->second,
450                               1.e2 * user->units->second); CHKERRQ(ierr);
451   ierr = TSSetFromOptions(*ts); CHKERRQ(ierr);
452   user->time = -1.0; // require all BCs and ctx to be updated
453   user->dt   = -1.0;
454   if (!app_ctx->cont_steps) { // print initial condition
455     if (!app_ctx->test_mode) {
456       ierr = TSMonitor_NS(*ts, 0, 0., *Q, user); CHKERRQ(ierr);
457     }
458   } else { // continue from time of last output
459     PetscReal   time;
460     PetscInt    count;
461     PetscViewer viewer;
462     char        file_path[PETSC_MAX_PATH_LEN];
463     ierr = PetscSNPrintf(file_path, sizeof file_path, "%s/ns-time.bin",
464                          app_ctx->output_dir); CHKERRQ(ierr);
465     ierr = PetscViewerBinaryOpen(comm, file_path, FILE_MODE_READ, &viewer);
466     CHKERRQ(ierr);
467     ierr = PetscViewerBinaryRead(viewer, &time, 1, &count, PETSC_REAL);
468     CHKERRQ(ierr);
469     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
470     ierr = TSSetTime(*ts, time * user->units->second); CHKERRQ(ierr);
471   }
472   if (!app_ctx->test_mode) {
473     ierr = TSMonitorSet(*ts, TSMonitor_NS, user, NULL); CHKERRQ(ierr);
474   }
475 
476   // Solve
477   PetscScalar start_time;
478   ierr = TSGetTime(*ts, &start_time); CHKERRQ(ierr);
479 
480   PetscPreLoadBegin(PETSC_FALSE, "Fluids Solve");
481   PetscCall(TSSetTime(*ts, start_time));
482   PetscCall(TSSetStepNumber(*ts, 0));
483   if (PetscPreLoadingOn) {
484     // LCOV_EXCL_START
485     SNES      snes;
486     Vec       Q_preload;
487     PetscReal rtol;
488     PetscCall(VecDuplicate(*Q, &Q_preload));
489     PetscCall(VecCopy(*Q, Q_preload));
490     PetscCall(TSGetSNES(*ts, &snes));
491     PetscCall(SNESGetTolerances(snes, NULL, &rtol, NULL, NULL, NULL));
492     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, .99, PETSC_DEFAULT,
493                                 PETSC_DEFAULT, PETSC_DEFAULT));
494     PetscCall(TSSetSolution(*ts, *Q));
495     PetscCall(TSStep(*ts));
496     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, rtol, PETSC_DEFAULT,
497                                 PETSC_DEFAULT, PETSC_DEFAULT));
498     PetscCall(VecDestroy(&Q_preload));
499     // LCOV_EXCL_STOP
500   } else {
501     ierr = PetscBarrier((PetscObject) *ts); CHKERRQ(ierr);
502     ierr = TSSolve(*ts, *Q); CHKERRQ(ierr);
503   }
504   PetscPreLoadEnd();
505 
506   PetscCall(TSGetSolveTime(*ts, &final_time));
507   *f_time = final_time;
508 
509   if (!app_ctx->test_mode) {
510     PetscLogEvent stage_id;
511     PetscStageLog stage_log;
512 
513     PetscCall(PetscLogStageGetId("Fluids Solve", &stage_id));
514     PetscCall(PetscLogGetStageLog(&stage_log));
515     PetscCall(PetscPrintf(PETSC_COMM_WORLD,
516                           "Time taken for solution (sec): %g\n",
517                           stage_log->stageInfo[stage_id].perfInfo.time));
518   }
519   PetscFunctionReturn(0);
520 }
521