xref: /libCEED/examples/fluids/src/setupts.c (revision 17b0d5c66e710b5b55588e66a2830e7df114beeb)
13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
377841947SLeila Ghaffari //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
577841947SLeila Ghaffari //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
777841947SLeila Ghaffari 
877841947SLeila Ghaffari /// @file
977841947SLeila Ghaffari /// Time-stepping functions for Navier-Stokes example using PETSc
1077841947SLeila Ghaffari 
1177841947SLeila Ghaffari #include "../navierstokes.h"
12ca69d878SAdeleke O. Bankole #include "../qfunctions/newtonian_state.h"
1377841947SLeila Ghaffari 
1477841947SLeila Ghaffari // Compute mass matrix for explicit scheme
152b730f8bSJeremy L Thompson PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm, CeedData ceed_data, Vec M) {
1677841947SLeila Ghaffari   Vec           M_loc;
1777841947SLeila Ghaffari   CeedQFunction qf_mass;
1877841947SLeila Ghaffari   CeedOperator  op_mass;
1977841947SLeila Ghaffari   CeedVector    m_ceed, ones_vec;
2077841947SLeila Ghaffari   CeedInt       num_comp_q, q_data_size;
2177841947SLeila Ghaffari   PetscFunctionBeginUser;
2277841947SLeila Ghaffari 
2377841947SLeila Ghaffari   // CEED Restriction
2477841947SLeila Ghaffari   CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q);
2577841947SLeila Ghaffari   CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &q_data_size);
2677841947SLeila Ghaffari   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &m_ceed, NULL);
2777841947SLeila Ghaffari   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &ones_vec, NULL);
2877841947SLeila Ghaffari   CeedVectorSetValue(ones_vec, 1.0);
2977841947SLeila Ghaffari 
3077841947SLeila Ghaffari   // CEED QFunction
31ef080ff9SJames Wright   PetscCall(CreateMassQFunction(ceed, num_comp_q, q_data_size, &qf_mass));
3277841947SLeila Ghaffari 
3377841947SLeila Ghaffari   // CEED Operator
3477841947SLeila Ghaffari   CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass);
35ef080ff9SJames Wright   CeedOperatorSetField(op_mass, "u", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
362b730f8bSJeremy L Thompson   CeedOperatorSetField(op_mass, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data);
372b730f8bSJeremy L Thompson   CeedOperatorSetField(op_mass, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
3877841947SLeila Ghaffari 
3977841947SLeila Ghaffari   // Place PETSc vector in CEED vector
4077841947SLeila Ghaffari   PetscMemType m_mem_type;
412b730f8bSJeremy L Thompson   PetscCall(DMGetLocalVector(dm, &M_loc));
42c798d55aSJames Wright   PetscCall(VecP2C(M_loc, &m_mem_type, m_ceed));
4377841947SLeila Ghaffari 
4477841947SLeila Ghaffari   // Apply CEED Operator
4577841947SLeila Ghaffari   CeedOperatorApply(op_mass, ones_vec, m_ceed, CEED_REQUEST_IMMEDIATE);
4677841947SLeila Ghaffari 
4777841947SLeila Ghaffari   // Restore vectors
48c798d55aSJames Wright   PetscCall(VecC2P(m_ceed, m_mem_type, M_loc));
4977841947SLeila Ghaffari 
5077841947SLeila Ghaffari   // Local-to-Global
512b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(M));
522b730f8bSJeremy L Thompson   PetscCall(DMLocalToGlobal(dm, M_loc, ADD_VALUES, M));
532b730f8bSJeremy L Thompson   PetscCall(DMRestoreLocalVector(dm, &M_loc));
5477841947SLeila Ghaffari 
5577841947SLeila Ghaffari   // Invert diagonally lumped mass vector for RHS function
562b730f8bSJeremy L Thompson   PetscCall(VecReciprocal(M));
5777841947SLeila Ghaffari 
5877841947SLeila Ghaffari   // Cleanup
5977841947SLeila Ghaffari   CeedVectorDestroy(&ones_vec);
6077841947SLeila Ghaffari   CeedVectorDestroy(&m_ceed);
6177841947SLeila Ghaffari   CeedQFunctionDestroy(&qf_mass);
6277841947SLeila Ghaffari   CeedOperatorDestroy(&op_mass);
6377841947SLeila Ghaffari 
6477841947SLeila Ghaffari   PetscFunctionReturn(0);
6577841947SLeila Ghaffari }
6677841947SLeila Ghaffari 
67a0b9a424SJames Wright // Insert Boundary values if it's a new time
68a0b9a424SJames Wright PetscErrorCode UpdateBoundaryValues(User user, Vec Q_loc, PetscReal t) {
69a0b9a424SJames Wright   PetscFunctionBeginUser;
70a0b9a424SJames Wright   if (user->time_bc_set != t) {
71a0b9a424SJames Wright     PetscCall(DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL));
72a0b9a424SJames Wright     user->time_bc_set = t;
73a0b9a424SJames Wright   }
74a0b9a424SJames Wright   PetscFunctionReturn(0);
75a0b9a424SJames Wright }
76a0b9a424SJames Wright 
77e42c1df6SJames Wright PetscErrorCode UpdateContextLabel(PetscScalar *previous_value, PetscScalar update_value, CeedOperator op, CeedContextFieldLabel label) {
78e42c1df6SJames Wright   PetscFunctionBeginUser;
79e42c1df6SJames Wright 
80e42c1df6SJames Wright   if (*previous_value != update_value) {
81e42c1df6SJames Wright     if (label) {
82*17b0d5c6SJeremy L Thompson       CeedOperatorSetContextDouble(op, label, &update_value);
83e42c1df6SJames Wright     }
84e42c1df6SJames Wright     *previous_value = update_value;
85e42c1df6SJames Wright   }
86e42c1df6SJames Wright   PetscFunctionReturn(0);
87e42c1df6SJames Wright }
88e42c1df6SJames Wright 
8977841947SLeila Ghaffari // RHS (Explicit time-stepper) function setup
9077841947SLeila Ghaffari //   This is the RHS of the ODE, given as u_t = G(t,u)
9177841947SLeila Ghaffari //   This function takes in a state vector Q and writes into G
9277841947SLeila Ghaffari PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) {
9377841947SLeila Ghaffari   User         user = *(User *)user_data;
94c798d55aSJames Wright   PetscScalar  dt;
955e82a6e1SJeremy L Thompson   Vec          Q_loc = user->Q_loc, G_loc;
9677841947SLeila Ghaffari   PetscMemType q_mem_type, g_mem_type;
9777841947SLeila Ghaffari   PetscFunctionBeginUser;
9877841947SLeila Ghaffari 
995e82a6e1SJeremy L Thompson   // Get local vector
1002b730f8bSJeremy L Thompson   PetscCall(DMGetLocalVector(user->dm, &G_loc));
1015e82a6e1SJeremy L Thompson 
1025e82a6e1SJeremy L Thompson   // Update time dependent data
103a0b9a424SJames Wright   PetscCall(UpdateBoundaryValues(user, Q_loc, t));
104e42c1df6SJames Wright   PetscCall(UpdateContextLabel(&user->time, t, user->op_rhs, user->phys->solution_time_label));
1052b730f8bSJeremy L Thompson   PetscCall(TSGetTimeStep(ts, &dt));
106e42c1df6SJames Wright   PetscCall(UpdateContextLabel(&user->dt, dt, user->op_rhs, user->phys->timestep_size_label));
10777841947SLeila Ghaffari 
10877841947SLeila Ghaffari   // Global-to-local
1092b730f8bSJeremy L Thompson   PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
11077841947SLeila Ghaffari 
11177841947SLeila Ghaffari   // Place PETSc vectors in CEED vectors
112c798d55aSJames Wright   PetscCall(VecReadP2C(Q_loc, &q_mem_type, user->q_ceed));
113c798d55aSJames Wright   PetscCall(VecP2C(G_loc, &g_mem_type, user->g_ceed));
11477841947SLeila Ghaffari 
11577841947SLeila Ghaffari   // Apply CEED operator
1162b730f8bSJeremy L Thompson   CeedOperatorApply(user->op_rhs, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE);
11777841947SLeila Ghaffari 
11877841947SLeila Ghaffari   // Restore vectors
119c798d55aSJames Wright   PetscCall(VecReadC2P(user->q_ceed, q_mem_type, Q_loc));
120c798d55aSJames Wright   PetscCall(VecC2P(user->g_ceed, g_mem_type, G_loc));
12177841947SLeila Ghaffari 
12277841947SLeila Ghaffari   // Local-to-Global
1232b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(G));
1242b730f8bSJeremy L Thompson   PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G));
12577841947SLeila Ghaffari 
12677841947SLeila Ghaffari   // Inverse of the lumped mass matrix (M is Minv)
1272b730f8bSJeremy L Thompson   PetscCall(VecPointwiseMult(G, G, user->M));
12877841947SLeila Ghaffari 
12977841947SLeila Ghaffari   // Restore vectors
1302b730f8bSJeremy L Thompson   PetscCall(DMRestoreLocalVector(user->dm, &G_loc));
13177841947SLeila Ghaffari 
13277841947SLeila Ghaffari   PetscFunctionReturn(0);
13377841947SLeila Ghaffari }
13477841947SLeila Ghaffari 
135ca69d878SAdeleke O. Bankole // Surface forces function setup
136ca69d878SAdeleke O. Bankole static PetscErrorCode Surface_Forces_NS(DM dm, Vec G_loc, PetscInt num_walls, const PetscInt walls[], PetscScalar *reaction_force) {
137ca69d878SAdeleke O. Bankole   DMLabel            face_label;
138ca69d878SAdeleke O. Bankole   const PetscScalar *g;
139d6734f85SAdeleke O. Bankole   PetscInt           dof, dim = 3;
140ca69d878SAdeleke O. Bankole   MPI_Comm           comm;
141d6734f85SAdeleke O. Bankole   PetscSection       s;
142ca69d878SAdeleke O. Bankole 
143ca69d878SAdeleke O. Bankole   PetscFunctionBeginUser;
144ca69d878SAdeleke O. Bankole   PetscCall(PetscArrayzero(reaction_force, num_walls * dim));
145ca69d878SAdeleke O. Bankole   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
146ca69d878SAdeleke O. Bankole   PetscCall(DMGetLabel(dm, "Face Sets", &face_label));
147ca69d878SAdeleke O. Bankole   PetscCall(VecGetArrayRead(G_loc, &g));
148ca69d878SAdeleke O. Bankole   for (PetscInt w = 0; w < num_walls; w++) {
149ca69d878SAdeleke O. Bankole     const PetscInt wall = walls[w];
150ca69d878SAdeleke O. Bankole     IS             wall_is;
151d6734f85SAdeleke O. Bankole     PetscCall(DMGetLocalSection(dm, &s));
152ca69d878SAdeleke O. Bankole     PetscCall(DMLabelGetStratumIS(face_label, wall, &wall_is));
153ca69d878SAdeleke O. Bankole     if (wall_is) {  // There exist such points on this process
154ca69d878SAdeleke O. Bankole       PetscInt        num_points;
155d6734f85SAdeleke O. Bankole       PetscInt        num_comp = 0;
156ca69d878SAdeleke O. Bankole       const PetscInt *points;
157d6734f85SAdeleke O. Bankole       PetscCall(PetscSectionGetFieldComponents(s, 0, &num_comp));
158ca69d878SAdeleke O. Bankole       PetscCall(ISGetSize(wall_is, &num_points));
159ca69d878SAdeleke O. Bankole       PetscCall(ISGetIndices(wall_is, &points));
160ca69d878SAdeleke O. Bankole       for (PetscInt i = 0; i < num_points; i++) {
161ca69d878SAdeleke O. Bankole         const PetscInt           p = points[i];
162ca69d878SAdeleke O. Bankole         const StateConservative *r;
163ca69d878SAdeleke O. Bankole         PetscCall(DMPlexPointLocalRead(dm, p, g, &r));
164d6734f85SAdeleke O. Bankole         PetscCall(PetscSectionGetDof(s, p, &dof));
165d6734f85SAdeleke O. Bankole         for (PetscInt node = 0; node < dof / num_comp; node++) {
166ca69d878SAdeleke O. Bankole           for (PetscInt j = 0; j < 3; j++) {
167d6734f85SAdeleke O. Bankole             reaction_force[w * dim + j] -= r[node].momentum[j];
168d6734f85SAdeleke O. Bankole           }
169ca69d878SAdeleke O. Bankole         }
170ca69d878SAdeleke O. Bankole       }
171ca69d878SAdeleke O. Bankole       PetscCall(ISRestoreIndices(wall_is, &points));
172ca69d878SAdeleke O. Bankole     }
173ca69d878SAdeleke O. Bankole     PetscCall(ISDestroy(&wall_is));
174ca69d878SAdeleke O. Bankole   }
175ca69d878SAdeleke O. Bankole   PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, reaction_force, dim * num_walls, MPIU_SCALAR, MPI_SUM, comm));
176ca69d878SAdeleke O. Bankole   //  Restore Vectors
177ca69d878SAdeleke O. Bankole   PetscCall(VecRestoreArrayRead(G_loc, &g));
178ca69d878SAdeleke O. Bankole 
179ca69d878SAdeleke O. Bankole   PetscFunctionReturn(0);
180ca69d878SAdeleke O. Bankole }
181ca69d878SAdeleke O. Bankole 
18277841947SLeila Ghaffari // Implicit time-stepper function setup
1832b730f8bSJeremy L Thompson PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) {
18477841947SLeila Ghaffari   User         user = *(User *)user_data;
185c798d55aSJames Wright   PetscScalar  dt;
1865e82a6e1SJeremy L Thompson   Vec          Q_loc = user->Q_loc, Q_dot_loc = user->Q_dot_loc, G_loc;
18777841947SLeila Ghaffari   PetscMemType q_mem_type, q_dot_mem_type, g_mem_type;
18877841947SLeila Ghaffari   PetscFunctionBeginUser;
18977841947SLeila Ghaffari 
1905e82a6e1SJeremy L Thompson   // Get local vectors
191ca69d878SAdeleke O. Bankole   PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
1925e82a6e1SJeremy L Thompson 
1935e82a6e1SJeremy L Thompson   // Update time dependent data
194a0b9a424SJames Wright   PetscCall(UpdateBoundaryValues(user, Q_loc, t));
195e42c1df6SJames Wright   PetscCall(UpdateContextLabel(&user->time, t, user->op_ifunction, user->phys->solution_time_label));
1962b730f8bSJeremy L Thompson   PetscCall(TSGetTimeStep(ts, &dt));
197e42c1df6SJames Wright   PetscCall(UpdateContextLabel(&user->dt, dt, user->op_ifunction, user->phys->timestep_size_label));
19877841947SLeila Ghaffari 
19977841947SLeila Ghaffari   // Global-to-local
200878eb0e7SJames Wright   PetscCall(DMGlobalToLocalBegin(user->dm, Q, INSERT_VALUES, Q_loc));
201878eb0e7SJames Wright   PetscCall(DMGlobalToLocalBegin(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
202878eb0e7SJames Wright   PetscCall(DMGlobalToLocalEnd(user->dm, Q, INSERT_VALUES, Q_loc));
203878eb0e7SJames Wright   PetscCall(DMGlobalToLocalEnd(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
20477841947SLeila Ghaffari 
20577841947SLeila Ghaffari   // Place PETSc vectors in CEED vectors
206c798d55aSJames Wright   PetscCall(VecReadP2C(Q_loc, &q_mem_type, user->q_ceed));
207c798d55aSJames Wright   PetscCall(VecReadP2C(Q_dot_loc, &q_dot_mem_type, user->q_dot_ceed));
208c798d55aSJames Wright   PetscCall(VecP2C(G_loc, &g_mem_type, user->g_ceed));
20977841947SLeila Ghaffari 
21077841947SLeila Ghaffari   // Apply CEED operator
2112b730f8bSJeremy L Thompson   CeedOperatorApply(user->op_ifunction, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE);
21277841947SLeila Ghaffari 
21377841947SLeila Ghaffari   // Restore vectors
214c798d55aSJames Wright   PetscCall(VecReadC2P(user->q_ceed, q_mem_type, Q_loc));
215c798d55aSJames Wright   PetscCall(VecReadC2P(user->q_dot_ceed, q_dot_mem_type, Q_dot_loc));
216c798d55aSJames Wright   PetscCall(VecC2P(user->g_ceed, g_mem_type, G_loc));
21777841947SLeila Ghaffari 
21877841947SLeila Ghaffari   // Local-to-Global
2192b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(G));
2202b730f8bSJeremy L Thompson   PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G));
22177841947SLeila Ghaffari 
22277841947SLeila Ghaffari   // Restore vectors
223ca69d878SAdeleke O. Bankole   PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
22477841947SLeila Ghaffari 
22577841947SLeila Ghaffari   PetscFunctionReturn(0);
22677841947SLeila Ghaffari }
22777841947SLeila Ghaffari 
228e334ad8fSJed Brown static PetscErrorCode MatMult_NS_IJacobian(Mat J, Vec Q, Vec G) {
229e334ad8fSJed Brown   User               user;
230e334ad8fSJed Brown   const PetscScalar *q;
231e334ad8fSJed Brown   PetscScalar       *g;
232e334ad8fSJed Brown   PetscMemType       q_mem_type, g_mem_type;
233e334ad8fSJed Brown   PetscFunctionBeginUser;
2342b730f8bSJeremy L Thompson   PetscCall(MatShellGetContext(J, &user));
2355e82a6e1SJeremy L Thompson   Vec Q_loc = user->Q_dot_loc,  // Note - Q_dot_loc has zero BCs
2365e82a6e1SJeremy L Thompson       G_loc;
2375e82a6e1SJeremy L Thompson 
238e334ad8fSJed Brown   // Get local vectors
2392b730f8bSJeremy L Thompson   PetscCall(DMGetLocalVector(user->dm, &G_loc));
240e334ad8fSJed Brown 
241e334ad8fSJed Brown   // Global-to-local
2422b730f8bSJeremy L Thompson   PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
243e334ad8fSJed Brown 
244e334ad8fSJed Brown   // Place PETSc vectors in CEED vectors
2452b730f8bSJeremy L Thompson   PetscCall(VecGetArrayReadAndMemType(Q_loc, &q, &q_mem_type));
2462b730f8bSJeremy L Thompson   PetscCall(VecGetArrayAndMemType(G_loc, &g, &g_mem_type));
2472b730f8bSJeremy L Thompson   CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER, (PetscScalar *)q);
248e334ad8fSJed Brown   CeedVectorSetArray(user->g_ceed, MemTypeP2C(g_mem_type), CEED_USE_POINTER, g);
249e334ad8fSJed Brown 
250e334ad8fSJed Brown   // Apply CEED operator
2512b730f8bSJeremy L Thompson   CeedOperatorApply(user->op_ijacobian, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE);
252e334ad8fSJed Brown 
253e334ad8fSJed Brown   // Restore vectors
254e334ad8fSJed Brown   CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL);
255e334ad8fSJed Brown   CeedVectorTakeArray(user->g_ceed, MemTypeP2C(g_mem_type), NULL);
2562b730f8bSJeremy L Thompson   PetscCall(VecRestoreArrayReadAndMemType(Q_loc, &q));
2572b730f8bSJeremy L Thompson   PetscCall(VecRestoreArrayAndMemType(G_loc, &g));
258e334ad8fSJed Brown 
259e334ad8fSJed Brown   // Local-to-Global
2602b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(G));
2612b730f8bSJeremy L Thompson   PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G));
262e334ad8fSJed Brown 
263e334ad8fSJed Brown   // Restore vectors
2642b730f8bSJeremy L Thompson   PetscCall(DMRestoreLocalVector(user->dm, &G_loc));
265e334ad8fSJed Brown   PetscFunctionReturn(0);
266e334ad8fSJed Brown }
267e334ad8fSJed Brown 
268e334ad8fSJed Brown PetscErrorCode MatGetDiagonal_NS_IJacobian(Mat A, Vec D) {
269e334ad8fSJed Brown   User         user;
270e334ad8fSJed Brown   Vec          D_loc;
271e334ad8fSJed Brown   PetscScalar *d;
272e334ad8fSJed Brown   PetscMemType mem_type;
273e334ad8fSJed Brown 
274e334ad8fSJed Brown   PetscFunctionBeginUser;
275e334ad8fSJed Brown   MatShellGetContext(A, &user);
276e334ad8fSJed Brown   PetscCall(DMGetLocalVector(user->dm, &D_loc));
277e334ad8fSJed Brown   PetscCall(VecGetArrayAndMemType(D_loc, &d, &mem_type));
278e334ad8fSJed Brown   CeedVectorSetArray(user->g_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, d);
2792b730f8bSJeremy L Thompson   CeedOperatorLinearAssembleDiagonal(user->op_ijacobian, user->g_ceed, CEED_REQUEST_IMMEDIATE);
280e334ad8fSJed Brown   CeedVectorTakeArray(user->g_ceed, MemTypeP2C(mem_type), NULL);
281e334ad8fSJed Brown   PetscCall(VecRestoreArrayAndMemType(D_loc, &d));
282e334ad8fSJed Brown   PetscCall(VecZeroEntries(D));
283e334ad8fSJed Brown   PetscCall(DMLocalToGlobal(user->dm, D_loc, ADD_VALUES, D));
284e334ad8fSJed Brown   PetscCall(DMRestoreLocalVector(user->dm, &D_loc));
285e334ad8fSJed Brown   VecViewFromOptions(D, NULL, "-diag_vec_view");
286e334ad8fSJed Brown   PetscFunctionReturn(0);
287e334ad8fSJed Brown }
288e334ad8fSJed Brown 
2892b730f8bSJeremy L Thompson static PetscErrorCode FormPreallocation(User user, PetscBool pbdiagonal, Mat J, CeedVector *coo_values) {
290544be873SJed Brown   PetscCount ncoo;
291544be873SJed Brown   PetscInt  *rows, *cols;
292544be873SJed Brown 
293544be873SJed Brown   PetscFunctionBeginUser;
294544be873SJed Brown   if (pbdiagonal) {
295544be873SJed Brown     CeedSize l_size;
296544be873SJed Brown     CeedOperatorGetActiveVectorLengths(user->op_ijacobian, &l_size, NULL);
297544be873SJed Brown     ncoo = l_size * 5;
298544be873SJed Brown     rows = malloc(ncoo * sizeof(rows[0]));
299544be873SJed Brown     cols = malloc(ncoo * sizeof(cols[0]));
300544be873SJed Brown     for (PetscCount n = 0; n < l_size / 5; n++) {
301544be873SJed Brown       for (PetscInt i = 0; i < 5; i++) {
302544be873SJed Brown         for (PetscInt j = 0; j < 5; j++) {
303544be873SJed Brown           rows[(n * 5 + i) * 5 + j] = n * 5 + i;
304544be873SJed Brown           cols[(n * 5 + i) * 5 + j] = n * 5 + j;
305544be873SJed Brown         }
306544be873SJed Brown       }
307544be873SJed Brown     }
308544be873SJed Brown   } else {
3092b730f8bSJeremy L Thompson     PetscCall(CeedOperatorLinearAssembleSymbolic(user->op_ijacobian, &ncoo, &rows, &cols));
310544be873SJed Brown   }
311544be873SJed Brown   PetscCall(MatSetPreallocationCOOLocal(J, ncoo, rows, cols));
312544be873SJed Brown   free(rows);
313544be873SJed Brown   free(cols);
314544be873SJed Brown   CeedVectorCreate(user->ceed, ncoo, coo_values);
315544be873SJed Brown   PetscFunctionReturn(0);
316544be873SJed Brown }
317544be873SJed Brown 
3182b730f8bSJeremy L Thompson static PetscErrorCode FormSetValues(User user, PetscBool pbdiagonal, Mat J, CeedVector coo_values) {
319544be873SJed Brown   CeedMemType        mem_type = CEED_MEM_HOST;
320544be873SJed Brown   const PetscScalar *values;
321544be873SJed Brown   MatType            mat_type;
322544be873SJed Brown 
323544be873SJed Brown   PetscFunctionBeginUser;
324544be873SJed Brown   PetscCall(MatGetType(J, &mat_type));
3252b730f8bSJeremy L Thompson   if (strstr(mat_type, "kokkos") || strstr(mat_type, "cusparse")) mem_type = CEED_MEM_DEVICE;
326544be873SJed Brown   if (user->app_ctx->pmat_pbdiagonal) {
3272b730f8bSJeremy L Thompson     CeedOperatorLinearAssemblePointBlockDiagonal(user->op_ijacobian, coo_values, CEED_REQUEST_IMMEDIATE);
328544be873SJed Brown   } else {
329544be873SJed Brown     CeedOperatorLinearAssemble(user->op_ijacobian, coo_values);
330544be873SJed Brown   }
331544be873SJed Brown   CeedVectorGetArrayRead(coo_values, mem_type, &values);
332544be873SJed Brown   PetscCall(MatSetValuesCOO(J, values, INSERT_VALUES));
333544be873SJed Brown   CeedVectorRestoreArrayRead(coo_values, &values);
334544be873SJed Brown   PetscFunctionReturn(0);
335544be873SJed Brown }
336544be873SJed Brown 
3372b730f8bSJeremy L Thompson PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, PetscReal shift, Mat J, Mat J_pre, void *user_data) {
338e334ad8fSJed Brown   User      user = *(User *)user_data;
339d6bf345cSJed Brown   PetscBool J_is_shell, J_is_mffd, J_pre_is_shell;
340e334ad8fSJed Brown   PetscFunctionBeginUser;
341*17b0d5c6SJeremy L Thompson   if (user->phys->ijacobian_time_shift_label) CeedOperatorSetContextDouble(user->op_ijacobian, user->phys->ijacobian_time_shift_label, &shift);
342d6bf345cSJed Brown   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATMFFD, &J_is_mffd));
343e334ad8fSJed Brown   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATSHELL, &J_is_shell));
3442b730f8bSJeremy L Thompson   PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATSHELL, &J_pre_is_shell));
345e334ad8fSJed Brown   if (!user->matrices_set_up) {
346e334ad8fSJed Brown     if (J_is_shell) {
347e334ad8fSJed Brown       PetscCall(MatShellSetContext(J, user));
3482b730f8bSJeremy L Thompson       PetscCall(MatShellSetOperation(J, MATOP_MULT, (void (*)(void))MatMult_NS_IJacobian));
3492b730f8bSJeremy L Thompson       PetscCall(MatShellSetOperation(J, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_NS_IJacobian));
350e334ad8fSJed Brown       PetscCall(MatSetUp(J));
351e334ad8fSJed Brown     }
352e334ad8fSJed Brown     if (!J_pre_is_shell) {
3532b730f8bSJeremy L Thompson       PetscCall(FormPreallocation(user, user->app_ctx->pmat_pbdiagonal, J_pre, &user->coo_values_pmat));
354544be873SJed Brown     }
355d6bf345cSJed Brown     if (J != J_pre && !J_is_shell && !J_is_mffd) {
356544be873SJed Brown       PetscCall(FormPreallocation(user, PETSC_FALSE, J, &user->coo_values_amat));
357544be873SJed Brown     }
358e334ad8fSJed Brown     user->matrices_set_up = true;
359e334ad8fSJed Brown   }
360e334ad8fSJed Brown   if (!J_pre_is_shell) {
3612b730f8bSJeremy L Thompson     PetscCall(FormSetValues(user, user->app_ctx->pmat_pbdiagonal, J_pre, user->coo_values_pmat));
362e334ad8fSJed Brown   }
363d6bf345cSJed Brown   if (user->coo_values_amat) {
364d6bf345cSJed Brown     PetscCall(FormSetValues(user, PETSC_FALSE, J, user->coo_values_amat));
365d6bf345cSJed Brown   } else if (J_is_mffd) {
366d6bf345cSJed Brown     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
367d6bf345cSJed Brown     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
368d6bf345cSJed Brown   }
369e334ad8fSJed Brown   PetscFunctionReturn(0);
370e334ad8fSJed Brown }
371e334ad8fSJed Brown 
3722b730f8bSJeremy L Thompson PetscErrorCode WriteOutput(User user, Vec Q, PetscInt step_no, PetscScalar time) {
37377841947SLeila Ghaffari   Vec         Q_loc;
37477841947SLeila Ghaffari   char        file_path[PETSC_MAX_PATH_LEN];
37577841947SLeila Ghaffari   PetscViewer viewer;
37677841947SLeila Ghaffari   PetscFunctionBeginUser;
37777841947SLeila Ghaffari 
37837cbb16aSJed Brown   if (user->app_ctx->checkpoint_vtk) {
37977841947SLeila Ghaffari     // Set up output
380f14660a4SJames Wright     PetscCall(DMGetLocalVector(user->dm, &Q_loc));
381f14660a4SJames Wright     PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec"));
382f14660a4SJames Wright     PetscCall(VecZeroEntries(Q_loc));
383f14660a4SJames Wright     PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
38477841947SLeila Ghaffari 
38577841947SLeila Ghaffari     // Output
38637cbb16aSJed Brown     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
387f14660a4SJames Wright 
3882b730f8bSJeremy L Thompson     PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path, FILE_MODE_WRITE, &viewer));
389f14660a4SJames Wright     PetscCall(VecView(Q_loc, viewer));
390f14660a4SJames Wright     PetscCall(PetscViewerDestroy(&viewer));
39177841947SLeila Ghaffari     if (user->dm_viz) {
39277841947SLeila Ghaffari       Vec         Q_refined, Q_refined_loc;
39377841947SLeila Ghaffari       char        file_path_refined[PETSC_MAX_PATH_LEN];
39477841947SLeila Ghaffari       PetscViewer viewer_refined;
39577841947SLeila Ghaffari 
396f14660a4SJames Wright       PetscCall(DMGetGlobalVector(user->dm_viz, &Q_refined));
397f14660a4SJames Wright       PetscCall(DMGetLocalVector(user->dm_viz, &Q_refined_loc));
398f14660a4SJames Wright       PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined"));
399f14660a4SJames Wright 
400f14660a4SJames Wright       PetscCall(MatInterpolate(user->interp_viz, Q, Q_refined));
401f14660a4SJames Wright       PetscCall(VecZeroEntries(Q_refined_loc));
4022b730f8bSJeremy L Thompson       PetscCall(DMGlobalToLocal(user->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc));
403f14660a4SJames Wright 
40437cbb16aSJed Brown       PetscCall(
40537cbb16aSJed Brown           PetscSNPrintf(file_path_refined, sizeof file_path_refined, "%s/nsrefined-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
406f14660a4SJames Wright 
4072b730f8bSJeremy L Thompson       PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined), file_path_refined, FILE_MODE_WRITE, &viewer_refined));
408f14660a4SJames Wright       PetscCall(VecView(Q_refined_loc, viewer_refined));
409f14660a4SJames Wright       PetscCall(DMRestoreLocalVector(user->dm_viz, &Q_refined_loc));
410f14660a4SJames Wright       PetscCall(DMRestoreGlobalVector(user->dm_viz, &Q_refined));
411f14660a4SJames Wright       PetscCall(PetscViewerDestroy(&viewer_refined));
41277841947SLeila Ghaffari     }
413f14660a4SJames Wright     PetscCall(DMRestoreLocalVector(user->dm, &Q_loc));
41437cbb16aSJed Brown   }
41577841947SLeila Ghaffari 
41677841947SLeila Ghaffari   // Save data in a binary file for continuation of simulations
41769293791SJames Wright   if (user->app_ctx->add_stepnum2bin) {
41837cbb16aSJed Brown     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", user->app_ctx->output_dir, step_no));
41969293791SJames Wright   } else {
4202b730f8bSJeremy L Thompson     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", user->app_ctx->output_dir));
42169293791SJames Wright   }
4222b730f8bSJeremy L Thompson   PetscCall(PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, &viewer));
423f14660a4SJames Wright 
4244de8550aSJed Brown   PetscInt token = FLUIDS_FILE_TOKEN;
4254de8550aSJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &token, 1, PETSC_INT));
4264de8550aSJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &step_no, 1, PETSC_INT));
4274de8550aSJed Brown   time /= user->units->second;  // Dimensionalize time back
4284de8550aSJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL));
429f14660a4SJames Wright   PetscCall(VecView(Q, viewer));
430f14660a4SJames Wright   PetscCall(PetscViewerDestroy(&viewer));
431f14660a4SJames Wright   PetscFunctionReturn(0);
432f14660a4SJames Wright }
433f14660a4SJames Wright 
434ca69d878SAdeleke O. Bankole // CSV Monitor
435ca69d878SAdeleke O. Bankole PetscErrorCode TSMonitor_WallForce(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
436ca69d878SAdeleke O. Bankole   User              user = ctx;
437ca69d878SAdeleke O. Bankole   Vec               G_loc;
438ca69d878SAdeleke O. Bankole   PetscInt          num_wall = user->app_ctx->wall_forces.num_wall, dim = 3;
439ca69d878SAdeleke O. Bankole   const PetscInt   *walls  = user->app_ctx->wall_forces.walls;
440ca69d878SAdeleke O. Bankole   PetscViewer       viewer = user->app_ctx->wall_forces.viewer;
441ca69d878SAdeleke O. Bankole   PetscViewerFormat format = user->app_ctx->wall_forces.viewer_format;
442ca69d878SAdeleke O. Bankole   PetscScalar      *reaction_force;
443ca69d878SAdeleke O. Bankole   PetscBool         iascii;
444ca69d878SAdeleke O. Bankole 
445ca69d878SAdeleke O. Bankole   PetscFunctionBeginUser;
446ca69d878SAdeleke O. Bankole   if (!viewer) PetscFunctionReturn(0);
447ca69d878SAdeleke O. Bankole   PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
448ca69d878SAdeleke O. Bankole   PetscCall(PetscMalloc1(num_wall * dim, &reaction_force));
449ca69d878SAdeleke O. Bankole   PetscCall(Surface_Forces_NS(user->dm, G_loc, num_wall, walls, reaction_force));
450ca69d878SAdeleke O. Bankole   PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
451ca69d878SAdeleke O. Bankole 
452ca69d878SAdeleke O. Bankole   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
453ca69d878SAdeleke O. Bankole 
454ca69d878SAdeleke O. Bankole   if (iascii) {
455ca69d878SAdeleke O. Bankole     if (format == PETSC_VIEWER_ASCII_CSV && !user->app_ctx->wall_forces.header_written) {
456ca69d878SAdeleke O. Bankole       PetscCall(PetscViewerASCIIPrintf(viewer, "Step,Time,Wall,ForceX,ForceY,ForceZ\n"));
457ca69d878SAdeleke O. Bankole       user->app_ctx->wall_forces.header_written = PETSC_TRUE;
458ca69d878SAdeleke O. Bankole     }
459ca69d878SAdeleke O. Bankole     for (PetscInt w = 0; w < num_wall; w++) {
460ca69d878SAdeleke O. Bankole       PetscInt wall = walls[w];
461ca69d878SAdeleke O. Bankole       if (format == PETSC_VIEWER_ASCII_CSV) {
462ca69d878SAdeleke O. Bankole         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT ",%g,%" PetscInt_FMT ",%g,%g,%g\n", step_no, time, wall,
463ca69d878SAdeleke O. Bankole                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
464ca69d878SAdeleke O. Bankole 
465ca69d878SAdeleke O. Bankole       } else {
466ca69d878SAdeleke O. Bankole         PetscCall(PetscViewerASCIIPrintf(viewer, "Wall %" PetscInt_FMT " Forces: Force_x = %12g, Force_y = %12g, Force_z = %12g\n", wall,
467ca69d878SAdeleke O. Bankole                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
468ca69d878SAdeleke O. Bankole       }
469ca69d878SAdeleke O. Bankole     }
470ca69d878SAdeleke O. Bankole   }
471ca69d878SAdeleke O. Bankole   PetscCall(PetscFree(reaction_force));
472ca69d878SAdeleke O. Bankole   PetscFunctionReturn(0);
473ca69d878SAdeleke O. Bankole }
474ca69d878SAdeleke O. Bankole 
475f14660a4SJames Wright // User provided TS Monitor
4762b730f8bSJeremy L Thompson PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
477f14660a4SJames Wright   User user = ctx;
478f14660a4SJames Wright   PetscFunctionBeginUser;
479f14660a4SJames Wright 
48037cbb16aSJed Brown   // Print every 'checkpoint_interval' steps
481894de27cSJames Wright   if (user->app_ctx->checkpoint_interval <= 0 || step_no % user->app_ctx->checkpoint_interval != 0 ||
482894de27cSJames Wright       (user->app_ctx->cont_steps == step_no && step_no != 0))
483894de27cSJames Wright     PetscFunctionReturn(0);
484f14660a4SJames Wright 
485f14660a4SJames Wright   PetscCall(WriteOutput(user, Q, step_no, time));
48677841947SLeila Ghaffari 
48777841947SLeila Ghaffari   PetscFunctionReturn(0);
48877841947SLeila Ghaffari }
48977841947SLeila Ghaffari 
49077841947SLeila Ghaffari // TS: Create, setup, and solve
4912b730f8bSJeremy L Thompson PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, Vec *Q, PetscScalar *f_time, TS *ts) {
49277841947SLeila Ghaffari   MPI_Comm    comm = user->comm;
49377841947SLeila Ghaffari   TSAdapt     adapt;
49477841947SLeila Ghaffari   PetscScalar final_time;
49577841947SLeila Ghaffari   PetscFunctionBeginUser;
49677841947SLeila Ghaffari 
4972b730f8bSJeremy L Thompson   PetscCall(TSCreate(comm, ts));
4982b730f8bSJeremy L Thompson   PetscCall(TSSetDM(*ts, dm));
49977841947SLeila Ghaffari   if (phys->implicit) {
5002b730f8bSJeremy L Thompson     PetscCall(TSSetType(*ts, TSBDF));
50177841947SLeila Ghaffari     if (user->op_ifunction) {
5022b730f8bSJeremy L Thompson       PetscCall(TSSetIFunction(*ts, NULL, IFunction_NS, &user));
50377841947SLeila Ghaffari     } else {  // Implicit integrators can fall back to using an RHSFunction
5042b730f8bSJeremy L Thompson       PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
50577841947SLeila Ghaffari     }
506e334ad8fSJed Brown     if (user->op_ijacobian) {
5072b730f8bSJeremy L Thompson       PetscCall(DMTSSetIJacobian(dm, FormIJacobian_NS, &user));
508544be873SJed Brown       if (app_ctx->amat_type) {
509544be873SJed Brown         Mat Pmat, Amat;
5102b730f8bSJeremy L Thompson         PetscCall(DMCreateMatrix(dm, &Pmat));
5112b730f8bSJeremy L Thompson         PetscCall(DMSetMatType(dm, app_ctx->amat_type));
5122b730f8bSJeremy L Thompson         PetscCall(DMCreateMatrix(dm, &Amat));
5132b730f8bSJeremy L Thompson         PetscCall(TSSetIJacobian(*ts, Amat, Pmat, NULL, NULL));
5142b730f8bSJeremy L Thompson         PetscCall(MatDestroy(&Amat));
5152b730f8bSJeremy L Thompson         PetscCall(MatDestroy(&Pmat));
516544be873SJed Brown       }
517e334ad8fSJed Brown     }
51877841947SLeila Ghaffari   } else {
5192b730f8bSJeremy L Thompson     if (!user->op_rhs) SETERRQ(comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction");
5202b730f8bSJeremy L Thompson     PetscCall(TSSetType(*ts, TSRK));
5212b730f8bSJeremy L Thompson     PetscCall(TSRKSetType(*ts, TSRK5F));
5222b730f8bSJeremy L Thompson     PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
52377841947SLeila Ghaffari   }
5242b730f8bSJeremy L Thompson   PetscCall(TSSetMaxTime(*ts, 500. * user->units->second));
5252b730f8bSJeremy L Thompson   PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER));
526cf7a0454SJed Brown   PetscCall(TSSetErrorIfStepFails(*ts, PETSC_FALSE));
5272b730f8bSJeremy L Thompson   PetscCall(TSSetTimeStep(*ts, 1.e-2 * user->units->second));
5288fb33541SJames Wright   if (app_ctx->test_type != TESTTYPE_NONE) {
5292b730f8bSJeremy L Thompson     PetscCall(TSSetMaxSteps(*ts, 10));
5302b730f8bSJeremy L Thompson   }
5312b730f8bSJeremy L Thompson   PetscCall(TSGetAdapt(*ts, &adapt));
5322b730f8bSJeremy L Thompson   PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * user->units->second, 1.e2 * user->units->second));
5332b730f8bSJeremy L Thompson   PetscCall(TSSetFromOptions(*ts));
534a0b9a424SJames Wright   user->time = user->time_bc_set = -1.0;  // require all BCs and ctx to be updated
5355e82a6e1SJeremy L Thompson   user->dt                       = -1.0;
53677841947SLeila Ghaffari   if (!app_ctx->cont_steps) {  // print initial condition
5378fb33541SJames Wright     if (app_ctx->test_type == TESTTYPE_NONE) {
5382b730f8bSJeremy L Thompson       PetscCall(TSMonitor_NS(*ts, 0, 0., *Q, user));
53977841947SLeila Ghaffari     }
54077841947SLeila Ghaffari   } else {  // continue from time of last output
54177841947SLeila Ghaffari     PetscInt    count;
54277841947SLeila Ghaffari     PetscViewer viewer;
5432b730f8bSJeremy L Thompson 
5444de8550aSJed Brown     if (app_ctx->cont_time <= 0) {  // Legacy files did not include step number and time
5452b730f8bSJeremy L Thompson       PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_time_file, FILE_MODE_READ, &viewer));
5464de8550aSJed Brown       PetscCall(PetscViewerBinaryRead(viewer, &app_ctx->cont_time, 1, &count, PETSC_REAL));
5472b730f8bSJeremy L Thompson       PetscCall(PetscViewerDestroy(&viewer));
5484de8550aSJed Brown       PetscCheck(app_ctx->cont_steps != -1, comm, PETSC_ERR_ARG_INCOMP,
5494de8550aSJed Brown                  "-continue step number not specified, but checkpoint file does not contain a step number (likely written by older code version)");
5504de8550aSJed Brown     }
5514de8550aSJed Brown     PetscCall(TSSetTime(*ts, app_ctx->cont_time * user->units->second));
552d8e0aecdSJed Brown     PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps));
55377841947SLeila Ghaffari   }
5548fb33541SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) {
5552b730f8bSJeremy L Thompson     PetscCall(TSMonitorSet(*ts, TSMonitor_NS, user, NULL));
5568fb33541SJames Wright   }
557ca69d878SAdeleke O. Bankole   if (app_ctx->wall_forces.viewer) {
558ca69d878SAdeleke O. Bankole     PetscCall(TSMonitorSet(*ts, TSMonitor_WallForce, user, NULL));
559ca69d878SAdeleke O. Bankole   }
560b7d66439SJames Wright   if (app_ctx->turb_spanstats_enable) {
561a175e481SJames Wright     PetscCall(TSMonitorSet(*ts, TSMonitor_Statistics, user, NULL));
562495b9769SJames Wright     CeedScalar previous_time = app_ctx->cont_time * user->units->second;
563*17b0d5c6SJeremy L Thompson     CeedOperatorSetContextDouble(user->spanstats.op_stats_collect, user->spanstats.previous_time_label, &previous_time);
564a175e481SJames Wright   }
56577841947SLeila Ghaffari 
56677841947SLeila Ghaffari   // Solve
567d8e0aecdSJed Brown   PetscReal start_time;
568d8e0aecdSJed Brown   PetscInt  start_step;
5692b730f8bSJeremy L Thompson   PetscCall(TSGetTime(*ts, &start_time));
570d8e0aecdSJed Brown   PetscCall(TSGetStepNumber(*ts, &start_step));
5717b39487dSJeremy L Thompson 
5727b39487dSJeremy L Thompson   PetscPreLoadBegin(PETSC_FALSE, "Fluids Solve");
5737b39487dSJeremy L Thompson   PetscCall(TSSetTime(*ts, start_time));
574d8e0aecdSJed Brown   PetscCall(TSSetStepNumber(*ts, start_step));
5757b39487dSJeremy L Thompson   if (PetscPreLoadingOn) {
5767b39487dSJeremy L Thompson     // LCOV_EXCL_START
5777b39487dSJeremy L Thompson     SNES      snes;
5787b39487dSJeremy L Thompson     Vec       Q_preload;
5797b39487dSJeremy L Thompson     PetscReal rtol;
5807b39487dSJeremy L Thompson     PetscCall(VecDuplicate(*Q, &Q_preload));
5817b39487dSJeremy L Thompson     PetscCall(VecCopy(*Q, Q_preload));
5827b39487dSJeremy L Thompson     PetscCall(TSGetSNES(*ts, &snes));
5837b39487dSJeremy L Thompson     PetscCall(SNESGetTolerances(snes, NULL, &rtol, NULL, NULL, NULL));
5842b730f8bSJeremy L Thompson     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, .99, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
5857b39487dSJeremy L Thompson     PetscCall(TSSetSolution(*ts, *Q));
5867b39487dSJeremy L Thompson     PetscCall(TSStep(*ts));
5872b730f8bSJeremy L Thompson     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, rtol, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
5887b39487dSJeremy L Thompson     PetscCall(VecDestroy(&Q_preload));
5897b39487dSJeremy L Thompson     // LCOV_EXCL_STOP
5907b39487dSJeremy L Thompson   } else {
5912b730f8bSJeremy L Thompson     PetscCall(PetscBarrier((PetscObject)*ts));
5922b730f8bSJeremy L Thompson     PetscCall(TSSolve(*ts, *Q));
5937b39487dSJeremy L Thompson   }
5947b39487dSJeremy L Thompson   PetscPreLoadEnd();
5957b39487dSJeremy L Thompson 
5967b39487dSJeremy L Thompson   PetscCall(TSGetSolveTime(*ts, &final_time));
59777841947SLeila Ghaffari   *f_time = final_time;
5987b39487dSJeremy L Thompson 
5998fb33541SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) {
600f14660a4SJames Wright     PetscInt step_no;
601f14660a4SJames Wright     PetscCall(TSGetStepNumber(*ts, &step_no));
602a175e481SJames Wright     if (user->app_ctx->checkpoint_interval > 0 || user->app_ctx->checkpoint_interval == -1) {
603f14660a4SJames Wright       PetscCall(WriteOutput(user, *Q, step_no, final_time));
604f14660a4SJames Wright     }
605f14660a4SJames Wright 
6067b39487dSJeremy L Thompson     PetscLogEvent stage_id;
6077b39487dSJeremy L Thompson     PetscStageLog stage_log;
6087b39487dSJeremy L Thompson 
6097b39487dSJeremy L Thompson     PetscCall(PetscLogStageGetId("Fluids Solve", &stage_id));
6107b39487dSJeremy L Thompson     PetscCall(PetscLogGetStageLog(&stage_log));
6112b730f8bSJeremy L Thompson     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time taken for solution (sec): %g\n", stage_log->stageInfo[stage_id].perfInfo.time));
61277841947SLeila Ghaffari   }
61377841947SLeila Ghaffari   PetscFunctionReturn(0);
61477841947SLeila Ghaffari }
615