xref: /honee/src/setupts.c (revision 2004e3acc5edde5448276a2cf66b233d507f627a)
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 
77c2cb7fc8SJames Wright PetscErrorCode UpdateContextLabel(PetscScalar *previous_value, PetscScalar update_value, CeedOperator op, CeedContextFieldLabel label) {
78c2cb7fc8SJames Wright   PetscFunctionBeginUser;
79c2cb7fc8SJames Wright 
80c2cb7fc8SJames Wright   if (*previous_value != update_value) {
81c2cb7fc8SJames Wright     if (label) {
82c2cb7fc8SJames Wright       CeedOperatorContextSetDouble(op, label, &update_value);
83c2cb7fc8SJames Wright     }
84c2cb7fc8SJames Wright     *previous_value = update_value;
85c2cb7fc8SJames Wright   }
86c2cb7fc8SJames Wright   PetscFunctionReturn(0);
87c2cb7fc8SJames Wright }
88c2cb7fc8SJames Wright 
89a515125bSLeila Ghaffari // RHS (Explicit time-stepper) function setup
90a515125bSLeila Ghaffari //   This is the RHS of the ODE, given as u_t = G(t,u)
91a515125bSLeila Ghaffari //   This function takes in a state vector Q and writes into G
92a515125bSLeila Ghaffari PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) {
93a515125bSLeila Ghaffari   User         user = *(User *)user_data;
94fd969b44SJames Wright   PetscScalar  dt;
95e2f84137SJeremy L Thompson   Vec          Q_loc = user->Q_loc, G_loc;
96a515125bSLeila Ghaffari   PetscMemType q_mem_type, g_mem_type;
97a515125bSLeila Ghaffari   PetscFunctionBeginUser;
98a515125bSLeila Ghaffari 
99e2f84137SJeremy L Thompson   // Get local vector
1002b916ea7SJeremy L Thompson   PetscCall(DMGetLocalVector(user->dm, &G_loc));
101e2f84137SJeremy L Thompson 
102e2f84137SJeremy L Thompson   // Update time dependent data
103c996854bSJames Wright   PetscCall(UpdateBoundaryValues(user, Q_loc, t));
104c2cb7fc8SJames Wright   PetscCall(UpdateContextLabel(&user->time, t, user->op_rhs, user->phys->solution_time_label));
1052b916ea7SJeremy L Thompson   PetscCall(TSGetTimeStep(ts, &dt));
106c2cb7fc8SJames Wright   PetscCall(UpdateContextLabel(&user->dt, dt, user->op_rhs, user->phys->timestep_size_label));
107a515125bSLeila Ghaffari 
108a515125bSLeila Ghaffari   // Global-to-local
1092b916ea7SJeremy L Thompson   PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
110a515125bSLeila Ghaffari 
111a515125bSLeila Ghaffari   // Place PETSc vectors in CEED vectors
112fd969b44SJames Wright   PetscCall(VecReadP2C(Q_loc, &q_mem_type, user->q_ceed));
113fd969b44SJames Wright   PetscCall(VecP2C(G_loc, &g_mem_type, user->g_ceed));
114a515125bSLeila Ghaffari 
115a515125bSLeila Ghaffari   // Apply CEED operator
1162b916ea7SJeremy L Thompson   CeedOperatorApply(user->op_rhs, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE);
117a515125bSLeila Ghaffari 
118a515125bSLeila Ghaffari   // Restore vectors
119fd969b44SJames Wright   PetscCall(VecReadC2P(user->q_ceed, q_mem_type, Q_loc));
120fd969b44SJames Wright   PetscCall(VecC2P(user->g_ceed, g_mem_type, G_loc));
121a515125bSLeila Ghaffari 
122a515125bSLeila Ghaffari   // Local-to-Global
1232b916ea7SJeremy L Thompson   PetscCall(VecZeroEntries(G));
1242b916ea7SJeremy L Thompson   PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G));
125a515125bSLeila Ghaffari 
126a515125bSLeila Ghaffari   // Inverse of the lumped mass matrix (M is Minv)
1272b916ea7SJeremy L Thompson   PetscCall(VecPointwiseMult(G, G, user->M));
128a515125bSLeila Ghaffari 
129a515125bSLeila Ghaffari   // Restore vectors
1302b916ea7SJeremy L Thompson   PetscCall(DMRestoreLocalVector(user->dm, &G_loc));
131a515125bSLeila Ghaffari 
132a515125bSLeila Ghaffari   PetscFunctionReturn(0);
133a515125bSLeila Ghaffari }
134a515125bSLeila Ghaffari 
135c5e9980aSAdeleke O. Bankole // Surface forces function setup
136c5e9980aSAdeleke O. Bankole static PetscErrorCode Surface_Forces_NS(DM dm, Vec G_loc, PetscInt num_walls, const PetscInt walls[], PetscScalar *reaction_force) {
137c5e9980aSAdeleke O. Bankole   DMLabel            face_label;
138c5e9980aSAdeleke O. Bankole   const PetscScalar *g;
139*2004e3acSAdeleke O. Bankole   PetscInt           dof, dim = 3;
140c5e9980aSAdeleke O. Bankole   MPI_Comm           comm;
141*2004e3acSAdeleke O. Bankole   PetscSection       s;
142c5e9980aSAdeleke O. Bankole 
143c5e9980aSAdeleke O. Bankole   PetscFunctionBeginUser;
144c5e9980aSAdeleke O. Bankole   PetscCall(PetscArrayzero(reaction_force, num_walls * dim));
145c5e9980aSAdeleke O. Bankole   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
146c5e9980aSAdeleke O. Bankole   PetscCall(DMGetLabel(dm, "Face Sets", &face_label));
147c5e9980aSAdeleke O. Bankole   PetscCall(VecGetArrayRead(G_loc, &g));
148c5e9980aSAdeleke O. Bankole   for (PetscInt w = 0; w < num_walls; w++) {
149c5e9980aSAdeleke O. Bankole     const PetscInt wall = walls[w];
150c5e9980aSAdeleke O. Bankole     IS             wall_is;
151*2004e3acSAdeleke O. Bankole     PetscCall(DMGetLocalSection(dm, &s));
152c5e9980aSAdeleke O. Bankole     PetscCall(DMLabelGetStratumIS(face_label, wall, &wall_is));
153c5e9980aSAdeleke O. Bankole     if (wall_is) {  // There exist such points on this process
154c5e9980aSAdeleke O. Bankole       PetscInt        num_points;
155*2004e3acSAdeleke O. Bankole       PetscInt        num_comp = 0;
156c5e9980aSAdeleke O. Bankole       const PetscInt *points;
157*2004e3acSAdeleke O. Bankole       PetscCall(PetscSectionGetFieldComponents(s, 0, &num_comp));
158c5e9980aSAdeleke O. Bankole       PetscCall(ISGetSize(wall_is, &num_points));
159c5e9980aSAdeleke O. Bankole       PetscCall(ISGetIndices(wall_is, &points));
160c5e9980aSAdeleke O. Bankole       for (PetscInt i = 0; i < num_points; i++) {
161c5e9980aSAdeleke O. Bankole         const PetscInt           p = points[i];
162c5e9980aSAdeleke O. Bankole         const StateConservative *r;
163c5e9980aSAdeleke O. Bankole         PetscCall(DMPlexPointLocalRead(dm, p, g, &r));
164*2004e3acSAdeleke O. Bankole         PetscCall(PetscSectionGetDof(s, p, &dof));
165*2004e3acSAdeleke O. Bankole         for (PetscInt node = 0; node < dof / num_comp; node++) {
166c5e9980aSAdeleke O. Bankole           for (PetscInt j = 0; j < 3; j++) {
167*2004e3acSAdeleke O. Bankole             reaction_force[w * dim + j] -= r[node].momentum[j];
168*2004e3acSAdeleke O. Bankole           }
169c5e9980aSAdeleke O. Bankole         }
170c5e9980aSAdeleke O. Bankole       }
171c5e9980aSAdeleke O. Bankole       PetscCall(ISRestoreIndices(wall_is, &points));
172c5e9980aSAdeleke O. Bankole     }
173c5e9980aSAdeleke O. Bankole     PetscCall(ISDestroy(&wall_is));
174c5e9980aSAdeleke O. Bankole   }
175c5e9980aSAdeleke O. Bankole   PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, reaction_force, dim * num_walls, MPIU_SCALAR, MPI_SUM, comm));
176c5e9980aSAdeleke O. Bankole   //  Restore Vectors
177c5e9980aSAdeleke O. Bankole   PetscCall(VecRestoreArrayRead(G_loc, &g));
178c5e9980aSAdeleke O. Bankole 
179c5e9980aSAdeleke O. Bankole   PetscFunctionReturn(0);
180c5e9980aSAdeleke O. Bankole }
181c5e9980aSAdeleke O. Bankole 
182a515125bSLeila Ghaffari // Implicit time-stepper function setup
1832b916ea7SJeremy L Thompson PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) {
184a515125bSLeila Ghaffari   User         user = *(User *)user_data;
185fd969b44SJames Wright   PetscScalar  dt;
186e2f84137SJeremy L Thompson   Vec          Q_loc = user->Q_loc, Q_dot_loc = user->Q_dot_loc, G_loc;
187a515125bSLeila Ghaffari   PetscMemType q_mem_type, q_dot_mem_type, g_mem_type;
188a515125bSLeila Ghaffari   PetscFunctionBeginUser;
189a515125bSLeila Ghaffari 
190e2f84137SJeremy L Thompson   // Get local vectors
191c5e9980aSAdeleke O. Bankole   PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
192e2f84137SJeremy L Thompson 
193e2f84137SJeremy L Thompson   // Update time dependent data
194c996854bSJames Wright   PetscCall(UpdateBoundaryValues(user, Q_loc, t));
195c2cb7fc8SJames Wright   PetscCall(UpdateContextLabel(&user->time, t, user->op_ifunction, user->phys->solution_time_label));
1962b916ea7SJeremy L Thompson   PetscCall(TSGetTimeStep(ts, &dt));
197c2cb7fc8SJames Wright   PetscCall(UpdateContextLabel(&user->dt, dt, user->op_ifunction, user->phys->timestep_size_label));
198a515125bSLeila Ghaffari 
199a515125bSLeila Ghaffari   // Global-to-local
20006108310SJames Wright   PetscCall(DMGlobalToLocalBegin(user->dm, Q, INSERT_VALUES, Q_loc));
20106108310SJames Wright   PetscCall(DMGlobalToLocalBegin(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
20206108310SJames Wright   PetscCall(DMGlobalToLocalEnd(user->dm, Q, INSERT_VALUES, Q_loc));
20306108310SJames Wright   PetscCall(DMGlobalToLocalEnd(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
204a515125bSLeila Ghaffari 
205a515125bSLeila Ghaffari   // Place PETSc vectors in CEED vectors
206fd969b44SJames Wright   PetscCall(VecReadP2C(Q_loc, &q_mem_type, user->q_ceed));
207fd969b44SJames Wright   PetscCall(VecReadP2C(Q_dot_loc, &q_dot_mem_type, user->q_dot_ceed));
208fd969b44SJames Wright   PetscCall(VecP2C(G_loc, &g_mem_type, user->g_ceed));
209a515125bSLeila Ghaffari 
210a515125bSLeila Ghaffari   // Apply CEED operator
2112b916ea7SJeremy L Thompson   CeedOperatorApply(user->op_ifunction, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE);
212a515125bSLeila Ghaffari 
213a515125bSLeila Ghaffari   // Restore vectors
214fd969b44SJames Wright   PetscCall(VecReadC2P(user->q_ceed, q_mem_type, Q_loc));
215fd969b44SJames Wright   PetscCall(VecReadC2P(user->q_dot_ceed, q_dot_mem_type, Q_dot_loc));
216fd969b44SJames Wright   PetscCall(VecC2P(user->g_ceed, g_mem_type, G_loc));
217a515125bSLeila Ghaffari 
218a515125bSLeila Ghaffari   // Local-to-Global
2192b916ea7SJeremy L Thompson   PetscCall(VecZeroEntries(G));
2202b916ea7SJeremy L Thompson   PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G));
221a515125bSLeila Ghaffari 
222a515125bSLeila Ghaffari   // Restore vectors
223c5e9980aSAdeleke O. Bankole   PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
224a515125bSLeila Ghaffari 
225a515125bSLeila Ghaffari   PetscFunctionReturn(0);
226a515125bSLeila Ghaffari }
227a515125bSLeila Ghaffari 
228f0b65372SJed Brown static PetscErrorCode MatMult_NS_IJacobian(Mat J, Vec Q, Vec G) {
229f0b65372SJed Brown   User               user;
230f0b65372SJed Brown   const PetscScalar *q;
231f0b65372SJed Brown   PetscScalar       *g;
232f0b65372SJed Brown   PetscMemType       q_mem_type, g_mem_type;
233f0b65372SJed Brown   PetscFunctionBeginUser;
2342b916ea7SJeremy L Thompson   PetscCall(MatShellGetContext(J, &user));
235e2f84137SJeremy L Thompson   Vec Q_loc = user->Q_dot_loc,  // Note - Q_dot_loc has zero BCs
236e2f84137SJeremy L Thompson       G_loc;
237e2f84137SJeremy L Thompson 
238f0b65372SJed Brown   // Get local vectors
2392b916ea7SJeremy L Thompson   PetscCall(DMGetLocalVector(user->dm, &G_loc));
240f0b65372SJed Brown 
241f0b65372SJed Brown   // Global-to-local
2422b916ea7SJeremy L Thompson   PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
243f0b65372SJed Brown 
244f0b65372SJed Brown   // Place PETSc vectors in CEED vectors
2452b916ea7SJeremy L Thompson   PetscCall(VecGetArrayReadAndMemType(Q_loc, &q, &q_mem_type));
2462b916ea7SJeremy L Thompson   PetscCall(VecGetArrayAndMemType(G_loc, &g, &g_mem_type));
2472b916ea7SJeremy L Thompson   CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER, (PetscScalar *)q);
248f0b65372SJed Brown   CeedVectorSetArray(user->g_ceed, MemTypeP2C(g_mem_type), CEED_USE_POINTER, g);
249f0b65372SJed Brown 
250f0b65372SJed Brown   // Apply CEED operator
2512b916ea7SJeremy L Thompson   CeedOperatorApply(user->op_ijacobian, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE);
252f0b65372SJed Brown 
253f0b65372SJed Brown   // Restore vectors
254f0b65372SJed Brown   CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL);
255f0b65372SJed Brown   CeedVectorTakeArray(user->g_ceed, MemTypeP2C(g_mem_type), NULL);
2562b916ea7SJeremy L Thompson   PetscCall(VecRestoreArrayReadAndMemType(Q_loc, &q));
2572b916ea7SJeremy L Thompson   PetscCall(VecRestoreArrayAndMemType(G_loc, &g));
258f0b65372SJed Brown 
259f0b65372SJed Brown   // Local-to-Global
2602b916ea7SJeremy L Thompson   PetscCall(VecZeroEntries(G));
2612b916ea7SJeremy L Thompson   PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G));
262f0b65372SJed Brown 
263f0b65372SJed Brown   // Restore vectors
2642b916ea7SJeremy L Thompson   PetscCall(DMRestoreLocalVector(user->dm, &G_loc));
265f0b65372SJed Brown   PetscFunctionReturn(0);
266f0b65372SJed Brown }
267f0b65372SJed Brown 
268f0b65372SJed Brown PetscErrorCode MatGetDiagonal_NS_IJacobian(Mat A, Vec D) {
269f0b65372SJed Brown   User         user;
270f0b65372SJed Brown   Vec          D_loc;
271f0b65372SJed Brown   PetscScalar *d;
272f0b65372SJed Brown   PetscMemType mem_type;
273f0b65372SJed Brown 
274f0b65372SJed Brown   PetscFunctionBeginUser;
275f0b65372SJed Brown   MatShellGetContext(A, &user);
276f0b65372SJed Brown   PetscCall(DMGetLocalVector(user->dm, &D_loc));
277f0b65372SJed Brown   PetscCall(VecGetArrayAndMemType(D_loc, &d, &mem_type));
278f0b65372SJed Brown   CeedVectorSetArray(user->g_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, d);
2792b916ea7SJeremy L Thompson   CeedOperatorLinearAssembleDiagonal(user->op_ijacobian, user->g_ceed, CEED_REQUEST_IMMEDIATE);
280f0b65372SJed Brown   CeedVectorTakeArray(user->g_ceed, MemTypeP2C(mem_type), NULL);
281f0b65372SJed Brown   PetscCall(VecRestoreArrayAndMemType(D_loc, &d));
282f0b65372SJed Brown   PetscCall(VecZeroEntries(D));
283f0b65372SJed Brown   PetscCall(DMLocalToGlobal(user->dm, D_loc, ADD_VALUES, D));
284f0b65372SJed Brown   PetscCall(DMRestoreLocalVector(user->dm, &D_loc));
285f0b65372SJed Brown   VecViewFromOptions(D, NULL, "-diag_vec_view");
286f0b65372SJed Brown   PetscFunctionReturn(0);
287f0b65372SJed Brown }
288f0b65372SJed Brown 
2892b916ea7SJeremy L Thompson static PetscErrorCode FormPreallocation(User user, PetscBool pbdiagonal, Mat J, CeedVector *coo_values) {
290b107fddaSJed Brown   PetscCount ncoo;
291b107fddaSJed Brown   PetscInt  *rows, *cols;
292b107fddaSJed Brown 
293b107fddaSJed Brown   PetscFunctionBeginUser;
294b107fddaSJed Brown   if (pbdiagonal) {
295b107fddaSJed Brown     CeedSize l_size;
296b107fddaSJed Brown     CeedOperatorGetActiveVectorLengths(user->op_ijacobian, &l_size, NULL);
297b107fddaSJed Brown     ncoo = l_size * 5;
298b107fddaSJed Brown     rows = malloc(ncoo * sizeof(rows[0]));
299b107fddaSJed Brown     cols = malloc(ncoo * sizeof(cols[0]));
300b107fddaSJed Brown     for (PetscCount n = 0; n < l_size / 5; n++) {
301b107fddaSJed Brown       for (PetscInt i = 0; i < 5; i++) {
302b107fddaSJed Brown         for (PetscInt j = 0; j < 5; j++) {
303b107fddaSJed Brown           rows[(n * 5 + i) * 5 + j] = n * 5 + i;
304b107fddaSJed Brown           cols[(n * 5 + i) * 5 + j] = n * 5 + j;
305b107fddaSJed Brown         }
306b107fddaSJed Brown       }
307b107fddaSJed Brown     }
308b107fddaSJed Brown   } else {
3092b916ea7SJeremy L Thompson     PetscCall(CeedOperatorLinearAssembleSymbolic(user->op_ijacobian, &ncoo, &rows, &cols));
310b107fddaSJed Brown   }
311b107fddaSJed Brown   PetscCall(MatSetPreallocationCOOLocal(J, ncoo, rows, cols));
312b107fddaSJed Brown   free(rows);
313b107fddaSJed Brown   free(cols);
314b107fddaSJed Brown   CeedVectorCreate(user->ceed, ncoo, coo_values);
315b107fddaSJed Brown   PetscFunctionReturn(0);
316b107fddaSJed Brown }
317b107fddaSJed Brown 
3182b916ea7SJeremy L Thompson static PetscErrorCode FormSetValues(User user, PetscBool pbdiagonal, Mat J, CeedVector coo_values) {
319b107fddaSJed Brown   CeedMemType        mem_type = CEED_MEM_HOST;
320b107fddaSJed Brown   const PetscScalar *values;
321b107fddaSJed Brown   MatType            mat_type;
322b107fddaSJed Brown 
323b107fddaSJed Brown   PetscFunctionBeginUser;
324b107fddaSJed Brown   PetscCall(MatGetType(J, &mat_type));
3252b916ea7SJeremy L Thompson   if (strstr(mat_type, "kokkos") || strstr(mat_type, "cusparse")) mem_type = CEED_MEM_DEVICE;
326b107fddaSJed Brown   if (user->app_ctx->pmat_pbdiagonal) {
3272b916ea7SJeremy L Thompson     CeedOperatorLinearAssemblePointBlockDiagonal(user->op_ijacobian, coo_values, CEED_REQUEST_IMMEDIATE);
328b107fddaSJed Brown   } else {
329b107fddaSJed Brown     CeedOperatorLinearAssemble(user->op_ijacobian, coo_values);
330b107fddaSJed Brown   }
331b107fddaSJed Brown   CeedVectorGetArrayRead(coo_values, mem_type, &values);
332b107fddaSJed Brown   PetscCall(MatSetValuesCOO(J, values, INSERT_VALUES));
333b107fddaSJed Brown   CeedVectorRestoreArrayRead(coo_values, &values);
334b107fddaSJed Brown   PetscFunctionReturn(0);
335b107fddaSJed Brown }
336b107fddaSJed Brown 
3372b916ea7SJeremy L Thompson PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, PetscReal shift, Mat J, Mat J_pre, void *user_data) {
338f0b65372SJed Brown   User      user = *(User *)user_data;
33904855949SJed Brown   PetscBool J_is_shell, J_is_mffd, J_pre_is_shell;
340f0b65372SJed Brown   PetscFunctionBeginUser;
3412b916ea7SJeremy L Thompson   if (user->phys->ijacobian_time_shift_label) CeedOperatorContextSetDouble(user->op_ijacobian, user->phys->ijacobian_time_shift_label, &shift);
34204855949SJed Brown   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATMFFD, &J_is_mffd));
343f0b65372SJed Brown   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATSHELL, &J_is_shell));
3442b916ea7SJeremy L Thompson   PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATSHELL, &J_pre_is_shell));
345f0b65372SJed Brown   if (!user->matrices_set_up) {
346f0b65372SJed Brown     if (J_is_shell) {
347f0b65372SJed Brown       PetscCall(MatShellSetContext(J, user));
3482b916ea7SJeremy L Thompson       PetscCall(MatShellSetOperation(J, MATOP_MULT, (void (*)(void))MatMult_NS_IJacobian));
3492b916ea7SJeremy L Thompson       PetscCall(MatShellSetOperation(J, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_NS_IJacobian));
350f0b65372SJed Brown       PetscCall(MatSetUp(J));
351f0b65372SJed Brown     }
352f0b65372SJed Brown     if (!J_pre_is_shell) {
3532b916ea7SJeremy L Thompson       PetscCall(FormPreallocation(user, user->app_ctx->pmat_pbdiagonal, J_pre, &user->coo_values_pmat));
354b107fddaSJed Brown     }
35504855949SJed Brown     if (J != J_pre && !J_is_shell && !J_is_mffd) {
356b107fddaSJed Brown       PetscCall(FormPreallocation(user, PETSC_FALSE, J, &user->coo_values_amat));
357b107fddaSJed Brown     }
358f0b65372SJed Brown     user->matrices_set_up = true;
359f0b65372SJed Brown   }
360f0b65372SJed Brown   if (!J_pre_is_shell) {
3612b916ea7SJeremy L Thompson     PetscCall(FormSetValues(user, user->app_ctx->pmat_pbdiagonal, J_pre, user->coo_values_pmat));
362f0b65372SJed Brown   }
36304855949SJed Brown   if (user->coo_values_amat) {
36404855949SJed Brown     PetscCall(FormSetValues(user, PETSC_FALSE, J, user->coo_values_amat));
36504855949SJed Brown   } else if (J_is_mffd) {
36604855949SJed Brown     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
36704855949SJed Brown     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
36804855949SJed Brown   }
369f0b65372SJed Brown   PetscFunctionReturn(0);
370f0b65372SJed Brown }
371f0b65372SJed Brown 
3722b916ea7SJeremy L Thompson PetscErrorCode WriteOutput(User user, Vec Q, PetscInt step_no, PetscScalar time) {
373a515125bSLeila Ghaffari   Vec         Q_loc;
374a515125bSLeila Ghaffari   char        file_path[PETSC_MAX_PATH_LEN];
375a515125bSLeila Ghaffari   PetscViewer viewer;
376a515125bSLeila Ghaffari   PetscFunctionBeginUser;
377a515125bSLeila Ghaffari 
378852e5969SJed Brown   if (user->app_ctx->checkpoint_vtk) {
379a515125bSLeila Ghaffari     // Set up output
3807538d537SJames Wright     PetscCall(DMGetLocalVector(user->dm, &Q_loc));
3817538d537SJames Wright     PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec"));
3827538d537SJames Wright     PetscCall(VecZeroEntries(Q_loc));
3837538d537SJames Wright     PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
384a515125bSLeila Ghaffari 
385a515125bSLeila Ghaffari     // Output
386852e5969SJed Brown     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
3877538d537SJames Wright 
3882b916ea7SJeremy L Thompson     PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path, FILE_MODE_WRITE, &viewer));
3897538d537SJames Wright     PetscCall(VecView(Q_loc, viewer));
3907538d537SJames Wright     PetscCall(PetscViewerDestroy(&viewer));
391a515125bSLeila Ghaffari     if (user->dm_viz) {
392a515125bSLeila Ghaffari       Vec         Q_refined, Q_refined_loc;
393a515125bSLeila Ghaffari       char        file_path_refined[PETSC_MAX_PATH_LEN];
394a515125bSLeila Ghaffari       PetscViewer viewer_refined;
395a515125bSLeila Ghaffari 
3967538d537SJames Wright       PetscCall(DMGetGlobalVector(user->dm_viz, &Q_refined));
3977538d537SJames Wright       PetscCall(DMGetLocalVector(user->dm_viz, &Q_refined_loc));
3987538d537SJames Wright       PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined"));
3997538d537SJames Wright 
4007538d537SJames Wright       PetscCall(MatInterpolate(user->interp_viz, Q, Q_refined));
4017538d537SJames Wright       PetscCall(VecZeroEntries(Q_refined_loc));
4022b916ea7SJeremy L Thompson       PetscCall(DMGlobalToLocal(user->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc));
4037538d537SJames Wright 
404852e5969SJed Brown       PetscCall(
405852e5969SJed Brown           PetscSNPrintf(file_path_refined, sizeof file_path_refined, "%s/nsrefined-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
4067538d537SJames Wright 
4072b916ea7SJeremy L Thompson       PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined), file_path_refined, FILE_MODE_WRITE, &viewer_refined));
4087538d537SJames Wright       PetscCall(VecView(Q_refined_loc, viewer_refined));
4097538d537SJames Wright       PetscCall(DMRestoreLocalVector(user->dm_viz, &Q_refined_loc));
4107538d537SJames Wright       PetscCall(DMRestoreGlobalVector(user->dm_viz, &Q_refined));
4117538d537SJames Wright       PetscCall(PetscViewerDestroy(&viewer_refined));
412a515125bSLeila Ghaffari     }
4137538d537SJames Wright     PetscCall(DMRestoreLocalVector(user->dm, &Q_loc));
414852e5969SJed Brown   }
415a515125bSLeila Ghaffari 
416a515125bSLeila Ghaffari   // Save data in a binary file for continuation of simulations
41791a36801SJames Wright   if (user->app_ctx->add_stepnum2bin) {
418852e5969SJed Brown     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", user->app_ctx->output_dir, step_no));
41991a36801SJames Wright   } else {
4202b916ea7SJeremy L Thompson     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", user->app_ctx->output_dir));
42191a36801SJames Wright   }
4222b916ea7SJeremy L Thompson   PetscCall(PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, &viewer));
4237538d537SJames Wright 
4249293eaa1SJed Brown   PetscInt token = FLUIDS_FILE_TOKEN;
4259293eaa1SJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &token, 1, PETSC_INT));
4269293eaa1SJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &step_no, 1, PETSC_INT));
4279293eaa1SJed Brown   time /= user->units->second;  // Dimensionalize time back
4289293eaa1SJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL));
4297538d537SJames Wright   PetscCall(VecView(Q, viewer));
4307538d537SJames Wright   PetscCall(PetscViewerDestroy(&viewer));
4317538d537SJames Wright   PetscFunctionReturn(0);
4327538d537SJames Wright }
4337538d537SJames Wright 
434c5e9980aSAdeleke O. Bankole // CSV Monitor
435c5e9980aSAdeleke O. Bankole PetscErrorCode TSMonitor_WallForce(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
436c5e9980aSAdeleke O. Bankole   User              user = ctx;
437c5e9980aSAdeleke O. Bankole   Vec               G_loc;
438c5e9980aSAdeleke O. Bankole   PetscInt          num_wall = user->app_ctx->wall_forces.num_wall, dim = 3;
439c5e9980aSAdeleke O. Bankole   const PetscInt   *walls  = user->app_ctx->wall_forces.walls;
440c5e9980aSAdeleke O. Bankole   PetscViewer       viewer = user->app_ctx->wall_forces.viewer;
441c5e9980aSAdeleke O. Bankole   PetscViewerFormat format = user->app_ctx->wall_forces.viewer_format;
442c5e9980aSAdeleke O. Bankole   PetscScalar      *reaction_force;
443c5e9980aSAdeleke O. Bankole   PetscBool         iascii;
444c5e9980aSAdeleke O. Bankole 
445c5e9980aSAdeleke O. Bankole   PetscFunctionBeginUser;
446c5e9980aSAdeleke O. Bankole   if (!viewer) PetscFunctionReturn(0);
447c5e9980aSAdeleke O. Bankole   PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
448c5e9980aSAdeleke O. Bankole   PetscCall(PetscMalloc1(num_wall * dim, &reaction_force));
449c5e9980aSAdeleke O. Bankole   PetscCall(Surface_Forces_NS(user->dm, G_loc, num_wall, walls, reaction_force));
450c5e9980aSAdeleke O. Bankole   PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
451c5e9980aSAdeleke O. Bankole 
452c5e9980aSAdeleke O. Bankole   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
453c5e9980aSAdeleke O. Bankole 
454c5e9980aSAdeleke O. Bankole   if (iascii) {
455c5e9980aSAdeleke O. Bankole     if (format == PETSC_VIEWER_ASCII_CSV && !user->app_ctx->wall_forces.header_written) {
456c5e9980aSAdeleke O. Bankole       PetscCall(PetscViewerASCIIPrintf(viewer, "Step,Time,Wall,ForceX,ForceY,ForceZ\n"));
457c5e9980aSAdeleke O. Bankole       user->app_ctx->wall_forces.header_written = PETSC_TRUE;
458c5e9980aSAdeleke O. Bankole     }
459c5e9980aSAdeleke O. Bankole     for (PetscInt w = 0; w < num_wall; w++) {
460c5e9980aSAdeleke O. Bankole       PetscInt wall = walls[w];
461c5e9980aSAdeleke O. Bankole       if (format == PETSC_VIEWER_ASCII_CSV) {
462c5e9980aSAdeleke O. Bankole         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT ",%g,%" PetscInt_FMT ",%g,%g,%g\n", step_no, time, wall,
463c5e9980aSAdeleke O. Bankole                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
464c5e9980aSAdeleke O. Bankole 
465c5e9980aSAdeleke O. Bankole       } else {
466c5e9980aSAdeleke O. Bankole         PetscCall(PetscViewerASCIIPrintf(viewer, "Wall %" PetscInt_FMT " Forces: Force_x = %12g, Force_y = %12g, Force_z = %12g\n", wall,
467c5e9980aSAdeleke O. Bankole                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
468c5e9980aSAdeleke O. Bankole       }
469c5e9980aSAdeleke O. Bankole     }
470c5e9980aSAdeleke O. Bankole   }
471c5e9980aSAdeleke O. Bankole   PetscCall(PetscFree(reaction_force));
472c5e9980aSAdeleke O. Bankole   PetscFunctionReturn(0);
473c5e9980aSAdeleke O. Bankole }
474c5e9980aSAdeleke O. Bankole 
4757538d537SJames Wright // User provided TS Monitor
4762b916ea7SJeremy L Thompson PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
4777538d537SJames Wright   User user = ctx;
4787538d537SJames Wright   PetscFunctionBeginUser;
4797538d537SJames Wright 
480852e5969SJed Brown   // Print every 'checkpoint_interval' steps
481c539088bSJames Wright   if (user->app_ctx->checkpoint_interval <= 0 || step_no % user->app_ctx->checkpoint_interval != 0 ||
482c539088bSJames Wright       (user->app_ctx->cont_steps == step_no && step_no != 0))
483c539088bSJames Wright     PetscFunctionReturn(0);
4847538d537SJames Wright 
4857538d537SJames Wright   PetscCall(WriteOutput(user, Q, step_no, time));
486a515125bSLeila Ghaffari 
487a515125bSLeila Ghaffari   PetscFunctionReturn(0);
488a515125bSLeila Ghaffari }
489a515125bSLeila Ghaffari 
490a515125bSLeila Ghaffari // TS: Create, setup, and solve
4912b916ea7SJeremy L Thompson PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, Vec *Q, PetscScalar *f_time, TS *ts) {
492a515125bSLeila Ghaffari   MPI_Comm    comm = user->comm;
493a515125bSLeila Ghaffari   TSAdapt     adapt;
494a515125bSLeila Ghaffari   PetscScalar final_time;
495a515125bSLeila Ghaffari   PetscFunctionBeginUser;
496a515125bSLeila Ghaffari 
4972b916ea7SJeremy L Thompson   PetscCall(TSCreate(comm, ts));
4982b916ea7SJeremy L Thompson   PetscCall(TSSetDM(*ts, dm));
499a515125bSLeila Ghaffari   if (phys->implicit) {
5002b916ea7SJeremy L Thompson     PetscCall(TSSetType(*ts, TSBDF));
501a515125bSLeila Ghaffari     if (user->op_ifunction) {
5022b916ea7SJeremy L Thompson       PetscCall(TSSetIFunction(*ts, NULL, IFunction_NS, &user));
503a515125bSLeila Ghaffari     } else {  // Implicit integrators can fall back to using an RHSFunction
5042b916ea7SJeremy L Thompson       PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
505a515125bSLeila Ghaffari     }
506f0b65372SJed Brown     if (user->op_ijacobian) {
5072b916ea7SJeremy L Thompson       PetscCall(DMTSSetIJacobian(dm, FormIJacobian_NS, &user));
508b107fddaSJed Brown       if (app_ctx->amat_type) {
509b107fddaSJed Brown         Mat Pmat, Amat;
5102b916ea7SJeremy L Thompson         PetscCall(DMCreateMatrix(dm, &Pmat));
5112b916ea7SJeremy L Thompson         PetscCall(DMSetMatType(dm, app_ctx->amat_type));
5122b916ea7SJeremy L Thompson         PetscCall(DMCreateMatrix(dm, &Amat));
5132b916ea7SJeremy L Thompson         PetscCall(TSSetIJacobian(*ts, Amat, Pmat, NULL, NULL));
5142b916ea7SJeremy L Thompson         PetscCall(MatDestroy(&Amat));
5152b916ea7SJeremy L Thompson         PetscCall(MatDestroy(&Pmat));
516b107fddaSJed Brown       }
517f0b65372SJed Brown     }
518a515125bSLeila Ghaffari   } else {
5192b916ea7SJeremy L Thompson     if (!user->op_rhs) SETERRQ(comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction");
5202b916ea7SJeremy L Thompson     PetscCall(TSSetType(*ts, TSRK));
5212b916ea7SJeremy L Thompson     PetscCall(TSRKSetType(*ts, TSRK5F));
5222b916ea7SJeremy L Thompson     PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
523a515125bSLeila Ghaffari   }
5242b916ea7SJeremy L Thompson   PetscCall(TSSetMaxTime(*ts, 500. * user->units->second));
5252b916ea7SJeremy L Thompson   PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER));
526f0784ed3SJed Brown   PetscCall(TSSetErrorIfStepFails(*ts, PETSC_FALSE));
5272b916ea7SJeremy L Thompson   PetscCall(TSSetTimeStep(*ts, 1.e-2 * user->units->second));
5280e1e9333SJames Wright   if (app_ctx->test_type != TESTTYPE_NONE) {
5292b916ea7SJeremy L Thompson     PetscCall(TSSetMaxSteps(*ts, 10));
5302b916ea7SJeremy L Thompson   }
5312b916ea7SJeremy L Thompson   PetscCall(TSGetAdapt(*ts, &adapt));
5322b916ea7SJeremy L Thompson   PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * user->units->second, 1.e2 * user->units->second));
5332b916ea7SJeremy L Thompson   PetscCall(TSSetFromOptions(*ts));
534c996854bSJames Wright   user->time = user->time_bc_set = -1.0;  // require all BCs and ctx to be updated
535e2f84137SJeremy L Thompson   user->dt                       = -1.0;
536a515125bSLeila Ghaffari   if (!app_ctx->cont_steps) {  // print initial condition
5370e1e9333SJames Wright     if (app_ctx->test_type == TESTTYPE_NONE) {
5382b916ea7SJeremy L Thompson       PetscCall(TSMonitor_NS(*ts, 0, 0., *Q, user));
539a515125bSLeila Ghaffari     }
540a515125bSLeila Ghaffari   } else {  // continue from time of last output
541a515125bSLeila Ghaffari     PetscInt    count;
542a515125bSLeila Ghaffari     PetscViewer viewer;
5432b916ea7SJeremy L Thompson 
5449293eaa1SJed Brown     if (app_ctx->cont_time <= 0) {  // Legacy files did not include step number and time
5452b916ea7SJeremy L Thompson       PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_time_file, FILE_MODE_READ, &viewer));
5469293eaa1SJed Brown       PetscCall(PetscViewerBinaryRead(viewer, &app_ctx->cont_time, 1, &count, PETSC_REAL));
5472b916ea7SJeremy L Thompson       PetscCall(PetscViewerDestroy(&viewer));
5489293eaa1SJed Brown       PetscCheck(app_ctx->cont_steps != -1, comm, PETSC_ERR_ARG_INCOMP,
5499293eaa1SJed Brown                  "-continue step number not specified, but checkpoint file does not contain a step number (likely written by older code version)");
5509293eaa1SJed Brown     }
5519293eaa1SJed Brown     PetscCall(TSSetTime(*ts, app_ctx->cont_time * user->units->second));
55274a6f4ddSJed Brown     PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps));
553a515125bSLeila Ghaffari   }
5540e1e9333SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) {
5552b916ea7SJeremy L Thompson     PetscCall(TSMonitorSet(*ts, TSMonitor_NS, user, NULL));
5560e1e9333SJames Wright   }
557c5e9980aSAdeleke O. Bankole   if (app_ctx->wall_forces.viewer) {
558c5e9980aSAdeleke O. Bankole     PetscCall(TSMonitorSet(*ts, TSMonitor_WallForce, user, NULL));
559c5e9980aSAdeleke O. Bankole   }
560c931fa59SJames Wright   if (app_ctx->turb_spanstats_enable) {
561b0488d1fSJames Wright     PetscCall(TSMonitorSet(*ts, TSMonitor_Statistics, user, NULL));
562b8daee98SJames Wright     CeedScalar previous_time = app_ctx->cont_time * user->units->second;
563b8daee98SJames Wright     CeedOperatorContextSetDouble(user->spanstats.op_stats_collect, user->spanstats.previous_time_label, &previous_time);
564b0488d1fSJames Wright   }
565a515125bSLeila Ghaffari 
566a515125bSLeila Ghaffari   // Solve
56774a6f4ddSJed Brown   PetscReal start_time;
56874a6f4ddSJed Brown   PetscInt  start_step;
5692b916ea7SJeremy L Thompson   PetscCall(TSGetTime(*ts, &start_time));
57074a6f4ddSJed Brown   PetscCall(TSGetStepNumber(*ts, &start_step));
57191982731SJeremy L Thompson 
57291982731SJeremy L Thompson   PetscPreLoadBegin(PETSC_FALSE, "Fluids Solve");
57391982731SJeremy L Thompson   PetscCall(TSSetTime(*ts, start_time));
57474a6f4ddSJed Brown   PetscCall(TSSetStepNumber(*ts, start_step));
57591982731SJeremy L Thompson   if (PetscPreLoadingOn) {
57691982731SJeremy L Thompson     // LCOV_EXCL_START
57791982731SJeremy L Thompson     SNES      snes;
57891982731SJeremy L Thompson     Vec       Q_preload;
57991982731SJeremy L Thompson     PetscReal rtol;
58091982731SJeremy L Thompson     PetscCall(VecDuplicate(*Q, &Q_preload));
58191982731SJeremy L Thompson     PetscCall(VecCopy(*Q, Q_preload));
58291982731SJeremy L Thompson     PetscCall(TSGetSNES(*ts, &snes));
58391982731SJeremy L Thompson     PetscCall(SNESGetTolerances(snes, NULL, &rtol, NULL, NULL, NULL));
5842b916ea7SJeremy L Thompson     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, .99, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
58591982731SJeremy L Thompson     PetscCall(TSSetSolution(*ts, *Q));
58691982731SJeremy L Thompson     PetscCall(TSStep(*ts));
5872b916ea7SJeremy L Thompson     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, rtol, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
58891982731SJeremy L Thompson     PetscCall(VecDestroy(&Q_preload));
58991982731SJeremy L Thompson     // LCOV_EXCL_STOP
59091982731SJeremy L Thompson   } else {
5912b916ea7SJeremy L Thompson     PetscCall(PetscBarrier((PetscObject)*ts));
5922b916ea7SJeremy L Thompson     PetscCall(TSSolve(*ts, *Q));
59391982731SJeremy L Thompson   }
59491982731SJeremy L Thompson   PetscPreLoadEnd();
59591982731SJeremy L Thompson 
59691982731SJeremy L Thompson   PetscCall(TSGetSolveTime(*ts, &final_time));
597a515125bSLeila Ghaffari   *f_time = final_time;
59891982731SJeremy L Thompson 
5990e1e9333SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) {
6007538d537SJames Wright     PetscInt step_no;
6017538d537SJames Wright     PetscCall(TSGetStepNumber(*ts, &step_no));
602b0488d1fSJames Wright     if (user->app_ctx->checkpoint_interval > 0 || user->app_ctx->checkpoint_interval == -1) {
6037538d537SJames Wright       PetscCall(WriteOutput(user, *Q, step_no, final_time));
6047538d537SJames Wright     }
6057538d537SJames Wright 
60691982731SJeremy L Thompson     PetscLogEvent stage_id;
60791982731SJeremy L Thompson     PetscStageLog stage_log;
60891982731SJeremy L Thompson 
60991982731SJeremy L Thompson     PetscCall(PetscLogStageGetId("Fluids Solve", &stage_id));
61091982731SJeremy L Thompson     PetscCall(PetscLogGetStageLog(&stage_log));
6112b916ea7SJeremy L Thompson     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time taken for solution (sec): %g\n", stage_log->stageInfo[stage_id].perfInfo.time));
612a515125bSLeila Ghaffari   }
613a515125bSLeila Ghaffari   PetscFunctionReturn(0);
614a515125bSLeila Ghaffari }
615