xref: /honee/src/setupts.c (revision c5e9980ad559552e5e02639665d11e7187691949)
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"
12*c5e9980aSAdeleke 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   CeedScalar  *m;
41a515125bSLeila Ghaffari   PetscMemType m_mem_type;
422b916ea7SJeremy L Thompson   PetscCall(DMGetLocalVector(dm, &M_loc));
432b916ea7SJeremy L Thompson   PetscCall(VecGetArrayAndMemType(M_loc, (PetscScalar **)&m, &m_mem_type));
44a515125bSLeila Ghaffari   CeedVectorSetArray(m_ceed, MemTypeP2C(m_mem_type), CEED_USE_POINTER, m);
45a515125bSLeila Ghaffari 
46a515125bSLeila Ghaffari   // Apply CEED Operator
47a515125bSLeila Ghaffari   CeedOperatorApply(op_mass, ones_vec, m_ceed, CEED_REQUEST_IMMEDIATE);
48a515125bSLeila Ghaffari 
49a515125bSLeila Ghaffari   // Restore vectors
50a515125bSLeila Ghaffari   CeedVectorTakeArray(m_ceed, MemTypeP2C(m_mem_type), NULL);
512b916ea7SJeremy L Thompson   PetscCall(VecRestoreArrayReadAndMemType(M_loc, (const PetscScalar **)&m));
52a515125bSLeila Ghaffari 
53a515125bSLeila Ghaffari   // Local-to-Global
542b916ea7SJeremy L Thompson   PetscCall(VecZeroEntries(M));
552b916ea7SJeremy L Thompson   PetscCall(DMLocalToGlobal(dm, M_loc, ADD_VALUES, M));
562b916ea7SJeremy L Thompson   PetscCall(DMRestoreLocalVector(dm, &M_loc));
57a515125bSLeila Ghaffari 
58a515125bSLeila Ghaffari   // Invert diagonally lumped mass vector for RHS function
592b916ea7SJeremy L Thompson   PetscCall(VecReciprocal(M));
60a515125bSLeila Ghaffari 
61a515125bSLeila Ghaffari   // Cleanup
62a515125bSLeila Ghaffari   CeedVectorDestroy(&ones_vec);
63a515125bSLeila Ghaffari   CeedVectorDestroy(&m_ceed);
64a515125bSLeila Ghaffari   CeedQFunctionDestroy(&qf_mass);
65a515125bSLeila Ghaffari   CeedOperatorDestroy(&op_mass);
66a515125bSLeila Ghaffari 
67a515125bSLeila Ghaffari   PetscFunctionReturn(0);
68a515125bSLeila Ghaffari }
69a515125bSLeila Ghaffari 
70c996854bSJames Wright // Insert Boundary values if it's a new time
71c996854bSJames Wright PetscErrorCode UpdateBoundaryValues(User user, Vec Q_loc, PetscReal t) {
72c996854bSJames Wright   PetscFunctionBeginUser;
73c996854bSJames Wright   if (user->time_bc_set != t) {
74c996854bSJames Wright     PetscCall(DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL));
75c996854bSJames Wright     user->time_bc_set = t;
76c996854bSJames Wright   }
77c996854bSJames Wright   PetscFunctionReturn(0);
78c996854bSJames Wright }
79c996854bSJames Wright 
80c2cb7fc8SJames Wright PetscErrorCode UpdateContextLabel(PetscScalar *previous_value, PetscScalar update_value, CeedOperator op, CeedContextFieldLabel label) {
81c2cb7fc8SJames Wright   PetscFunctionBeginUser;
82c2cb7fc8SJames Wright 
83c2cb7fc8SJames Wright   if (*previous_value != update_value) {
84c2cb7fc8SJames Wright     if (label) {
85c2cb7fc8SJames Wright       CeedOperatorContextSetDouble(op, label, &update_value);
86c2cb7fc8SJames Wright     }
87c2cb7fc8SJames Wright     *previous_value = update_value;
88c2cb7fc8SJames Wright   }
89c2cb7fc8SJames Wright   PetscFunctionReturn(0);
90c2cb7fc8SJames Wright }
91c2cb7fc8SJames Wright 
92a515125bSLeila Ghaffari // RHS (Explicit time-stepper) function setup
93a515125bSLeila Ghaffari //   This is the RHS of the ODE, given as u_t = G(t,u)
94a515125bSLeila Ghaffari //   This function takes in a state vector Q and writes into G
95a515125bSLeila Ghaffari PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) {
96a515125bSLeila Ghaffari   User         user = *(User *)user_data;
97c2cb7fc8SJames Wright   PetscScalar *q, *g, dt;
98e2f84137SJeremy L Thompson   Vec          Q_loc = user->Q_loc, G_loc;
99a515125bSLeila Ghaffari   PetscMemType q_mem_type, g_mem_type;
100a515125bSLeila Ghaffari   PetscFunctionBeginUser;
101a515125bSLeila Ghaffari 
102e2f84137SJeremy L Thompson   // Get local vector
1032b916ea7SJeremy L Thompson   PetscCall(DMGetLocalVector(user->dm, &G_loc));
104e2f84137SJeremy L Thompson 
105e2f84137SJeremy L Thompson   // Update time dependent data
106c996854bSJames Wright   PetscCall(UpdateBoundaryValues(user, Q_loc, t));
107c2cb7fc8SJames Wright   PetscCall(UpdateContextLabel(&user->time, t, user->op_rhs, user->phys->solution_time_label));
1082b916ea7SJeremy L Thompson   PetscCall(TSGetTimeStep(ts, &dt));
109c2cb7fc8SJames Wright   PetscCall(UpdateContextLabel(&user->dt, dt, user->op_rhs, user->phys->timestep_size_label));
110a515125bSLeila Ghaffari 
111a515125bSLeila Ghaffari   // Global-to-local
1122b916ea7SJeremy L Thompson   PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
113a515125bSLeila Ghaffari 
114a515125bSLeila Ghaffari   // Place PETSc vectors in CEED vectors
1152b916ea7SJeremy L Thompson   PetscCall(VecGetArrayReadAndMemType(Q_loc, (const PetscScalar **)&q, &q_mem_type));
1162b916ea7SJeremy L Thompson   PetscCall(VecGetArrayAndMemType(G_loc, &g, &g_mem_type));
117a515125bSLeila Ghaffari   CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER, q);
118a515125bSLeila Ghaffari   CeedVectorSetArray(user->g_ceed, MemTypeP2C(g_mem_type), CEED_USE_POINTER, g);
119a515125bSLeila Ghaffari 
120a515125bSLeila Ghaffari   // Apply CEED operator
1212b916ea7SJeremy L Thompson   CeedOperatorApply(user->op_rhs, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE);
122a515125bSLeila Ghaffari 
123a515125bSLeila Ghaffari   // Restore vectors
124a515125bSLeila Ghaffari   CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL);
125a515125bSLeila Ghaffari   CeedVectorTakeArray(user->g_ceed, MemTypeP2C(g_mem_type), NULL);
1262b916ea7SJeremy L Thompson   PetscCall(VecRestoreArrayReadAndMemType(Q_loc, (const PetscScalar **)&q));
1272b916ea7SJeremy L Thompson   PetscCall(VecRestoreArrayAndMemType(G_loc, &g));
128a515125bSLeila Ghaffari 
129a515125bSLeila Ghaffari   // Local-to-Global
1302b916ea7SJeremy L Thompson   PetscCall(VecZeroEntries(G));
1312b916ea7SJeremy L Thompson   PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G));
132a515125bSLeila Ghaffari 
133a515125bSLeila Ghaffari   // Inverse of the lumped mass matrix (M is Minv)
1342b916ea7SJeremy L Thompson   PetscCall(VecPointwiseMult(G, G, user->M));
135a515125bSLeila Ghaffari 
136a515125bSLeila Ghaffari   // Restore vectors
1372b916ea7SJeremy L Thompson   PetscCall(DMRestoreLocalVector(user->dm, &G_loc));
138a515125bSLeila Ghaffari 
139a515125bSLeila Ghaffari   PetscFunctionReturn(0);
140a515125bSLeila Ghaffari }
141a515125bSLeila Ghaffari 
142*c5e9980aSAdeleke O. Bankole // Surface forces function setup
143*c5e9980aSAdeleke O. Bankole static PetscErrorCode Surface_Forces_NS(DM dm, Vec G_loc, PetscInt num_walls, const PetscInt walls[], PetscScalar *reaction_force) {
144*c5e9980aSAdeleke O. Bankole   DMLabel            face_label;
145*c5e9980aSAdeleke O. Bankole   const PetscScalar *g;
146*c5e9980aSAdeleke O. Bankole   PetscInt           dim = 3;
147*c5e9980aSAdeleke O. Bankole   MPI_Comm           comm;
148*c5e9980aSAdeleke O. Bankole 
149*c5e9980aSAdeleke O. Bankole   PetscFunctionBeginUser;
150*c5e9980aSAdeleke O. Bankole   PetscCall(PetscArrayzero(reaction_force, num_walls * dim));
151*c5e9980aSAdeleke O. Bankole   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
152*c5e9980aSAdeleke O. Bankole   PetscCall(DMGetLabel(dm, "Face Sets", &face_label));
153*c5e9980aSAdeleke O. Bankole   PetscCall(VecGetArrayRead(G_loc, &g));
154*c5e9980aSAdeleke O. Bankole   for (PetscInt w = 0; w < num_walls; w++) {
155*c5e9980aSAdeleke O. Bankole     const PetscInt wall = walls[w];
156*c5e9980aSAdeleke O. Bankole     IS             wall_is;
157*c5e9980aSAdeleke O. Bankole     PetscCall(DMLabelGetStratumIS(face_label, wall, &wall_is));
158*c5e9980aSAdeleke O. Bankole     if (wall_is) {  // There exist such points on this process
159*c5e9980aSAdeleke O. Bankole       PetscInt        num_points;
160*c5e9980aSAdeleke O. Bankole       const PetscInt *points;
161*c5e9980aSAdeleke O. Bankole       PetscCall(ISGetSize(wall_is, &num_points));
162*c5e9980aSAdeleke O. Bankole       PetscCall(ISGetIndices(wall_is, &points));
163*c5e9980aSAdeleke O. Bankole       for (PetscInt i = 0; i < num_points; i++) {
164*c5e9980aSAdeleke O. Bankole         const PetscInt           p = points[i];
165*c5e9980aSAdeleke O. Bankole         const StateConservative *r;
166*c5e9980aSAdeleke O. Bankole         PetscCall(DMPlexPointLocalRead(dm, p, g, &r));
167*c5e9980aSAdeleke O. Bankole         if (!r) continue;
168*c5e9980aSAdeleke O. Bankole         for (PetscInt j = 0; j < 3; j++) {
169*c5e9980aSAdeleke O. Bankole           reaction_force[w * dim + j] -= r->momentum[j];
170*c5e9980aSAdeleke O. Bankole         }
171*c5e9980aSAdeleke O. Bankole       }
172*c5e9980aSAdeleke O. Bankole       PetscCall(ISRestoreIndices(wall_is, &points));
173*c5e9980aSAdeleke O. Bankole     }
174*c5e9980aSAdeleke O. Bankole     PetscCall(ISDestroy(&wall_is));
175*c5e9980aSAdeleke O. Bankole   }
176*c5e9980aSAdeleke O. Bankole   PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, reaction_force, dim * num_walls, MPIU_SCALAR, MPI_SUM, comm));
177*c5e9980aSAdeleke O. Bankole   //  Restore Vectors
178*c5e9980aSAdeleke O. Bankole   PetscCall(VecRestoreArrayRead(G_loc, &g));
179*c5e9980aSAdeleke O. Bankole 
180*c5e9980aSAdeleke O. Bankole   PetscFunctionReturn(0);
181*c5e9980aSAdeleke O. Bankole }
182*c5e9980aSAdeleke O. Bankole 
183a515125bSLeila Ghaffari // Implicit time-stepper function setup
1842b916ea7SJeremy L Thompson PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) {
185a515125bSLeila Ghaffari   User               user = *(User *)user_data;
186a515125bSLeila Ghaffari   const PetscScalar *q, *q_dot;
187c2cb7fc8SJames Wright   PetscScalar       *g, dt;
188e2f84137SJeremy L Thompson   Vec                Q_loc = user->Q_loc, Q_dot_loc = user->Q_dot_loc, G_loc;
189a515125bSLeila Ghaffari   PetscMemType       q_mem_type, q_dot_mem_type, g_mem_type;
190a515125bSLeila Ghaffari   PetscFunctionBeginUser;
191a515125bSLeila Ghaffari 
192e2f84137SJeremy L Thompson   // Get local vectors
193*c5e9980aSAdeleke O. Bankole   PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
194e2f84137SJeremy L Thompson 
195e2f84137SJeremy L Thompson   // Update time dependent data
196c996854bSJames Wright   PetscCall(UpdateBoundaryValues(user, Q_loc, t));
197c2cb7fc8SJames Wright   PetscCall(UpdateContextLabel(&user->time, t, user->op_ifunction, user->phys->solution_time_label));
1982b916ea7SJeremy L Thompson   PetscCall(TSGetTimeStep(ts, &dt));
199c2cb7fc8SJames Wright   PetscCall(UpdateContextLabel(&user->dt, dt, user->op_ifunction, user->phys->timestep_size_label));
200a515125bSLeila Ghaffari 
201a515125bSLeila Ghaffari   // Global-to-local
20206108310SJames Wright   PetscCall(DMGlobalToLocalBegin(user->dm, Q, INSERT_VALUES, Q_loc));
20306108310SJames Wright   PetscCall(DMGlobalToLocalBegin(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
20406108310SJames Wright   PetscCall(DMGlobalToLocalEnd(user->dm, Q, INSERT_VALUES, Q_loc));
20506108310SJames Wright   PetscCall(DMGlobalToLocalEnd(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
206a515125bSLeila Ghaffari 
207a515125bSLeila Ghaffari   // Place PETSc vectors in CEED vectors
2082b916ea7SJeremy L Thompson   PetscCall(VecGetArrayReadAndMemType(Q_loc, &q, &q_mem_type));
2092b916ea7SJeremy L Thompson   PetscCall(VecGetArrayReadAndMemType(Q_dot_loc, &q_dot, &q_dot_mem_type));
2102b916ea7SJeremy L Thompson   PetscCall(VecGetArrayAndMemType(G_loc, &g, &g_mem_type));
2112b916ea7SJeremy L Thompson   CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER, (PetscScalar *)q);
2122b916ea7SJeremy L Thompson   CeedVectorSetArray(user->q_dot_ceed, MemTypeP2C(q_dot_mem_type), CEED_USE_POINTER, (PetscScalar *)q_dot);
213a515125bSLeila Ghaffari   CeedVectorSetArray(user->g_ceed, MemTypeP2C(g_mem_type), CEED_USE_POINTER, g);
214a515125bSLeila Ghaffari 
215a515125bSLeila Ghaffari   // Apply CEED operator
2162b916ea7SJeremy L Thompson   CeedOperatorApply(user->op_ifunction, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE);
217a515125bSLeila Ghaffari 
218a515125bSLeila Ghaffari   // Restore vectors
219a515125bSLeila Ghaffari   CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL);
220a515125bSLeila Ghaffari   CeedVectorTakeArray(user->q_dot_ceed, MemTypeP2C(q_dot_mem_type), NULL);
221a515125bSLeila Ghaffari   CeedVectorTakeArray(user->g_ceed, MemTypeP2C(g_mem_type), NULL);
2222b916ea7SJeremy L Thompson   PetscCall(VecRestoreArrayReadAndMemType(Q_loc, &q));
2232b916ea7SJeremy L Thompson   PetscCall(VecRestoreArrayReadAndMemType(Q_dot_loc, &q_dot));
2242b916ea7SJeremy L Thompson   PetscCall(VecRestoreArrayAndMemType(G_loc, &g));
225a515125bSLeila Ghaffari 
226a515125bSLeila Ghaffari   // Local-to-Global
2272b916ea7SJeremy L Thompson   PetscCall(VecZeroEntries(G));
2282b916ea7SJeremy L Thompson   PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G));
229a515125bSLeila Ghaffari 
230a515125bSLeila Ghaffari   // Restore vectors
231*c5e9980aSAdeleke O. Bankole   PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
232a515125bSLeila Ghaffari 
233a515125bSLeila Ghaffari   PetscFunctionReturn(0);
234a515125bSLeila Ghaffari }
235a515125bSLeila Ghaffari 
236f0b65372SJed Brown static PetscErrorCode MatMult_NS_IJacobian(Mat J, Vec Q, Vec G) {
237f0b65372SJed Brown   User               user;
238f0b65372SJed Brown   const PetscScalar *q;
239f0b65372SJed Brown   PetscScalar       *g;
240f0b65372SJed Brown   PetscMemType       q_mem_type, g_mem_type;
241f0b65372SJed Brown   PetscFunctionBeginUser;
2422b916ea7SJeremy L Thompson   PetscCall(MatShellGetContext(J, &user));
243e2f84137SJeremy L Thompson   Vec Q_loc = user->Q_dot_loc,  // Note - Q_dot_loc has zero BCs
244e2f84137SJeremy L Thompson       G_loc;
245e2f84137SJeremy L Thompson 
246f0b65372SJed Brown   // Get local vectors
2472b916ea7SJeremy L Thompson   PetscCall(DMGetLocalVector(user->dm, &G_loc));
248f0b65372SJed Brown 
249f0b65372SJed Brown   // Global-to-local
2502b916ea7SJeremy L Thompson   PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
251f0b65372SJed Brown 
252f0b65372SJed Brown   // Place PETSc vectors in CEED vectors
2532b916ea7SJeremy L Thompson   PetscCall(VecGetArrayReadAndMemType(Q_loc, &q, &q_mem_type));
2542b916ea7SJeremy L Thompson   PetscCall(VecGetArrayAndMemType(G_loc, &g, &g_mem_type));
2552b916ea7SJeremy L Thompson   CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER, (PetscScalar *)q);
256f0b65372SJed Brown   CeedVectorSetArray(user->g_ceed, MemTypeP2C(g_mem_type), CEED_USE_POINTER, g);
257f0b65372SJed Brown 
258f0b65372SJed Brown   // Apply CEED operator
2592b916ea7SJeremy L Thompson   CeedOperatorApply(user->op_ijacobian, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE);
260f0b65372SJed Brown 
261f0b65372SJed Brown   // Restore vectors
262f0b65372SJed Brown   CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL);
263f0b65372SJed Brown   CeedVectorTakeArray(user->g_ceed, MemTypeP2C(g_mem_type), NULL);
2642b916ea7SJeremy L Thompson   PetscCall(VecRestoreArrayReadAndMemType(Q_loc, &q));
2652b916ea7SJeremy L Thompson   PetscCall(VecRestoreArrayAndMemType(G_loc, &g));
266f0b65372SJed Brown 
267f0b65372SJed Brown   // Local-to-Global
2682b916ea7SJeremy L Thompson   PetscCall(VecZeroEntries(G));
2692b916ea7SJeremy L Thompson   PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G));
270f0b65372SJed Brown 
271f0b65372SJed Brown   // Restore vectors
2722b916ea7SJeremy L Thompson   PetscCall(DMRestoreLocalVector(user->dm, &G_loc));
273f0b65372SJed Brown   PetscFunctionReturn(0);
274f0b65372SJed Brown }
275f0b65372SJed Brown 
276f0b65372SJed Brown PetscErrorCode MatGetDiagonal_NS_IJacobian(Mat A, Vec D) {
277f0b65372SJed Brown   User         user;
278f0b65372SJed Brown   Vec          D_loc;
279f0b65372SJed Brown   PetscScalar *d;
280f0b65372SJed Brown   PetscMemType mem_type;
281f0b65372SJed Brown 
282f0b65372SJed Brown   PetscFunctionBeginUser;
283f0b65372SJed Brown   MatShellGetContext(A, &user);
284f0b65372SJed Brown   PetscCall(DMGetLocalVector(user->dm, &D_loc));
285f0b65372SJed Brown   PetscCall(VecGetArrayAndMemType(D_loc, &d, &mem_type));
286f0b65372SJed Brown   CeedVectorSetArray(user->g_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, d);
2872b916ea7SJeremy L Thompson   CeedOperatorLinearAssembleDiagonal(user->op_ijacobian, user->g_ceed, CEED_REQUEST_IMMEDIATE);
288f0b65372SJed Brown   CeedVectorTakeArray(user->g_ceed, MemTypeP2C(mem_type), NULL);
289f0b65372SJed Brown   PetscCall(VecRestoreArrayAndMemType(D_loc, &d));
290f0b65372SJed Brown   PetscCall(VecZeroEntries(D));
291f0b65372SJed Brown   PetscCall(DMLocalToGlobal(user->dm, D_loc, ADD_VALUES, D));
292f0b65372SJed Brown   PetscCall(DMRestoreLocalVector(user->dm, &D_loc));
293f0b65372SJed Brown   VecViewFromOptions(D, NULL, "-diag_vec_view");
294f0b65372SJed Brown   PetscFunctionReturn(0);
295f0b65372SJed Brown }
296f0b65372SJed Brown 
2972b916ea7SJeremy L Thompson static PetscErrorCode FormPreallocation(User user, PetscBool pbdiagonal, Mat J, CeedVector *coo_values) {
298b107fddaSJed Brown   PetscCount ncoo;
299b107fddaSJed Brown   PetscInt  *rows, *cols;
300b107fddaSJed Brown 
301b107fddaSJed Brown   PetscFunctionBeginUser;
302b107fddaSJed Brown   if (pbdiagonal) {
303b107fddaSJed Brown     CeedSize l_size;
304b107fddaSJed Brown     CeedOperatorGetActiveVectorLengths(user->op_ijacobian, &l_size, NULL);
305b107fddaSJed Brown     ncoo = l_size * 5;
306b107fddaSJed Brown     rows = malloc(ncoo * sizeof(rows[0]));
307b107fddaSJed Brown     cols = malloc(ncoo * sizeof(cols[0]));
308b107fddaSJed Brown     for (PetscCount n = 0; n < l_size / 5; n++) {
309b107fddaSJed Brown       for (PetscInt i = 0; i < 5; i++) {
310b107fddaSJed Brown         for (PetscInt j = 0; j < 5; j++) {
311b107fddaSJed Brown           rows[(n * 5 + i) * 5 + j] = n * 5 + i;
312b107fddaSJed Brown           cols[(n * 5 + i) * 5 + j] = n * 5 + j;
313b107fddaSJed Brown         }
314b107fddaSJed Brown       }
315b107fddaSJed Brown     }
316b107fddaSJed Brown   } else {
3172b916ea7SJeremy L Thompson     PetscCall(CeedOperatorLinearAssembleSymbolic(user->op_ijacobian, &ncoo, &rows, &cols));
318b107fddaSJed Brown   }
319b107fddaSJed Brown   PetscCall(MatSetPreallocationCOOLocal(J, ncoo, rows, cols));
320b107fddaSJed Brown   free(rows);
321b107fddaSJed Brown   free(cols);
322b107fddaSJed Brown   CeedVectorCreate(user->ceed, ncoo, coo_values);
323b107fddaSJed Brown   PetscFunctionReturn(0);
324b107fddaSJed Brown }
325b107fddaSJed Brown 
3262b916ea7SJeremy L Thompson static PetscErrorCode FormSetValues(User user, PetscBool pbdiagonal, Mat J, CeedVector coo_values) {
327b107fddaSJed Brown   CeedMemType        mem_type = CEED_MEM_HOST;
328b107fddaSJed Brown   const PetscScalar *values;
329b107fddaSJed Brown   MatType            mat_type;
330b107fddaSJed Brown 
331b107fddaSJed Brown   PetscFunctionBeginUser;
332b107fddaSJed Brown   PetscCall(MatGetType(J, &mat_type));
3332b916ea7SJeremy L Thompson   if (strstr(mat_type, "kokkos") || strstr(mat_type, "cusparse")) mem_type = CEED_MEM_DEVICE;
334b107fddaSJed Brown   if (user->app_ctx->pmat_pbdiagonal) {
3352b916ea7SJeremy L Thompson     CeedOperatorLinearAssemblePointBlockDiagonal(user->op_ijacobian, coo_values, CEED_REQUEST_IMMEDIATE);
336b107fddaSJed Brown   } else {
337b107fddaSJed Brown     CeedOperatorLinearAssemble(user->op_ijacobian, coo_values);
338b107fddaSJed Brown   }
339b107fddaSJed Brown   CeedVectorGetArrayRead(coo_values, mem_type, &values);
340b107fddaSJed Brown   PetscCall(MatSetValuesCOO(J, values, INSERT_VALUES));
341b107fddaSJed Brown   CeedVectorRestoreArrayRead(coo_values, &values);
342b107fddaSJed Brown   PetscFunctionReturn(0);
343b107fddaSJed Brown }
344b107fddaSJed Brown 
3452b916ea7SJeremy L Thompson PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, PetscReal shift, Mat J, Mat J_pre, void *user_data) {
346f0b65372SJed Brown   User      user = *(User *)user_data;
34704855949SJed Brown   PetscBool J_is_shell, J_is_mffd, J_pre_is_shell;
348f0b65372SJed Brown   PetscFunctionBeginUser;
3492b916ea7SJeremy L Thompson   if (user->phys->ijacobian_time_shift_label) CeedOperatorContextSetDouble(user->op_ijacobian, user->phys->ijacobian_time_shift_label, &shift);
35004855949SJed Brown   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATMFFD, &J_is_mffd));
351f0b65372SJed Brown   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATSHELL, &J_is_shell));
3522b916ea7SJeremy L Thompson   PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATSHELL, &J_pre_is_shell));
353f0b65372SJed Brown   if (!user->matrices_set_up) {
354f0b65372SJed Brown     if (J_is_shell) {
355f0b65372SJed Brown       PetscCall(MatShellSetContext(J, user));
3562b916ea7SJeremy L Thompson       PetscCall(MatShellSetOperation(J, MATOP_MULT, (void (*)(void))MatMult_NS_IJacobian));
3572b916ea7SJeremy L Thompson       PetscCall(MatShellSetOperation(J, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_NS_IJacobian));
358f0b65372SJed Brown       PetscCall(MatSetUp(J));
359f0b65372SJed Brown     }
360f0b65372SJed Brown     if (!J_pre_is_shell) {
3612b916ea7SJeremy L Thompson       PetscCall(FormPreallocation(user, user->app_ctx->pmat_pbdiagonal, J_pre, &user->coo_values_pmat));
362b107fddaSJed Brown     }
36304855949SJed Brown     if (J != J_pre && !J_is_shell && !J_is_mffd) {
364b107fddaSJed Brown       PetscCall(FormPreallocation(user, PETSC_FALSE, J, &user->coo_values_amat));
365b107fddaSJed Brown     }
366f0b65372SJed Brown     user->matrices_set_up = true;
367f0b65372SJed Brown   }
368f0b65372SJed Brown   if (!J_pre_is_shell) {
3692b916ea7SJeremy L Thompson     PetscCall(FormSetValues(user, user->app_ctx->pmat_pbdiagonal, J_pre, user->coo_values_pmat));
370f0b65372SJed Brown   }
37104855949SJed Brown   if (user->coo_values_amat) {
37204855949SJed Brown     PetscCall(FormSetValues(user, PETSC_FALSE, J, user->coo_values_amat));
37304855949SJed Brown   } else if (J_is_mffd) {
37404855949SJed Brown     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
37504855949SJed Brown     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
37604855949SJed Brown   }
377f0b65372SJed Brown   PetscFunctionReturn(0);
378f0b65372SJed Brown }
379f0b65372SJed Brown 
3802b916ea7SJeremy L Thompson PetscErrorCode WriteOutput(User user, Vec Q, PetscInt step_no, PetscScalar time) {
381a515125bSLeila Ghaffari   Vec         Q_loc;
382a515125bSLeila Ghaffari   char        file_path[PETSC_MAX_PATH_LEN];
383a515125bSLeila Ghaffari   PetscViewer viewer;
384a515125bSLeila Ghaffari   PetscFunctionBeginUser;
385a515125bSLeila Ghaffari 
386852e5969SJed Brown   if (user->app_ctx->checkpoint_vtk) {
387a515125bSLeila Ghaffari     // Set up output
3887538d537SJames Wright     PetscCall(DMGetLocalVector(user->dm, &Q_loc));
3897538d537SJames Wright     PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec"));
3907538d537SJames Wright     PetscCall(VecZeroEntries(Q_loc));
3917538d537SJames Wright     PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
392a515125bSLeila Ghaffari 
393a515125bSLeila Ghaffari     // Output
394852e5969SJed Brown     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
3957538d537SJames Wright 
3962b916ea7SJeremy L Thompson     PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path, FILE_MODE_WRITE, &viewer));
3977538d537SJames Wright     PetscCall(VecView(Q_loc, viewer));
3987538d537SJames Wright     PetscCall(PetscViewerDestroy(&viewer));
399a515125bSLeila Ghaffari     if (user->dm_viz) {
400a515125bSLeila Ghaffari       Vec         Q_refined, Q_refined_loc;
401a515125bSLeila Ghaffari       char        file_path_refined[PETSC_MAX_PATH_LEN];
402a515125bSLeila Ghaffari       PetscViewer viewer_refined;
403a515125bSLeila Ghaffari 
4047538d537SJames Wright       PetscCall(DMGetGlobalVector(user->dm_viz, &Q_refined));
4057538d537SJames Wright       PetscCall(DMGetLocalVector(user->dm_viz, &Q_refined_loc));
4067538d537SJames Wright       PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined"));
4077538d537SJames Wright 
4087538d537SJames Wright       PetscCall(MatInterpolate(user->interp_viz, Q, Q_refined));
4097538d537SJames Wright       PetscCall(VecZeroEntries(Q_refined_loc));
4102b916ea7SJeremy L Thompson       PetscCall(DMGlobalToLocal(user->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc));
4117538d537SJames Wright 
412852e5969SJed Brown       PetscCall(
413852e5969SJed Brown           PetscSNPrintf(file_path_refined, sizeof file_path_refined, "%s/nsrefined-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
4147538d537SJames Wright 
4152b916ea7SJeremy L Thompson       PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined), file_path_refined, FILE_MODE_WRITE, &viewer_refined));
4167538d537SJames Wright       PetscCall(VecView(Q_refined_loc, viewer_refined));
4177538d537SJames Wright       PetscCall(DMRestoreLocalVector(user->dm_viz, &Q_refined_loc));
4187538d537SJames Wright       PetscCall(DMRestoreGlobalVector(user->dm_viz, &Q_refined));
4197538d537SJames Wright       PetscCall(PetscViewerDestroy(&viewer_refined));
420a515125bSLeila Ghaffari     }
4217538d537SJames Wright     PetscCall(DMRestoreLocalVector(user->dm, &Q_loc));
422852e5969SJed Brown   }
423a515125bSLeila Ghaffari 
424a515125bSLeila Ghaffari   // Save data in a binary file for continuation of simulations
42591a36801SJames Wright   if (user->app_ctx->add_stepnum2bin) {
426852e5969SJed Brown     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", user->app_ctx->output_dir, step_no));
42791a36801SJames Wright   } else {
4282b916ea7SJeremy L Thompson     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", user->app_ctx->output_dir));
42991a36801SJames Wright   }
4302b916ea7SJeremy L Thompson   PetscCall(PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, &viewer));
4317538d537SJames Wright 
4329293eaa1SJed Brown   PetscInt token = FLUIDS_FILE_TOKEN;
4339293eaa1SJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &token, 1, PETSC_INT));
4349293eaa1SJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &step_no, 1, PETSC_INT));
4359293eaa1SJed Brown   time /= user->units->second;  // Dimensionalize time back
4369293eaa1SJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL));
4377538d537SJames Wright   PetscCall(VecView(Q, viewer));
4387538d537SJames Wright   PetscCall(PetscViewerDestroy(&viewer));
4397538d537SJames Wright   PetscFunctionReturn(0);
4407538d537SJames Wright }
4417538d537SJames Wright 
442*c5e9980aSAdeleke O. Bankole // CSV Monitor
443*c5e9980aSAdeleke O. Bankole PetscErrorCode TSMonitor_WallForce(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
444*c5e9980aSAdeleke O. Bankole   User              user = ctx;
445*c5e9980aSAdeleke O. Bankole   Vec               G_loc;
446*c5e9980aSAdeleke O. Bankole   PetscInt          num_wall = user->app_ctx->wall_forces.num_wall, dim = 3;
447*c5e9980aSAdeleke O. Bankole   const PetscInt   *walls  = user->app_ctx->wall_forces.walls;
448*c5e9980aSAdeleke O. Bankole   PetscViewer       viewer = user->app_ctx->wall_forces.viewer;
449*c5e9980aSAdeleke O. Bankole   PetscViewerFormat format = user->app_ctx->wall_forces.viewer_format;
450*c5e9980aSAdeleke O. Bankole   PetscScalar      *reaction_force;
451*c5e9980aSAdeleke O. Bankole   PetscBool         iascii;
452*c5e9980aSAdeleke O. Bankole 
453*c5e9980aSAdeleke O. Bankole   PetscFunctionBeginUser;
454*c5e9980aSAdeleke O. Bankole   if (!viewer) PetscFunctionReturn(0);
455*c5e9980aSAdeleke O. Bankole   PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
456*c5e9980aSAdeleke O. Bankole   PetscCall(PetscMalloc1(num_wall * dim, &reaction_force));
457*c5e9980aSAdeleke O. Bankole   PetscCall(Surface_Forces_NS(user->dm, G_loc, num_wall, walls, reaction_force));
458*c5e9980aSAdeleke O. Bankole   PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
459*c5e9980aSAdeleke O. Bankole 
460*c5e9980aSAdeleke O. Bankole   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
461*c5e9980aSAdeleke O. Bankole 
462*c5e9980aSAdeleke O. Bankole   if (iascii) {
463*c5e9980aSAdeleke O. Bankole     if (format == PETSC_VIEWER_ASCII_CSV && !user->app_ctx->wall_forces.header_written) {
464*c5e9980aSAdeleke O. Bankole       PetscCall(PetscViewerASCIIPrintf(viewer, "Step,Time,Wall,ForceX,ForceY,ForceZ\n"));
465*c5e9980aSAdeleke O. Bankole       user->app_ctx->wall_forces.header_written = PETSC_TRUE;
466*c5e9980aSAdeleke O. Bankole     }
467*c5e9980aSAdeleke O. Bankole     for (PetscInt w = 0; w < num_wall; w++) {
468*c5e9980aSAdeleke O. Bankole       PetscInt wall = walls[w];
469*c5e9980aSAdeleke O. Bankole       if (format == PETSC_VIEWER_ASCII_CSV) {
470*c5e9980aSAdeleke O. Bankole         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT ",%g,%" PetscInt_FMT ",%g,%g,%g\n", step_no, time, wall,
471*c5e9980aSAdeleke O. Bankole                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
472*c5e9980aSAdeleke O. Bankole 
473*c5e9980aSAdeleke O. Bankole       } else {
474*c5e9980aSAdeleke O. Bankole         PetscCall(PetscViewerASCIIPrintf(viewer, "Wall %" PetscInt_FMT " Forces: Force_x = %12g, Force_y = %12g, Force_z = %12g\n", wall,
475*c5e9980aSAdeleke O. Bankole                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
476*c5e9980aSAdeleke O. Bankole       }
477*c5e9980aSAdeleke O. Bankole     }
478*c5e9980aSAdeleke O. Bankole   }
479*c5e9980aSAdeleke O. Bankole   PetscCall(PetscFree(reaction_force));
480*c5e9980aSAdeleke O. Bankole   PetscFunctionReturn(0);
481*c5e9980aSAdeleke O. Bankole }
482*c5e9980aSAdeleke O. Bankole 
4837538d537SJames Wright // User provided TS Monitor
4842b916ea7SJeremy L Thompson PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
4857538d537SJames Wright   User user = ctx;
4867538d537SJames Wright   PetscFunctionBeginUser;
4877538d537SJames Wright 
488852e5969SJed Brown   // Print every 'checkpoint_interval' steps
489c539088bSJames Wright   if (user->app_ctx->checkpoint_interval <= 0 || step_no % user->app_ctx->checkpoint_interval != 0 ||
490c539088bSJames Wright       (user->app_ctx->cont_steps == step_no && step_no != 0))
491c539088bSJames Wright     PetscFunctionReturn(0);
4927538d537SJames Wright 
4937538d537SJames Wright   PetscCall(WriteOutput(user, Q, step_no, time));
494a515125bSLeila Ghaffari 
495a515125bSLeila Ghaffari   PetscFunctionReturn(0);
496a515125bSLeila Ghaffari }
497a515125bSLeila Ghaffari 
498a515125bSLeila Ghaffari // TS: Create, setup, and solve
4992b916ea7SJeremy L Thompson PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, Vec *Q, PetscScalar *f_time, TS *ts) {
500a515125bSLeila Ghaffari   MPI_Comm    comm = user->comm;
501a515125bSLeila Ghaffari   TSAdapt     adapt;
502a515125bSLeila Ghaffari   PetscScalar final_time;
503a515125bSLeila Ghaffari   PetscFunctionBeginUser;
504a515125bSLeila Ghaffari 
5052b916ea7SJeremy L Thompson   PetscCall(TSCreate(comm, ts));
5062b916ea7SJeremy L Thompson   PetscCall(TSSetDM(*ts, dm));
507a515125bSLeila Ghaffari   if (phys->implicit) {
5082b916ea7SJeremy L Thompson     PetscCall(TSSetType(*ts, TSBDF));
509a515125bSLeila Ghaffari     if (user->op_ifunction) {
5102b916ea7SJeremy L Thompson       PetscCall(TSSetIFunction(*ts, NULL, IFunction_NS, &user));
511a515125bSLeila Ghaffari     } else {  // Implicit integrators can fall back to using an RHSFunction
5122b916ea7SJeremy L Thompson       PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
513a515125bSLeila Ghaffari     }
514f0b65372SJed Brown     if (user->op_ijacobian) {
5152b916ea7SJeremy L Thompson       PetscCall(DMTSSetIJacobian(dm, FormIJacobian_NS, &user));
516b107fddaSJed Brown       if (app_ctx->amat_type) {
517b107fddaSJed Brown         Mat Pmat, Amat;
5182b916ea7SJeremy L Thompson         PetscCall(DMCreateMatrix(dm, &Pmat));
5192b916ea7SJeremy L Thompson         PetscCall(DMSetMatType(dm, app_ctx->amat_type));
5202b916ea7SJeremy L Thompson         PetscCall(DMCreateMatrix(dm, &Amat));
5212b916ea7SJeremy L Thompson         PetscCall(TSSetIJacobian(*ts, Amat, Pmat, NULL, NULL));
5222b916ea7SJeremy L Thompson         PetscCall(MatDestroy(&Amat));
5232b916ea7SJeremy L Thompson         PetscCall(MatDestroy(&Pmat));
524b107fddaSJed Brown       }
525f0b65372SJed Brown     }
526a515125bSLeila Ghaffari   } else {
5272b916ea7SJeremy L Thompson     if (!user->op_rhs) SETERRQ(comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction");
5282b916ea7SJeremy L Thompson     PetscCall(TSSetType(*ts, TSRK));
5292b916ea7SJeremy L Thompson     PetscCall(TSRKSetType(*ts, TSRK5F));
5302b916ea7SJeremy L Thompson     PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
531a515125bSLeila Ghaffari   }
5322b916ea7SJeremy L Thompson   PetscCall(TSSetMaxTime(*ts, 500. * user->units->second));
5332b916ea7SJeremy L Thompson   PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER));
534f0784ed3SJed Brown   PetscCall(TSSetErrorIfStepFails(*ts, PETSC_FALSE));
5352b916ea7SJeremy L Thompson   PetscCall(TSSetTimeStep(*ts, 1.e-2 * user->units->second));
5360e1e9333SJames Wright   if (app_ctx->test_type != TESTTYPE_NONE) {
5372b916ea7SJeremy L Thompson     PetscCall(TSSetMaxSteps(*ts, 10));
5382b916ea7SJeremy L Thompson   }
5392b916ea7SJeremy L Thompson   PetscCall(TSGetAdapt(*ts, &adapt));
5402b916ea7SJeremy L Thompson   PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * user->units->second, 1.e2 * user->units->second));
5412b916ea7SJeremy L Thompson   PetscCall(TSSetFromOptions(*ts));
542c996854bSJames Wright   user->time = user->time_bc_set = -1.0;  // require all BCs and ctx to be updated
543e2f84137SJeremy L Thompson   user->dt                       = -1.0;
544a515125bSLeila Ghaffari   if (!app_ctx->cont_steps) {  // print initial condition
5450e1e9333SJames Wright     if (app_ctx->test_type == TESTTYPE_NONE) {
5462b916ea7SJeremy L Thompson       PetscCall(TSMonitor_NS(*ts, 0, 0., *Q, user));
547a515125bSLeila Ghaffari     }
548a515125bSLeila Ghaffari   } else {  // continue from time of last output
549a515125bSLeila Ghaffari     PetscInt    count;
550a515125bSLeila Ghaffari     PetscViewer viewer;
5512b916ea7SJeremy L Thompson 
5529293eaa1SJed Brown     if (app_ctx->cont_time <= 0) {  // Legacy files did not include step number and time
5532b916ea7SJeremy L Thompson       PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_time_file, FILE_MODE_READ, &viewer));
5549293eaa1SJed Brown       PetscCall(PetscViewerBinaryRead(viewer, &app_ctx->cont_time, 1, &count, PETSC_REAL));
5552b916ea7SJeremy L Thompson       PetscCall(PetscViewerDestroy(&viewer));
5569293eaa1SJed Brown       PetscCheck(app_ctx->cont_steps != -1, comm, PETSC_ERR_ARG_INCOMP,
5579293eaa1SJed Brown                  "-continue step number not specified, but checkpoint file does not contain a step number (likely written by older code version)");
5589293eaa1SJed Brown     }
5599293eaa1SJed Brown     PetscCall(TSSetTime(*ts, app_ctx->cont_time * user->units->second));
56074a6f4ddSJed Brown     PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps));
561a515125bSLeila Ghaffari   }
5620e1e9333SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) {
5632b916ea7SJeremy L Thompson     PetscCall(TSMonitorSet(*ts, TSMonitor_NS, user, NULL));
5640e1e9333SJames Wright   }
565*c5e9980aSAdeleke O. Bankole   if (app_ctx->wall_forces.viewer) {
566*c5e9980aSAdeleke O. Bankole     PetscCall(TSMonitorSet(*ts, TSMonitor_WallForce, user, NULL));
567*c5e9980aSAdeleke O. Bankole   }
568c931fa59SJames Wright   if (app_ctx->turb_spanstats_enable) {
569b0488d1fSJames Wright     PetscCall(TSMonitorSet(*ts, TSMonitor_Statistics, user, NULL));
570b8daee98SJames Wright     CeedScalar previous_time = app_ctx->cont_time * user->units->second;
571b8daee98SJames Wright     CeedOperatorContextSetDouble(user->spanstats.op_stats_collect, user->spanstats.previous_time_label, &previous_time);
572b0488d1fSJames Wright   }
573a515125bSLeila Ghaffari 
574a515125bSLeila Ghaffari   // Solve
57574a6f4ddSJed Brown   PetscReal start_time;
57674a6f4ddSJed Brown   PetscInt  start_step;
5772b916ea7SJeremy L Thompson   PetscCall(TSGetTime(*ts, &start_time));
57874a6f4ddSJed Brown   PetscCall(TSGetStepNumber(*ts, &start_step));
57991982731SJeremy L Thompson 
58091982731SJeremy L Thompson   PetscPreLoadBegin(PETSC_FALSE, "Fluids Solve");
58191982731SJeremy L Thompson   PetscCall(TSSetTime(*ts, start_time));
58274a6f4ddSJed Brown   PetscCall(TSSetStepNumber(*ts, start_step));
58391982731SJeremy L Thompson   if (PetscPreLoadingOn) {
58491982731SJeremy L Thompson     // LCOV_EXCL_START
58591982731SJeremy L Thompson     SNES      snes;
58691982731SJeremy L Thompson     Vec       Q_preload;
58791982731SJeremy L Thompson     PetscReal rtol;
58891982731SJeremy L Thompson     PetscCall(VecDuplicate(*Q, &Q_preload));
58991982731SJeremy L Thompson     PetscCall(VecCopy(*Q, Q_preload));
59091982731SJeremy L Thompson     PetscCall(TSGetSNES(*ts, &snes));
59191982731SJeremy L Thompson     PetscCall(SNESGetTolerances(snes, NULL, &rtol, NULL, NULL, NULL));
5922b916ea7SJeremy L Thompson     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, .99, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
59391982731SJeremy L Thompson     PetscCall(TSSetSolution(*ts, *Q));
59491982731SJeremy L Thompson     PetscCall(TSStep(*ts));
5952b916ea7SJeremy L Thompson     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, rtol, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
59691982731SJeremy L Thompson     PetscCall(VecDestroy(&Q_preload));
59791982731SJeremy L Thompson     // LCOV_EXCL_STOP
59891982731SJeremy L Thompson   } else {
5992b916ea7SJeremy L Thompson     PetscCall(PetscBarrier((PetscObject)*ts));
6002b916ea7SJeremy L Thompson     PetscCall(TSSolve(*ts, *Q));
60191982731SJeremy L Thompson   }
60291982731SJeremy L Thompson   PetscPreLoadEnd();
60391982731SJeremy L Thompson 
60491982731SJeremy L Thompson   PetscCall(TSGetSolveTime(*ts, &final_time));
605a515125bSLeila Ghaffari   *f_time = final_time;
60691982731SJeremy L Thompson 
6070e1e9333SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) {
6087538d537SJames Wright     PetscInt step_no;
6097538d537SJames Wright     PetscCall(TSGetStepNumber(*ts, &step_no));
610b0488d1fSJames Wright     if (user->app_ctx->checkpoint_interval > 0 || user->app_ctx->checkpoint_interval == -1) {
6117538d537SJames Wright       PetscCall(WriteOutput(user, *Q, step_no, final_time));
6127538d537SJames Wright     }
6137538d537SJames Wright 
61491982731SJeremy L Thompson     PetscLogEvent stage_id;
61591982731SJeremy L Thompson     PetscStageLog stage_log;
61691982731SJeremy L Thompson 
61791982731SJeremy L Thompson     PetscCall(PetscLogStageGetId("Fluids Solve", &stage_id));
61891982731SJeremy L Thompson     PetscCall(PetscLogGetStageLog(&stage_log));
6192b916ea7SJeremy L Thompson     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time taken for solution (sec): %g\n", stage_log->stageInfo[stage_id].perfInfo.time));
620a515125bSLeila Ghaffari   }
621a515125bSLeila Ghaffari   PetscFunctionReturn(0);
622a515125bSLeila Ghaffari }
623