xref: /honee/src/setupts.c (revision 91f639d21ae2340098eda83a53dbe75d172c5677)
1727da7e7SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2727da7e7SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3a515125bSLeila Ghaffari //
4727da7e7SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5a515125bSLeila Ghaffari //
6727da7e7SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
7a515125bSLeila Ghaffari 
8a515125bSLeila Ghaffari /// @file
9a515125bSLeila Ghaffari /// Time-stepping functions for Navier-Stokes example using PETSc
10a515125bSLeila Ghaffari 
11a515125bSLeila Ghaffari #include "../navierstokes.h"
12c5e9980aSAdeleke O. Bankole #include "../qfunctions/newtonian_state.h"
13a515125bSLeila Ghaffari 
14a515125bSLeila Ghaffari // Compute mass matrix for explicit scheme
152b916ea7SJeremy L Thompson PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm, CeedData ceed_data, Vec M) {
16a515125bSLeila Ghaffari   Vec           M_loc;
17a515125bSLeila Ghaffari   CeedQFunction qf_mass;
18a515125bSLeila Ghaffari   CeedOperator  op_mass;
19a515125bSLeila Ghaffari   CeedVector    m_ceed, ones_vec;
20a515125bSLeila Ghaffari   CeedInt       num_comp_q, q_data_size;
21a515125bSLeila Ghaffari   PetscFunctionBeginUser;
22a515125bSLeila Ghaffari 
23a515125bSLeila Ghaffari   // CEED Restriction
24a515125bSLeila Ghaffari   CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q);
25a515125bSLeila Ghaffari   CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &q_data_size);
26a515125bSLeila Ghaffari   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &m_ceed, NULL);
27a515125bSLeila Ghaffari   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &ones_vec, NULL);
28a515125bSLeila Ghaffari   CeedVectorSetValue(ones_vec, 1.0);
29a515125bSLeila Ghaffari 
30a515125bSLeila Ghaffari   // CEED QFunction
319f59f36eSJames Wright   PetscCall(CreateMassQFunction(ceed, num_comp_q, q_data_size, &qf_mass));
32a515125bSLeila Ghaffari 
33a515125bSLeila Ghaffari   // CEED Operator
34a515125bSLeila Ghaffari   CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass);
359f59f36eSJames Wright   CeedOperatorSetField(op_mass, "u", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
362b916ea7SJeremy L Thompson   CeedOperatorSetField(op_mass, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data);
372b916ea7SJeremy L Thompson   CeedOperatorSetField(op_mass, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
38a515125bSLeila Ghaffari 
39a515125bSLeila Ghaffari   // Place PETSc vector in CEED vector
40a515125bSLeila Ghaffari   PetscMemType m_mem_type;
412b916ea7SJeremy L Thompson   PetscCall(DMGetLocalVector(dm, &M_loc));
42fd969b44SJames Wright   PetscCall(VecP2C(M_loc, &m_mem_type, m_ceed));
43a515125bSLeila Ghaffari 
44a515125bSLeila Ghaffari   // Apply CEED Operator
45a515125bSLeila Ghaffari   CeedOperatorApply(op_mass, ones_vec, m_ceed, CEED_REQUEST_IMMEDIATE);
46a515125bSLeila Ghaffari 
47a515125bSLeila Ghaffari   // Restore vectors
48fd969b44SJames Wright   PetscCall(VecC2P(m_ceed, m_mem_type, M_loc));
49a515125bSLeila Ghaffari 
50a515125bSLeila Ghaffari   // Local-to-Global
512b916ea7SJeremy L Thompson   PetscCall(VecZeroEntries(M));
522b916ea7SJeremy L Thompson   PetscCall(DMLocalToGlobal(dm, M_loc, ADD_VALUES, M));
532b916ea7SJeremy L Thompson   PetscCall(DMRestoreLocalVector(dm, &M_loc));
54a515125bSLeila Ghaffari 
55a515125bSLeila Ghaffari   // Invert diagonally lumped mass vector for RHS function
562b916ea7SJeremy L Thompson   PetscCall(VecReciprocal(M));
57a515125bSLeila Ghaffari 
58a515125bSLeila Ghaffari   // Cleanup
59a515125bSLeila Ghaffari   CeedVectorDestroy(&ones_vec);
60a515125bSLeila Ghaffari   CeedVectorDestroy(&m_ceed);
61a515125bSLeila Ghaffari   CeedQFunctionDestroy(&qf_mass);
62a515125bSLeila Ghaffari   CeedOperatorDestroy(&op_mass);
63a515125bSLeila Ghaffari 
64a515125bSLeila Ghaffari   PetscFunctionReturn(0);
65a515125bSLeila Ghaffari }
66a515125bSLeila Ghaffari 
67c996854bSJames Wright // Insert Boundary values if it's a new time
68c996854bSJames Wright PetscErrorCode UpdateBoundaryValues(User user, Vec Q_loc, PetscReal t) {
69c996854bSJames Wright   PetscFunctionBeginUser;
70c996854bSJames Wright   if (user->time_bc_set != t) {
71c996854bSJames Wright     PetscCall(DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL));
72c996854bSJames Wright     user->time_bc_set = t;
73c996854bSJames Wright   }
74c996854bSJames Wright   PetscFunctionReturn(0);
75c996854bSJames Wright }
76c996854bSJames Wright 
77*91f639d2SJames Wright // @brief Update the context label value to new value if necessary.
78*91f639d2SJames Wright // @note This only supports labels with scalar label values (ie. not arrays)
79*91f639d2SJames Wright PetscErrorCode UpdateContextLabel(MPI_Comm comm, PetscScalar update_value, CeedOperator op, CeedContextFieldLabel label) {
80*91f639d2SJames Wright   PetscScalar label_value;
81c2cb7fc8SJames Wright 
82*91f639d2SJames Wright   PetscFunctionBeginUser;
83*91f639d2SJames Wright   PetscCheck(label, comm, PETSC_ERR_ARG_BADPTR, "Label should be non-NULL");
84*91f639d2SJames Wright 
85*91f639d2SJames Wright   {
86*91f639d2SJames Wright     size_t             num_elements;
87*91f639d2SJames Wright     const PetscScalar *label_values;
88*91f639d2SJames Wright     CeedOperatorGetContextDoubleRead(op, label, &num_elements, &label_values);
89*91f639d2SJames Wright     PetscCheck(num_elements == 1, comm, PETSC_ERR_SUP, "%s does not support labels with more than 1 value. Label has %zu values", __func__,
90*91f639d2SJames Wright                num_elements);
91*91f639d2SJames Wright     label_value = *label_values;
92*91f639d2SJames Wright     CeedOperatorRestoreContextDoubleRead(op, label, &label_values);
93c2cb7fc8SJames Wright   }
94*91f639d2SJames Wright 
95*91f639d2SJames Wright   if (label_value != update_value) {
96*91f639d2SJames Wright     CeedOperatorSetContextDouble(op, label, &update_value);
97c2cb7fc8SJames Wright   }
98c2cb7fc8SJames Wright   PetscFunctionReturn(0);
99c2cb7fc8SJames Wright }
100c2cb7fc8SJames Wright 
101a515125bSLeila Ghaffari // RHS (Explicit time-stepper) function setup
102a515125bSLeila Ghaffari //   This is the RHS of the ODE, given as u_t = G(t,u)
103a515125bSLeila Ghaffari //   This function takes in a state vector Q and writes into G
104a515125bSLeila Ghaffari PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) {
105a515125bSLeila Ghaffari   User         user = *(User *)user_data;
106*91f639d2SJames Wright   MPI_Comm     comm = PetscObjectComm((PetscObject)ts);
107fd969b44SJames Wright   PetscScalar  dt;
108e2f84137SJeremy L Thompson   Vec          Q_loc = user->Q_loc, G_loc;
109a515125bSLeila Ghaffari   PetscMemType q_mem_type, g_mem_type;
110a515125bSLeila Ghaffari   PetscFunctionBeginUser;
111a515125bSLeila Ghaffari 
112e2f84137SJeremy L Thompson   // Get local vector
1132b916ea7SJeremy L Thompson   PetscCall(DMGetLocalVector(user->dm, &G_loc));
114e2f84137SJeremy L Thompson 
115e2f84137SJeremy L Thompson   // Update time dependent data
116c996854bSJames Wright   PetscCall(UpdateBoundaryValues(user, Q_loc, t));
117*91f639d2SJames Wright   if (user->phys->solution_time_label) PetscCall(UpdateContextLabel(comm, t, user->op_rhs, user->phys->solution_time_label));
1182b916ea7SJeremy L Thompson   PetscCall(TSGetTimeStep(ts, &dt));
119*91f639d2SJames Wright   if (user->phys->timestep_size_label) PetscCall(UpdateContextLabel(comm, dt, user->op_rhs, user->phys->timestep_size_label));
120a515125bSLeila Ghaffari 
121a515125bSLeila Ghaffari   // Global-to-local
1222b916ea7SJeremy L Thompson   PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
123a515125bSLeila Ghaffari 
124a515125bSLeila Ghaffari   // Place PETSc vectors in CEED vectors
125fd969b44SJames Wright   PetscCall(VecReadP2C(Q_loc, &q_mem_type, user->q_ceed));
126fd969b44SJames Wright   PetscCall(VecP2C(G_loc, &g_mem_type, user->g_ceed));
127a515125bSLeila Ghaffari 
128a515125bSLeila Ghaffari   // Apply CEED operator
1292b916ea7SJeremy L Thompson   CeedOperatorApply(user->op_rhs, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE);
130a515125bSLeila Ghaffari 
131a515125bSLeila Ghaffari   // Restore vectors
132fd969b44SJames Wright   PetscCall(VecReadC2P(user->q_ceed, q_mem_type, Q_loc));
133fd969b44SJames Wright   PetscCall(VecC2P(user->g_ceed, g_mem_type, G_loc));
134a515125bSLeila Ghaffari 
135a515125bSLeila Ghaffari   // Local-to-Global
1362b916ea7SJeremy L Thompson   PetscCall(VecZeroEntries(G));
1372b916ea7SJeremy L Thompson   PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G));
138a515125bSLeila Ghaffari 
139a515125bSLeila Ghaffari   // Inverse of the lumped mass matrix (M is Minv)
1402b916ea7SJeremy L Thompson   PetscCall(VecPointwiseMult(G, G, user->M));
141a515125bSLeila Ghaffari 
142a515125bSLeila Ghaffari   // Restore vectors
1432b916ea7SJeremy L Thompson   PetscCall(DMRestoreLocalVector(user->dm, &G_loc));
144a515125bSLeila Ghaffari 
145a515125bSLeila Ghaffari   PetscFunctionReturn(0);
146a515125bSLeila Ghaffari }
147a515125bSLeila Ghaffari 
148c5e9980aSAdeleke O. Bankole // Surface forces function setup
149c5e9980aSAdeleke O. Bankole static PetscErrorCode Surface_Forces_NS(DM dm, Vec G_loc, PetscInt num_walls, const PetscInt walls[], PetscScalar *reaction_force) {
150c5e9980aSAdeleke O. Bankole   DMLabel            face_label;
151c5e9980aSAdeleke O. Bankole   const PetscScalar *g;
1522004e3acSAdeleke O. Bankole   PetscInt           dof, dim = 3;
153c5e9980aSAdeleke O. Bankole   MPI_Comm           comm;
1542004e3acSAdeleke O. Bankole   PetscSection       s;
155c5e9980aSAdeleke O. Bankole 
156c5e9980aSAdeleke O. Bankole   PetscFunctionBeginUser;
157c5e9980aSAdeleke O. Bankole   PetscCall(PetscArrayzero(reaction_force, num_walls * dim));
158c5e9980aSAdeleke O. Bankole   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
159c5e9980aSAdeleke O. Bankole   PetscCall(DMGetLabel(dm, "Face Sets", &face_label));
160c5e9980aSAdeleke O. Bankole   PetscCall(VecGetArrayRead(G_loc, &g));
161c5e9980aSAdeleke O. Bankole   for (PetscInt w = 0; w < num_walls; w++) {
162c5e9980aSAdeleke O. Bankole     const PetscInt wall = walls[w];
163c5e9980aSAdeleke O. Bankole     IS             wall_is;
1642004e3acSAdeleke O. Bankole     PetscCall(DMGetLocalSection(dm, &s));
165c5e9980aSAdeleke O. Bankole     PetscCall(DMLabelGetStratumIS(face_label, wall, &wall_is));
166c5e9980aSAdeleke O. Bankole     if (wall_is) {  // There exist such points on this process
167c5e9980aSAdeleke O. Bankole       PetscInt        num_points;
1682004e3acSAdeleke O. Bankole       PetscInt        num_comp = 0;
169c5e9980aSAdeleke O. Bankole       const PetscInt *points;
1702004e3acSAdeleke O. Bankole       PetscCall(PetscSectionGetFieldComponents(s, 0, &num_comp));
171c5e9980aSAdeleke O. Bankole       PetscCall(ISGetSize(wall_is, &num_points));
172c5e9980aSAdeleke O. Bankole       PetscCall(ISGetIndices(wall_is, &points));
173c5e9980aSAdeleke O. Bankole       for (PetscInt i = 0; i < num_points; i++) {
174c5e9980aSAdeleke O. Bankole         const PetscInt           p = points[i];
175c5e9980aSAdeleke O. Bankole         const StateConservative *r;
176c5e9980aSAdeleke O. Bankole         PetscCall(DMPlexPointLocalRead(dm, p, g, &r));
1772004e3acSAdeleke O. Bankole         PetscCall(PetscSectionGetDof(s, p, &dof));
1782004e3acSAdeleke O. Bankole         for (PetscInt node = 0; node < dof / num_comp; node++) {
179c5e9980aSAdeleke O. Bankole           for (PetscInt j = 0; j < 3; j++) {
1802004e3acSAdeleke O. Bankole             reaction_force[w * dim + j] -= r[node].momentum[j];
1812004e3acSAdeleke O. Bankole           }
182c5e9980aSAdeleke O. Bankole         }
183c5e9980aSAdeleke O. Bankole       }
184c5e9980aSAdeleke O. Bankole       PetscCall(ISRestoreIndices(wall_is, &points));
185c5e9980aSAdeleke O. Bankole     }
186c5e9980aSAdeleke O. Bankole     PetscCall(ISDestroy(&wall_is));
187c5e9980aSAdeleke O. Bankole   }
188c5e9980aSAdeleke O. Bankole   PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, reaction_force, dim * num_walls, MPIU_SCALAR, MPI_SUM, comm));
189c5e9980aSAdeleke O. Bankole   //  Restore Vectors
190c5e9980aSAdeleke O. Bankole   PetscCall(VecRestoreArrayRead(G_loc, &g));
191c5e9980aSAdeleke O. Bankole 
192c5e9980aSAdeleke O. Bankole   PetscFunctionReturn(0);
193c5e9980aSAdeleke O. Bankole }
194c5e9980aSAdeleke O. Bankole 
195a515125bSLeila Ghaffari // Implicit time-stepper function setup
1962b916ea7SJeremy L Thompson PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) {
197a515125bSLeila Ghaffari   User         user = *(User *)user_data;
198*91f639d2SJames Wright   MPI_Comm     comm = PetscObjectComm((PetscObject)ts);
199fd969b44SJames Wright   PetscScalar  dt;
200e2f84137SJeremy L Thompson   Vec          Q_loc = user->Q_loc, Q_dot_loc = user->Q_dot_loc, G_loc;
201a515125bSLeila Ghaffari   PetscMemType q_mem_type, q_dot_mem_type, g_mem_type;
202a515125bSLeila Ghaffari   PetscFunctionBeginUser;
203a515125bSLeila Ghaffari 
204e2f84137SJeremy L Thompson   // Get local vectors
205c5e9980aSAdeleke O. Bankole   PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
206e2f84137SJeremy L Thompson 
207e2f84137SJeremy L Thompson   // Update time dependent data
208c996854bSJames Wright   PetscCall(UpdateBoundaryValues(user, Q_loc, t));
209*91f639d2SJames Wright   if (user->phys->solution_time_label) PetscCall(UpdateContextLabel(comm, t, user->op_ifunction, user->phys->solution_time_label));
2102b916ea7SJeremy L Thompson   PetscCall(TSGetTimeStep(ts, &dt));
211*91f639d2SJames Wright   if (user->phys->timestep_size_label) PetscCall(UpdateContextLabel(comm, dt, user->op_ifunction, user->phys->timestep_size_label));
212a515125bSLeila Ghaffari 
213a515125bSLeila Ghaffari   // Global-to-local
21406108310SJames Wright   PetscCall(DMGlobalToLocalBegin(user->dm, Q, INSERT_VALUES, Q_loc));
21506108310SJames Wright   PetscCall(DMGlobalToLocalBegin(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
21606108310SJames Wright   PetscCall(DMGlobalToLocalEnd(user->dm, Q, INSERT_VALUES, Q_loc));
21706108310SJames Wright   PetscCall(DMGlobalToLocalEnd(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
218a515125bSLeila Ghaffari 
219a515125bSLeila Ghaffari   // Place PETSc vectors in CEED vectors
220fd969b44SJames Wright   PetscCall(VecReadP2C(Q_loc, &q_mem_type, user->q_ceed));
221fd969b44SJames Wright   PetscCall(VecReadP2C(Q_dot_loc, &q_dot_mem_type, user->q_dot_ceed));
222fd969b44SJames Wright   PetscCall(VecP2C(G_loc, &g_mem_type, user->g_ceed));
223a515125bSLeila Ghaffari 
224a515125bSLeila Ghaffari   // Apply CEED operator
2252b916ea7SJeremy L Thompson   CeedOperatorApply(user->op_ifunction, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE);
226a515125bSLeila Ghaffari 
227a515125bSLeila Ghaffari   // Restore vectors
228fd969b44SJames Wright   PetscCall(VecReadC2P(user->q_ceed, q_mem_type, Q_loc));
229fd969b44SJames Wright   PetscCall(VecReadC2P(user->q_dot_ceed, q_dot_mem_type, Q_dot_loc));
230fd969b44SJames Wright   PetscCall(VecC2P(user->g_ceed, g_mem_type, G_loc));
231a515125bSLeila Ghaffari 
232a515125bSLeila Ghaffari   // Local-to-Global
2332b916ea7SJeremy L Thompson   PetscCall(VecZeroEntries(G));
2342b916ea7SJeremy L Thompson   PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G));
235a515125bSLeila Ghaffari 
236a515125bSLeila Ghaffari   // Restore vectors
237c5e9980aSAdeleke O. Bankole   PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
238a515125bSLeila Ghaffari 
239a515125bSLeila Ghaffari   PetscFunctionReturn(0);
240a515125bSLeila Ghaffari }
241a515125bSLeila Ghaffari 
242f0b65372SJed Brown static PetscErrorCode MatMult_NS_IJacobian(Mat J, Vec Q, Vec G) {
243f0b65372SJed Brown   User               user;
244f0b65372SJed Brown   const PetscScalar *q;
245f0b65372SJed Brown   PetscScalar       *g;
246f0b65372SJed Brown   PetscMemType       q_mem_type, g_mem_type;
247f0b65372SJed Brown   PetscFunctionBeginUser;
2482b916ea7SJeremy L Thompson   PetscCall(MatShellGetContext(J, &user));
249e2f84137SJeremy L Thompson   Vec Q_loc = user->Q_dot_loc,  // Note - Q_dot_loc has zero BCs
250e2f84137SJeremy L Thompson       G_loc;
251e2f84137SJeremy L Thompson 
252f0b65372SJed Brown   // Get local vectors
2532b916ea7SJeremy L Thompson   PetscCall(DMGetLocalVector(user->dm, &G_loc));
254f0b65372SJed Brown 
255f0b65372SJed Brown   // Global-to-local
2562b916ea7SJeremy L Thompson   PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
257f0b65372SJed Brown 
258f0b65372SJed Brown   // Place PETSc vectors in CEED vectors
2592b916ea7SJeremy L Thompson   PetscCall(VecGetArrayReadAndMemType(Q_loc, &q, &q_mem_type));
2602b916ea7SJeremy L Thompson   PetscCall(VecGetArrayAndMemType(G_loc, &g, &g_mem_type));
2612b916ea7SJeremy L Thompson   CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER, (PetscScalar *)q);
262f0b65372SJed Brown   CeedVectorSetArray(user->g_ceed, MemTypeP2C(g_mem_type), CEED_USE_POINTER, g);
263f0b65372SJed Brown 
264f0b65372SJed Brown   // Apply CEED operator
2652b916ea7SJeremy L Thompson   CeedOperatorApply(user->op_ijacobian, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE);
266f0b65372SJed Brown 
267f0b65372SJed Brown   // Restore vectors
268f0b65372SJed Brown   CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL);
269f0b65372SJed Brown   CeedVectorTakeArray(user->g_ceed, MemTypeP2C(g_mem_type), NULL);
2702b916ea7SJeremy L Thompson   PetscCall(VecRestoreArrayReadAndMemType(Q_loc, &q));
2712b916ea7SJeremy L Thompson   PetscCall(VecRestoreArrayAndMemType(G_loc, &g));
272f0b65372SJed Brown 
273f0b65372SJed Brown   // Local-to-Global
2742b916ea7SJeremy L Thompson   PetscCall(VecZeroEntries(G));
2752b916ea7SJeremy L Thompson   PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G));
276f0b65372SJed Brown 
277f0b65372SJed Brown   // Restore vectors
2782b916ea7SJeremy L Thompson   PetscCall(DMRestoreLocalVector(user->dm, &G_loc));
279f0b65372SJed Brown   PetscFunctionReturn(0);
280f0b65372SJed Brown }
281f0b65372SJed Brown 
282f0b65372SJed Brown PetscErrorCode MatGetDiagonal_NS_IJacobian(Mat A, Vec D) {
283f0b65372SJed Brown   User         user;
284f0b65372SJed Brown   Vec          D_loc;
285f0b65372SJed Brown   PetscScalar *d;
286f0b65372SJed Brown   PetscMemType mem_type;
287f0b65372SJed Brown 
288f0b65372SJed Brown   PetscFunctionBeginUser;
289f0b65372SJed Brown   MatShellGetContext(A, &user);
290f0b65372SJed Brown   PetscCall(DMGetLocalVector(user->dm, &D_loc));
291f0b65372SJed Brown   PetscCall(VecGetArrayAndMemType(D_loc, &d, &mem_type));
292f0b65372SJed Brown   CeedVectorSetArray(user->g_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, d);
2932b916ea7SJeremy L Thompson   CeedOperatorLinearAssembleDiagonal(user->op_ijacobian, user->g_ceed, CEED_REQUEST_IMMEDIATE);
294f0b65372SJed Brown   CeedVectorTakeArray(user->g_ceed, MemTypeP2C(mem_type), NULL);
295f0b65372SJed Brown   PetscCall(VecRestoreArrayAndMemType(D_loc, &d));
296f0b65372SJed Brown   PetscCall(VecZeroEntries(D));
297f0b65372SJed Brown   PetscCall(DMLocalToGlobal(user->dm, D_loc, ADD_VALUES, D));
298f0b65372SJed Brown   PetscCall(DMRestoreLocalVector(user->dm, &D_loc));
299f0b65372SJed Brown   VecViewFromOptions(D, NULL, "-diag_vec_view");
300f0b65372SJed Brown   PetscFunctionReturn(0);
301f0b65372SJed Brown }
302f0b65372SJed Brown 
3032b916ea7SJeremy L Thompson static PetscErrorCode FormPreallocation(User user, PetscBool pbdiagonal, Mat J, CeedVector *coo_values) {
304b107fddaSJed Brown   PetscCount ncoo;
305b107fddaSJed Brown   PetscInt  *rows, *cols;
306b107fddaSJed Brown 
307b107fddaSJed Brown   PetscFunctionBeginUser;
308b107fddaSJed Brown   if (pbdiagonal) {
309b107fddaSJed Brown     CeedSize l_size;
310b107fddaSJed Brown     CeedOperatorGetActiveVectorLengths(user->op_ijacobian, &l_size, NULL);
311b107fddaSJed Brown     ncoo = l_size * 5;
312b107fddaSJed Brown     rows = malloc(ncoo * sizeof(rows[0]));
313b107fddaSJed Brown     cols = malloc(ncoo * sizeof(cols[0]));
314b107fddaSJed Brown     for (PetscCount n = 0; n < l_size / 5; n++) {
315b107fddaSJed Brown       for (PetscInt i = 0; i < 5; i++) {
316b107fddaSJed Brown         for (PetscInt j = 0; j < 5; j++) {
317b107fddaSJed Brown           rows[(n * 5 + i) * 5 + j] = n * 5 + i;
318b107fddaSJed Brown           cols[(n * 5 + i) * 5 + j] = n * 5 + j;
319b107fddaSJed Brown         }
320b107fddaSJed Brown       }
321b107fddaSJed Brown     }
322b107fddaSJed Brown   } else {
3232b916ea7SJeremy L Thompson     PetscCall(CeedOperatorLinearAssembleSymbolic(user->op_ijacobian, &ncoo, &rows, &cols));
324b107fddaSJed Brown   }
325b107fddaSJed Brown   PetscCall(MatSetPreallocationCOOLocal(J, ncoo, rows, cols));
326b107fddaSJed Brown   free(rows);
327b107fddaSJed Brown   free(cols);
328b107fddaSJed Brown   CeedVectorCreate(user->ceed, ncoo, coo_values);
329b107fddaSJed Brown   PetscFunctionReturn(0);
330b107fddaSJed Brown }
331b107fddaSJed Brown 
3322b916ea7SJeremy L Thompson static PetscErrorCode FormSetValues(User user, PetscBool pbdiagonal, Mat J, CeedVector coo_values) {
333b107fddaSJed Brown   CeedMemType        mem_type = CEED_MEM_HOST;
334b107fddaSJed Brown   const PetscScalar *values;
335b107fddaSJed Brown   MatType            mat_type;
336b107fddaSJed Brown 
337b107fddaSJed Brown   PetscFunctionBeginUser;
338b107fddaSJed Brown   PetscCall(MatGetType(J, &mat_type));
3392b916ea7SJeremy L Thompson   if (strstr(mat_type, "kokkos") || strstr(mat_type, "cusparse")) mem_type = CEED_MEM_DEVICE;
340b107fddaSJed Brown   if (user->app_ctx->pmat_pbdiagonal) {
3412b916ea7SJeremy L Thompson     CeedOperatorLinearAssemblePointBlockDiagonal(user->op_ijacobian, coo_values, CEED_REQUEST_IMMEDIATE);
342b107fddaSJed Brown   } else {
343b107fddaSJed Brown     CeedOperatorLinearAssemble(user->op_ijacobian, coo_values);
344b107fddaSJed Brown   }
345b107fddaSJed Brown   CeedVectorGetArrayRead(coo_values, mem_type, &values);
346b107fddaSJed Brown   PetscCall(MatSetValuesCOO(J, values, INSERT_VALUES));
347b107fddaSJed Brown   CeedVectorRestoreArrayRead(coo_values, &values);
348b107fddaSJed Brown   PetscFunctionReturn(0);
349b107fddaSJed Brown }
350b107fddaSJed Brown 
3512b916ea7SJeremy L Thompson PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, PetscReal shift, Mat J, Mat J_pre, void *user_data) {
352f0b65372SJed Brown   User      user = *(User *)user_data;
35304855949SJed Brown   PetscBool J_is_shell, J_is_mffd, J_pre_is_shell;
354f0b65372SJed Brown   PetscFunctionBeginUser;
355a5f46a7bSJeremy L Thompson   if (user->phys->ijacobian_time_shift_label) CeedOperatorSetContextDouble(user->op_ijacobian, user->phys->ijacobian_time_shift_label, &shift);
35604855949SJed Brown   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATMFFD, &J_is_mffd));
357f0b65372SJed Brown   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATSHELL, &J_is_shell));
3582b916ea7SJeremy L Thompson   PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATSHELL, &J_pre_is_shell));
359f0b65372SJed Brown   if (!user->matrices_set_up) {
360f0b65372SJed Brown     if (J_is_shell) {
361f0b65372SJed Brown       PetscCall(MatShellSetContext(J, user));
3622b916ea7SJeremy L Thompson       PetscCall(MatShellSetOperation(J, MATOP_MULT, (void (*)(void))MatMult_NS_IJacobian));
3632b916ea7SJeremy L Thompson       PetscCall(MatShellSetOperation(J, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_NS_IJacobian));
364f0b65372SJed Brown       PetscCall(MatSetUp(J));
365f0b65372SJed Brown     }
366f0b65372SJed Brown     if (!J_pre_is_shell) {
3672b916ea7SJeremy L Thompson       PetscCall(FormPreallocation(user, user->app_ctx->pmat_pbdiagonal, J_pre, &user->coo_values_pmat));
368b107fddaSJed Brown     }
36904855949SJed Brown     if (J != J_pre && !J_is_shell && !J_is_mffd) {
370b107fddaSJed Brown       PetscCall(FormPreallocation(user, PETSC_FALSE, J, &user->coo_values_amat));
371b107fddaSJed Brown     }
372f0b65372SJed Brown     user->matrices_set_up = true;
373f0b65372SJed Brown   }
374f0b65372SJed Brown   if (!J_pre_is_shell) {
3752b916ea7SJeremy L Thompson     PetscCall(FormSetValues(user, user->app_ctx->pmat_pbdiagonal, J_pre, user->coo_values_pmat));
376f0b65372SJed Brown   }
37704855949SJed Brown   if (user->coo_values_amat) {
37804855949SJed Brown     PetscCall(FormSetValues(user, PETSC_FALSE, J, user->coo_values_amat));
37904855949SJed Brown   } else if (J_is_mffd) {
38004855949SJed Brown     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
38104855949SJed Brown     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
38204855949SJed Brown   }
383f0b65372SJed Brown   PetscFunctionReturn(0);
384f0b65372SJed Brown }
385f0b65372SJed Brown 
3862b916ea7SJeremy L Thompson PetscErrorCode WriteOutput(User user, Vec Q, PetscInt step_no, PetscScalar time) {
387a515125bSLeila Ghaffari   Vec         Q_loc;
388a515125bSLeila Ghaffari   char        file_path[PETSC_MAX_PATH_LEN];
389a515125bSLeila Ghaffari   PetscViewer viewer;
390a515125bSLeila Ghaffari   PetscFunctionBeginUser;
391a515125bSLeila Ghaffari 
392852e5969SJed Brown   if (user->app_ctx->checkpoint_vtk) {
393a515125bSLeila Ghaffari     // Set up output
3947538d537SJames Wright     PetscCall(DMGetLocalVector(user->dm, &Q_loc));
3957538d537SJames Wright     PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec"));
3967538d537SJames Wright     PetscCall(VecZeroEntries(Q_loc));
3977538d537SJames Wright     PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
398a515125bSLeila Ghaffari 
399a515125bSLeila Ghaffari     // Output
400852e5969SJed Brown     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
4017538d537SJames Wright 
4022b916ea7SJeremy L Thompson     PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path, FILE_MODE_WRITE, &viewer));
4037538d537SJames Wright     PetscCall(VecView(Q_loc, viewer));
4047538d537SJames Wright     PetscCall(PetscViewerDestroy(&viewer));
405a515125bSLeila Ghaffari     if (user->dm_viz) {
406a515125bSLeila Ghaffari       Vec         Q_refined, Q_refined_loc;
407a515125bSLeila Ghaffari       char        file_path_refined[PETSC_MAX_PATH_LEN];
408a515125bSLeila Ghaffari       PetscViewer viewer_refined;
409a515125bSLeila Ghaffari 
4107538d537SJames Wright       PetscCall(DMGetGlobalVector(user->dm_viz, &Q_refined));
4117538d537SJames Wright       PetscCall(DMGetLocalVector(user->dm_viz, &Q_refined_loc));
4127538d537SJames Wright       PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined"));
4137538d537SJames Wright 
4147538d537SJames Wright       PetscCall(MatInterpolate(user->interp_viz, Q, Q_refined));
4157538d537SJames Wright       PetscCall(VecZeroEntries(Q_refined_loc));
4162b916ea7SJeremy L Thompson       PetscCall(DMGlobalToLocal(user->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc));
4177538d537SJames Wright 
418852e5969SJed Brown       PetscCall(
419852e5969SJed Brown           PetscSNPrintf(file_path_refined, sizeof file_path_refined, "%s/nsrefined-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
4207538d537SJames Wright 
4212b916ea7SJeremy L Thompson       PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined), file_path_refined, FILE_MODE_WRITE, &viewer_refined));
4227538d537SJames Wright       PetscCall(VecView(Q_refined_loc, viewer_refined));
4237538d537SJames Wright       PetscCall(DMRestoreLocalVector(user->dm_viz, &Q_refined_loc));
4247538d537SJames Wright       PetscCall(DMRestoreGlobalVector(user->dm_viz, &Q_refined));
4257538d537SJames Wright       PetscCall(PetscViewerDestroy(&viewer_refined));
426a515125bSLeila Ghaffari     }
4277538d537SJames Wright     PetscCall(DMRestoreLocalVector(user->dm, &Q_loc));
428852e5969SJed Brown   }
429a515125bSLeila Ghaffari 
430a515125bSLeila Ghaffari   // Save data in a binary file for continuation of simulations
43191a36801SJames Wright   if (user->app_ctx->add_stepnum2bin) {
432852e5969SJed Brown     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", user->app_ctx->output_dir, step_no));
43391a36801SJames Wright   } else {
4342b916ea7SJeremy L Thompson     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", user->app_ctx->output_dir));
43591a36801SJames Wright   }
4362b916ea7SJeremy L Thompson   PetscCall(PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, &viewer));
4377538d537SJames Wright 
4389293eaa1SJed Brown   PetscInt token = FLUIDS_FILE_TOKEN;
4399293eaa1SJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &token, 1, PETSC_INT));
4409293eaa1SJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &step_no, 1, PETSC_INT));
4419293eaa1SJed Brown   time /= user->units->second;  // Dimensionalize time back
4429293eaa1SJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL));
4437538d537SJames Wright   PetscCall(VecView(Q, viewer));
4447538d537SJames Wright   PetscCall(PetscViewerDestroy(&viewer));
4457538d537SJames Wright   PetscFunctionReturn(0);
4467538d537SJames Wright }
4477538d537SJames Wright 
448c5e9980aSAdeleke O. Bankole // CSV Monitor
449c5e9980aSAdeleke O. Bankole PetscErrorCode TSMonitor_WallForce(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
450c5e9980aSAdeleke O. Bankole   User              user = ctx;
451c5e9980aSAdeleke O. Bankole   Vec               G_loc;
452c5e9980aSAdeleke O. Bankole   PetscInt          num_wall = user->app_ctx->wall_forces.num_wall, dim = 3;
453c5e9980aSAdeleke O. Bankole   const PetscInt   *walls  = user->app_ctx->wall_forces.walls;
454c5e9980aSAdeleke O. Bankole   PetscViewer       viewer = user->app_ctx->wall_forces.viewer;
455c5e9980aSAdeleke O. Bankole   PetscViewerFormat format = user->app_ctx->wall_forces.viewer_format;
456c5e9980aSAdeleke O. Bankole   PetscScalar      *reaction_force;
457c5e9980aSAdeleke O. Bankole   PetscBool         iascii;
458c5e9980aSAdeleke O. Bankole 
459c5e9980aSAdeleke O. Bankole   PetscFunctionBeginUser;
460c5e9980aSAdeleke O. Bankole   if (!viewer) PetscFunctionReturn(0);
461c5e9980aSAdeleke O. Bankole   PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
462c5e9980aSAdeleke O. Bankole   PetscCall(PetscMalloc1(num_wall * dim, &reaction_force));
463c5e9980aSAdeleke O. Bankole   PetscCall(Surface_Forces_NS(user->dm, G_loc, num_wall, walls, reaction_force));
464c5e9980aSAdeleke O. Bankole   PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
465c5e9980aSAdeleke O. Bankole 
466c5e9980aSAdeleke O. Bankole   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
467c5e9980aSAdeleke O. Bankole 
468c5e9980aSAdeleke O. Bankole   if (iascii) {
469c5e9980aSAdeleke O. Bankole     if (format == PETSC_VIEWER_ASCII_CSV && !user->app_ctx->wall_forces.header_written) {
470c5e9980aSAdeleke O. Bankole       PetscCall(PetscViewerASCIIPrintf(viewer, "Step,Time,Wall,ForceX,ForceY,ForceZ\n"));
471c5e9980aSAdeleke O. Bankole       user->app_ctx->wall_forces.header_written = PETSC_TRUE;
472c5e9980aSAdeleke O. Bankole     }
473c5e9980aSAdeleke O. Bankole     for (PetscInt w = 0; w < num_wall; w++) {
474c5e9980aSAdeleke O. Bankole       PetscInt wall = walls[w];
475c5e9980aSAdeleke O. Bankole       if (format == PETSC_VIEWER_ASCII_CSV) {
476c5e9980aSAdeleke O. Bankole         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT ",%g,%" PetscInt_FMT ",%g,%g,%g\n", step_no, time, wall,
477c5e9980aSAdeleke O. Bankole                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
478c5e9980aSAdeleke O. Bankole 
479c5e9980aSAdeleke O. Bankole       } else {
480c5e9980aSAdeleke O. Bankole         PetscCall(PetscViewerASCIIPrintf(viewer, "Wall %" PetscInt_FMT " Forces: Force_x = %12g, Force_y = %12g, Force_z = %12g\n", wall,
481c5e9980aSAdeleke O. Bankole                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
482c5e9980aSAdeleke O. Bankole       }
483c5e9980aSAdeleke O. Bankole     }
484c5e9980aSAdeleke O. Bankole   }
485c5e9980aSAdeleke O. Bankole   PetscCall(PetscFree(reaction_force));
486c5e9980aSAdeleke O. Bankole   PetscFunctionReturn(0);
487c5e9980aSAdeleke O. Bankole }
488c5e9980aSAdeleke O. Bankole 
4897538d537SJames Wright // User provided TS Monitor
4902b916ea7SJeremy L Thompson PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
4917538d537SJames Wright   User user = ctx;
4927538d537SJames Wright   PetscFunctionBeginUser;
4937538d537SJames Wright 
494852e5969SJed Brown   // Print every 'checkpoint_interval' steps
495c539088bSJames Wright   if (user->app_ctx->checkpoint_interval <= 0 || step_no % user->app_ctx->checkpoint_interval != 0 ||
496c539088bSJames Wright       (user->app_ctx->cont_steps == step_no && step_no != 0))
497c539088bSJames Wright     PetscFunctionReturn(0);
4987538d537SJames Wright 
4997538d537SJames Wright   PetscCall(WriteOutput(user, Q, step_no, time));
500a515125bSLeila Ghaffari 
501a515125bSLeila Ghaffari   PetscFunctionReturn(0);
502a515125bSLeila Ghaffari }
503a515125bSLeila Ghaffari 
504a515125bSLeila Ghaffari // TS: Create, setup, and solve
5052b916ea7SJeremy L Thompson PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, Vec *Q, PetscScalar *f_time, TS *ts) {
506a515125bSLeila Ghaffari   MPI_Comm    comm = user->comm;
507a515125bSLeila Ghaffari   TSAdapt     adapt;
508a515125bSLeila Ghaffari   PetscScalar final_time;
509a515125bSLeila Ghaffari   PetscFunctionBeginUser;
510a515125bSLeila Ghaffari 
5112b916ea7SJeremy L Thompson   PetscCall(TSCreate(comm, ts));
5122b916ea7SJeremy L Thompson   PetscCall(TSSetDM(*ts, dm));
513a515125bSLeila Ghaffari   if (phys->implicit) {
5142b916ea7SJeremy L Thompson     PetscCall(TSSetType(*ts, TSBDF));
515a515125bSLeila Ghaffari     if (user->op_ifunction) {
5162b916ea7SJeremy L Thompson       PetscCall(TSSetIFunction(*ts, NULL, IFunction_NS, &user));
517a515125bSLeila Ghaffari     } else {  // Implicit integrators can fall back to using an RHSFunction
5182b916ea7SJeremy L Thompson       PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
519a515125bSLeila Ghaffari     }
520f0b65372SJed Brown     if (user->op_ijacobian) {
5212b916ea7SJeremy L Thompson       PetscCall(DMTSSetIJacobian(dm, FormIJacobian_NS, &user));
522b107fddaSJed Brown       if (app_ctx->amat_type) {
523b107fddaSJed Brown         Mat Pmat, Amat;
5242b916ea7SJeremy L Thompson         PetscCall(DMCreateMatrix(dm, &Pmat));
5252b916ea7SJeremy L Thompson         PetscCall(DMSetMatType(dm, app_ctx->amat_type));
5262b916ea7SJeremy L Thompson         PetscCall(DMCreateMatrix(dm, &Amat));
5272b916ea7SJeremy L Thompson         PetscCall(TSSetIJacobian(*ts, Amat, Pmat, NULL, NULL));
5282b916ea7SJeremy L Thompson         PetscCall(MatDestroy(&Amat));
5292b916ea7SJeremy L Thompson         PetscCall(MatDestroy(&Pmat));
530b107fddaSJed Brown       }
531f0b65372SJed Brown     }
532a515125bSLeila Ghaffari   } else {
5332b916ea7SJeremy L Thompson     if (!user->op_rhs) SETERRQ(comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction");
5342b916ea7SJeremy L Thompson     PetscCall(TSSetType(*ts, TSRK));
5352b916ea7SJeremy L Thompson     PetscCall(TSRKSetType(*ts, TSRK5F));
5362b916ea7SJeremy L Thompson     PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
537a515125bSLeila Ghaffari   }
5382b916ea7SJeremy L Thompson   PetscCall(TSSetMaxTime(*ts, 500. * user->units->second));
5392b916ea7SJeremy L Thompson   PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER));
540f0784ed3SJed Brown   PetscCall(TSSetErrorIfStepFails(*ts, PETSC_FALSE));
5412b916ea7SJeremy L Thompson   PetscCall(TSSetTimeStep(*ts, 1.e-2 * user->units->second));
5420e1e9333SJames Wright   if (app_ctx->test_type != TESTTYPE_NONE) {
5432b916ea7SJeremy L Thompson     PetscCall(TSSetMaxSteps(*ts, 10));
5442b916ea7SJeremy L Thompson   }
5452b916ea7SJeremy L Thompson   PetscCall(TSGetAdapt(*ts, &adapt));
5462b916ea7SJeremy L Thompson   PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * user->units->second, 1.e2 * user->units->second));
5472b916ea7SJeremy L Thompson   PetscCall(TSSetFromOptions(*ts));
548*91f639d2SJames Wright   user->time_bc_set = -1.0;    // require all BCs be updated
549a515125bSLeila Ghaffari   if (!app_ctx->cont_steps) {  // print initial condition
5500e1e9333SJames Wright     if (app_ctx->test_type == TESTTYPE_NONE) {
5512b916ea7SJeremy L Thompson       PetscCall(TSMonitor_NS(*ts, 0, 0., *Q, user));
552a515125bSLeila Ghaffari     }
553a515125bSLeila Ghaffari   } else {  // continue from time of last output
554a515125bSLeila Ghaffari     PetscInt    count;
555a515125bSLeila Ghaffari     PetscViewer viewer;
5562b916ea7SJeremy L Thompson 
5579293eaa1SJed Brown     if (app_ctx->cont_time <= 0) {  // Legacy files did not include step number and time
5582b916ea7SJeremy L Thompson       PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_time_file, FILE_MODE_READ, &viewer));
5599293eaa1SJed Brown       PetscCall(PetscViewerBinaryRead(viewer, &app_ctx->cont_time, 1, &count, PETSC_REAL));
5602b916ea7SJeremy L Thompson       PetscCall(PetscViewerDestroy(&viewer));
5619293eaa1SJed Brown       PetscCheck(app_ctx->cont_steps != -1, comm, PETSC_ERR_ARG_INCOMP,
5629293eaa1SJed Brown                  "-continue step number not specified, but checkpoint file does not contain a step number (likely written by older code version)");
5639293eaa1SJed Brown     }
5649293eaa1SJed Brown     PetscCall(TSSetTime(*ts, app_ctx->cont_time * user->units->second));
56574a6f4ddSJed Brown     PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps));
566a515125bSLeila Ghaffari   }
5670e1e9333SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) {
5682b916ea7SJeremy L Thompson     PetscCall(TSMonitorSet(*ts, TSMonitor_NS, user, NULL));
5690e1e9333SJames Wright   }
570c5e9980aSAdeleke O. Bankole   if (app_ctx->wall_forces.viewer) {
571c5e9980aSAdeleke O. Bankole     PetscCall(TSMonitorSet(*ts, TSMonitor_WallForce, user, NULL));
572c5e9980aSAdeleke O. Bankole   }
573c931fa59SJames Wright   if (app_ctx->turb_spanstats_enable) {
574b0488d1fSJames Wright     PetscCall(TSMonitorSet(*ts, TSMonitor_Statistics, user, NULL));
575b8daee98SJames Wright     CeedScalar previous_time = app_ctx->cont_time * user->units->second;
576a5f46a7bSJeremy L Thompson     CeedOperatorSetContextDouble(user->spanstats.op_stats_collect, user->spanstats.previous_time_label, &previous_time);
577b0488d1fSJames Wright   }
578a515125bSLeila Ghaffari 
579a515125bSLeila Ghaffari   // Solve
58074a6f4ddSJed Brown   PetscReal start_time;
58174a6f4ddSJed Brown   PetscInt  start_step;
5822b916ea7SJeremy L Thompson   PetscCall(TSGetTime(*ts, &start_time));
58374a6f4ddSJed Brown   PetscCall(TSGetStepNumber(*ts, &start_step));
58491982731SJeremy L Thompson 
58591982731SJeremy L Thompson   PetscPreLoadBegin(PETSC_FALSE, "Fluids Solve");
58691982731SJeremy L Thompson   PetscCall(TSSetTime(*ts, start_time));
58774a6f4ddSJed Brown   PetscCall(TSSetStepNumber(*ts, start_step));
58891982731SJeremy L Thompson   if (PetscPreLoadingOn) {
58991982731SJeremy L Thompson     // LCOV_EXCL_START
59091982731SJeremy L Thompson     SNES      snes;
59191982731SJeremy L Thompson     Vec       Q_preload;
59291982731SJeremy L Thompson     PetscReal rtol;
59391982731SJeremy L Thompson     PetscCall(VecDuplicate(*Q, &Q_preload));
59491982731SJeremy L Thompson     PetscCall(VecCopy(*Q, Q_preload));
59591982731SJeremy L Thompson     PetscCall(TSGetSNES(*ts, &snes));
59691982731SJeremy L Thompson     PetscCall(SNESGetTolerances(snes, NULL, &rtol, NULL, NULL, NULL));
5972b916ea7SJeremy L Thompson     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, .99, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
59891982731SJeremy L Thompson     PetscCall(TSSetSolution(*ts, *Q));
59991982731SJeremy L Thompson     PetscCall(TSStep(*ts));
6002b916ea7SJeremy L Thompson     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, rtol, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
60191982731SJeremy L Thompson     PetscCall(VecDestroy(&Q_preload));
60291982731SJeremy L Thompson     // LCOV_EXCL_STOP
60391982731SJeremy L Thompson   } else {
6042b916ea7SJeremy L Thompson     PetscCall(PetscBarrier((PetscObject)*ts));
6052b916ea7SJeremy L Thompson     PetscCall(TSSolve(*ts, *Q));
60691982731SJeremy L Thompson   }
60791982731SJeremy L Thompson   PetscPreLoadEnd();
60891982731SJeremy L Thompson 
60991982731SJeremy L Thompson   PetscCall(TSGetSolveTime(*ts, &final_time));
610a515125bSLeila Ghaffari   *f_time = final_time;
61191982731SJeremy L Thompson 
6120e1e9333SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) {
6137538d537SJames Wright     PetscInt step_no;
6147538d537SJames Wright     PetscCall(TSGetStepNumber(*ts, &step_no));
615b0488d1fSJames Wright     if (user->app_ctx->checkpoint_interval > 0 || user->app_ctx->checkpoint_interval == -1) {
6167538d537SJames Wright       PetscCall(WriteOutput(user, *Q, step_no, final_time));
6177538d537SJames Wright     }
6187538d537SJames Wright 
61991982731SJeremy L Thompson     PetscLogEvent stage_id;
62091982731SJeremy L Thompson     PetscStageLog stage_log;
62191982731SJeremy L Thompson 
62291982731SJeremy L Thompson     PetscCall(PetscLogStageGetId("Fluids Solve", &stage_id));
62391982731SJeremy L Thompson     PetscCall(PetscLogGetStageLog(&stage_log));
6242b916ea7SJeremy L Thompson     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time taken for solution (sec): %g\n", stage_log->stageInfo[stage_id].perfInfo.time));
625a515125bSLeila Ghaffari   }
626a515125bSLeila Ghaffari   PetscFunctionReturn(0);
627a515125bSLeila Ghaffari }
628