xref: /libCEED/examples/fluids/src/setupts.c (revision 1563739518ce9a7db5e4395aacf3f48f57a2512c)
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 
77*15637395SJames Wright // @brief Update the context label value to new value if necessary.
78*15637395SJames Wright // @note This only supports labels with scalar label values (ie. not arrays)
79*15637395SJames Wright PetscErrorCode UpdateContextLabel(MPI_Comm comm, PetscScalar update_value, CeedOperator op, CeedContextFieldLabel label) {
80*15637395SJames Wright   PetscScalar label_value;
81e42c1df6SJames Wright 
82*15637395SJames Wright   PetscFunctionBeginUser;
83*15637395SJames Wright   PetscCheck(label, comm, PETSC_ERR_ARG_BADPTR, "Label should be non-NULL");
84*15637395SJames Wright 
85*15637395SJames Wright   {
86*15637395SJames Wright     size_t             num_elements;
87*15637395SJames Wright     const PetscScalar *label_values;
88*15637395SJames Wright     CeedOperatorGetContextDoubleRead(op, label, &num_elements, &label_values);
89*15637395SJames Wright     PetscCheck(num_elements == 1, comm, PETSC_ERR_SUP, "%s does not support labels with more than 1 value. Label has %zu values", __func__,
90*15637395SJames Wright                num_elements);
91*15637395SJames Wright     label_value = *label_values;
92*15637395SJames Wright     CeedOperatorRestoreContextDoubleRead(op, label, &label_values);
93e42c1df6SJames Wright   }
94*15637395SJames Wright 
95*15637395SJames Wright   if (label_value != update_value) {
96*15637395SJames Wright     CeedOperatorSetContextDouble(op, label, &update_value);
97e42c1df6SJames Wright   }
98e42c1df6SJames Wright   PetscFunctionReturn(0);
99e42c1df6SJames Wright }
100e42c1df6SJames Wright 
10177841947SLeila Ghaffari // RHS (Explicit time-stepper) function setup
10277841947SLeila Ghaffari //   This is the RHS of the ODE, given as u_t = G(t,u)
10377841947SLeila Ghaffari //   This function takes in a state vector Q and writes into G
10477841947SLeila Ghaffari PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) {
10577841947SLeila Ghaffari   User         user = *(User *)user_data;
106*15637395SJames Wright   MPI_Comm     comm = PetscObjectComm((PetscObject)ts);
107c798d55aSJames Wright   PetscScalar  dt;
1085e82a6e1SJeremy L Thompson   Vec          Q_loc = user->Q_loc, G_loc;
10977841947SLeila Ghaffari   PetscMemType q_mem_type, g_mem_type;
11077841947SLeila Ghaffari   PetscFunctionBeginUser;
11177841947SLeila Ghaffari 
1125e82a6e1SJeremy L Thompson   // Get local vector
1132b730f8bSJeremy L Thompson   PetscCall(DMGetLocalVector(user->dm, &G_loc));
1145e82a6e1SJeremy L Thompson 
1155e82a6e1SJeremy L Thompson   // Update time dependent data
116a0b9a424SJames Wright   PetscCall(UpdateBoundaryValues(user, Q_loc, t));
117*15637395SJames Wright   if (user->phys->solution_time_label) PetscCall(UpdateContextLabel(comm, t, user->op_rhs, user->phys->solution_time_label));
1182b730f8bSJeremy L Thompson   PetscCall(TSGetTimeStep(ts, &dt));
119*15637395SJames Wright   if (user->phys->timestep_size_label) PetscCall(UpdateContextLabel(comm, dt, user->op_rhs, user->phys->timestep_size_label));
12077841947SLeila Ghaffari 
12177841947SLeila Ghaffari   // Global-to-local
1222b730f8bSJeremy L Thompson   PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
12377841947SLeila Ghaffari 
12477841947SLeila Ghaffari   // Place PETSc vectors in CEED vectors
125c798d55aSJames Wright   PetscCall(VecReadP2C(Q_loc, &q_mem_type, user->q_ceed));
126c798d55aSJames Wright   PetscCall(VecP2C(G_loc, &g_mem_type, user->g_ceed));
12777841947SLeila Ghaffari 
12877841947SLeila Ghaffari   // Apply CEED operator
1292b730f8bSJeremy L Thompson   CeedOperatorApply(user->op_rhs, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE);
13077841947SLeila Ghaffari 
13177841947SLeila Ghaffari   // Restore vectors
132c798d55aSJames Wright   PetscCall(VecReadC2P(user->q_ceed, q_mem_type, Q_loc));
133c798d55aSJames Wright   PetscCall(VecC2P(user->g_ceed, g_mem_type, G_loc));
13477841947SLeila Ghaffari 
13577841947SLeila Ghaffari   // Local-to-Global
1362b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(G));
1372b730f8bSJeremy L Thompson   PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G));
13877841947SLeila Ghaffari 
13977841947SLeila Ghaffari   // Inverse of the lumped mass matrix (M is Minv)
1402b730f8bSJeremy L Thompson   PetscCall(VecPointwiseMult(G, G, user->M));
14177841947SLeila Ghaffari 
14277841947SLeila Ghaffari   // Restore vectors
1432b730f8bSJeremy L Thompson   PetscCall(DMRestoreLocalVector(user->dm, &G_loc));
14477841947SLeila Ghaffari 
14577841947SLeila Ghaffari   PetscFunctionReturn(0);
14677841947SLeila Ghaffari }
14777841947SLeila Ghaffari 
148ca69d878SAdeleke O. Bankole // Surface forces function setup
149ca69d878SAdeleke O. Bankole static PetscErrorCode Surface_Forces_NS(DM dm, Vec G_loc, PetscInt num_walls, const PetscInt walls[], PetscScalar *reaction_force) {
150ca69d878SAdeleke O. Bankole   DMLabel            face_label;
151ca69d878SAdeleke O. Bankole   const PetscScalar *g;
152d6734f85SAdeleke O. Bankole   PetscInt           dof, dim = 3;
153ca69d878SAdeleke O. Bankole   MPI_Comm           comm;
154d6734f85SAdeleke O. Bankole   PetscSection       s;
155ca69d878SAdeleke O. Bankole 
156ca69d878SAdeleke O. Bankole   PetscFunctionBeginUser;
157ca69d878SAdeleke O. Bankole   PetscCall(PetscArrayzero(reaction_force, num_walls * dim));
158ca69d878SAdeleke O. Bankole   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
159ca69d878SAdeleke O. Bankole   PetscCall(DMGetLabel(dm, "Face Sets", &face_label));
160ca69d878SAdeleke O. Bankole   PetscCall(VecGetArrayRead(G_loc, &g));
161ca69d878SAdeleke O. Bankole   for (PetscInt w = 0; w < num_walls; w++) {
162ca69d878SAdeleke O. Bankole     const PetscInt wall = walls[w];
163ca69d878SAdeleke O. Bankole     IS             wall_is;
164d6734f85SAdeleke O. Bankole     PetscCall(DMGetLocalSection(dm, &s));
165ca69d878SAdeleke O. Bankole     PetscCall(DMLabelGetStratumIS(face_label, wall, &wall_is));
166ca69d878SAdeleke O. Bankole     if (wall_is) {  // There exist such points on this process
167ca69d878SAdeleke O. Bankole       PetscInt        num_points;
168d6734f85SAdeleke O. Bankole       PetscInt        num_comp = 0;
169ca69d878SAdeleke O. Bankole       const PetscInt *points;
170d6734f85SAdeleke O. Bankole       PetscCall(PetscSectionGetFieldComponents(s, 0, &num_comp));
171ca69d878SAdeleke O. Bankole       PetscCall(ISGetSize(wall_is, &num_points));
172ca69d878SAdeleke O. Bankole       PetscCall(ISGetIndices(wall_is, &points));
173ca69d878SAdeleke O. Bankole       for (PetscInt i = 0; i < num_points; i++) {
174ca69d878SAdeleke O. Bankole         const PetscInt           p = points[i];
175ca69d878SAdeleke O. Bankole         const StateConservative *r;
176ca69d878SAdeleke O. Bankole         PetscCall(DMPlexPointLocalRead(dm, p, g, &r));
177d6734f85SAdeleke O. Bankole         PetscCall(PetscSectionGetDof(s, p, &dof));
178d6734f85SAdeleke O. Bankole         for (PetscInt node = 0; node < dof / num_comp; node++) {
179ca69d878SAdeleke O. Bankole           for (PetscInt j = 0; j < 3; j++) {
180d6734f85SAdeleke O. Bankole             reaction_force[w * dim + j] -= r[node].momentum[j];
181d6734f85SAdeleke O. Bankole           }
182ca69d878SAdeleke O. Bankole         }
183ca69d878SAdeleke O. Bankole       }
184ca69d878SAdeleke O. Bankole       PetscCall(ISRestoreIndices(wall_is, &points));
185ca69d878SAdeleke O. Bankole     }
186ca69d878SAdeleke O. Bankole     PetscCall(ISDestroy(&wall_is));
187ca69d878SAdeleke O. Bankole   }
188ca69d878SAdeleke O. Bankole   PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, reaction_force, dim * num_walls, MPIU_SCALAR, MPI_SUM, comm));
189ca69d878SAdeleke O. Bankole   //  Restore Vectors
190ca69d878SAdeleke O. Bankole   PetscCall(VecRestoreArrayRead(G_loc, &g));
191ca69d878SAdeleke O. Bankole 
192ca69d878SAdeleke O. Bankole   PetscFunctionReturn(0);
193ca69d878SAdeleke O. Bankole }
194ca69d878SAdeleke O. Bankole 
19577841947SLeila Ghaffari // Implicit time-stepper function setup
1962b730f8bSJeremy L Thompson PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) {
19777841947SLeila Ghaffari   User         user = *(User *)user_data;
198*15637395SJames Wright   MPI_Comm     comm = PetscObjectComm((PetscObject)ts);
199c798d55aSJames Wright   PetscScalar  dt;
2005e82a6e1SJeremy L Thompson   Vec          Q_loc = user->Q_loc, Q_dot_loc = user->Q_dot_loc, G_loc;
20177841947SLeila Ghaffari   PetscMemType q_mem_type, q_dot_mem_type, g_mem_type;
20277841947SLeila Ghaffari   PetscFunctionBeginUser;
20377841947SLeila Ghaffari 
2045e82a6e1SJeremy L Thompson   // Get local vectors
205ca69d878SAdeleke O. Bankole   PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
2065e82a6e1SJeremy L Thompson 
2075e82a6e1SJeremy L Thompson   // Update time dependent data
208a0b9a424SJames Wright   PetscCall(UpdateBoundaryValues(user, Q_loc, t));
209*15637395SJames Wright   if (user->phys->solution_time_label) PetscCall(UpdateContextLabel(comm, t, user->op_ifunction, user->phys->solution_time_label));
2102b730f8bSJeremy L Thompson   PetscCall(TSGetTimeStep(ts, &dt));
211*15637395SJames Wright   if (user->phys->timestep_size_label) PetscCall(UpdateContextLabel(comm, dt, user->op_ifunction, user->phys->timestep_size_label));
21277841947SLeila Ghaffari 
21377841947SLeila Ghaffari   // Global-to-local
214878eb0e7SJames Wright   PetscCall(DMGlobalToLocalBegin(user->dm, Q, INSERT_VALUES, Q_loc));
215878eb0e7SJames Wright   PetscCall(DMGlobalToLocalBegin(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
216878eb0e7SJames Wright   PetscCall(DMGlobalToLocalEnd(user->dm, Q, INSERT_VALUES, Q_loc));
217878eb0e7SJames Wright   PetscCall(DMGlobalToLocalEnd(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
21877841947SLeila Ghaffari 
21977841947SLeila Ghaffari   // Place PETSc vectors in CEED vectors
220c798d55aSJames Wright   PetscCall(VecReadP2C(Q_loc, &q_mem_type, user->q_ceed));
221c798d55aSJames Wright   PetscCall(VecReadP2C(Q_dot_loc, &q_dot_mem_type, user->q_dot_ceed));
222c798d55aSJames Wright   PetscCall(VecP2C(G_loc, &g_mem_type, user->g_ceed));
22377841947SLeila Ghaffari 
22477841947SLeila Ghaffari   // Apply CEED operator
2252b730f8bSJeremy L Thompson   CeedOperatorApply(user->op_ifunction, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE);
22677841947SLeila Ghaffari 
22777841947SLeila Ghaffari   // Restore vectors
228c798d55aSJames Wright   PetscCall(VecReadC2P(user->q_ceed, q_mem_type, Q_loc));
229c798d55aSJames Wright   PetscCall(VecReadC2P(user->q_dot_ceed, q_dot_mem_type, Q_dot_loc));
230c798d55aSJames Wright   PetscCall(VecC2P(user->g_ceed, g_mem_type, G_loc));
23177841947SLeila Ghaffari 
23277841947SLeila Ghaffari   // Local-to-Global
2332b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(G));
2342b730f8bSJeremy L Thompson   PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G));
23577841947SLeila Ghaffari 
23677841947SLeila Ghaffari   // Restore vectors
237ca69d878SAdeleke O. Bankole   PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
23877841947SLeila Ghaffari 
23977841947SLeila Ghaffari   PetscFunctionReturn(0);
24077841947SLeila Ghaffari }
24177841947SLeila Ghaffari 
242e334ad8fSJed Brown static PetscErrorCode MatMult_NS_IJacobian(Mat J, Vec Q, Vec G) {
243e334ad8fSJed Brown   User               user;
244e334ad8fSJed Brown   const PetscScalar *q;
245e334ad8fSJed Brown   PetscScalar       *g;
246e334ad8fSJed Brown   PetscMemType       q_mem_type, g_mem_type;
247e334ad8fSJed Brown   PetscFunctionBeginUser;
2482b730f8bSJeremy L Thompson   PetscCall(MatShellGetContext(J, &user));
2495e82a6e1SJeremy L Thompson   Vec Q_loc = user->Q_dot_loc,  // Note - Q_dot_loc has zero BCs
2505e82a6e1SJeremy L Thompson       G_loc;
2515e82a6e1SJeremy L Thompson 
252e334ad8fSJed Brown   // Get local vectors
2532b730f8bSJeremy L Thompson   PetscCall(DMGetLocalVector(user->dm, &G_loc));
254e334ad8fSJed Brown 
255e334ad8fSJed Brown   // Global-to-local
2562b730f8bSJeremy L Thompson   PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
257e334ad8fSJed Brown 
258e334ad8fSJed Brown   // Place PETSc vectors in CEED vectors
2592b730f8bSJeremy L Thompson   PetscCall(VecGetArrayReadAndMemType(Q_loc, &q, &q_mem_type));
2602b730f8bSJeremy L Thompson   PetscCall(VecGetArrayAndMemType(G_loc, &g, &g_mem_type));
2612b730f8bSJeremy L Thompson   CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER, (PetscScalar *)q);
262e334ad8fSJed Brown   CeedVectorSetArray(user->g_ceed, MemTypeP2C(g_mem_type), CEED_USE_POINTER, g);
263e334ad8fSJed Brown 
264e334ad8fSJed Brown   // Apply CEED operator
2652b730f8bSJeremy L Thompson   CeedOperatorApply(user->op_ijacobian, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE);
266e334ad8fSJed Brown 
267e334ad8fSJed Brown   // Restore vectors
268e334ad8fSJed Brown   CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL);
269e334ad8fSJed Brown   CeedVectorTakeArray(user->g_ceed, MemTypeP2C(g_mem_type), NULL);
2702b730f8bSJeremy L Thompson   PetscCall(VecRestoreArrayReadAndMemType(Q_loc, &q));
2712b730f8bSJeremy L Thompson   PetscCall(VecRestoreArrayAndMemType(G_loc, &g));
272e334ad8fSJed Brown 
273e334ad8fSJed Brown   // Local-to-Global
2742b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(G));
2752b730f8bSJeremy L Thompson   PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G));
276e334ad8fSJed Brown 
277e334ad8fSJed Brown   // Restore vectors
2782b730f8bSJeremy L Thompson   PetscCall(DMRestoreLocalVector(user->dm, &G_loc));
279e334ad8fSJed Brown   PetscFunctionReturn(0);
280e334ad8fSJed Brown }
281e334ad8fSJed Brown 
282e334ad8fSJed Brown PetscErrorCode MatGetDiagonal_NS_IJacobian(Mat A, Vec D) {
283e334ad8fSJed Brown   User         user;
284e334ad8fSJed Brown   Vec          D_loc;
285e334ad8fSJed Brown   PetscScalar *d;
286e334ad8fSJed Brown   PetscMemType mem_type;
287e334ad8fSJed Brown 
288e334ad8fSJed Brown   PetscFunctionBeginUser;
289e334ad8fSJed Brown   MatShellGetContext(A, &user);
290e334ad8fSJed Brown   PetscCall(DMGetLocalVector(user->dm, &D_loc));
291e334ad8fSJed Brown   PetscCall(VecGetArrayAndMemType(D_loc, &d, &mem_type));
292e334ad8fSJed Brown   CeedVectorSetArray(user->g_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, d);
2932b730f8bSJeremy L Thompson   CeedOperatorLinearAssembleDiagonal(user->op_ijacobian, user->g_ceed, CEED_REQUEST_IMMEDIATE);
294e334ad8fSJed Brown   CeedVectorTakeArray(user->g_ceed, MemTypeP2C(mem_type), NULL);
295e334ad8fSJed Brown   PetscCall(VecRestoreArrayAndMemType(D_loc, &d));
296e334ad8fSJed Brown   PetscCall(VecZeroEntries(D));
297e334ad8fSJed Brown   PetscCall(DMLocalToGlobal(user->dm, D_loc, ADD_VALUES, D));
298e334ad8fSJed Brown   PetscCall(DMRestoreLocalVector(user->dm, &D_loc));
299e334ad8fSJed Brown   VecViewFromOptions(D, NULL, "-diag_vec_view");
300e334ad8fSJed Brown   PetscFunctionReturn(0);
301e334ad8fSJed Brown }
302e334ad8fSJed Brown 
3032b730f8bSJeremy L Thompson static PetscErrorCode FormPreallocation(User user, PetscBool pbdiagonal, Mat J, CeedVector *coo_values) {
304544be873SJed Brown   PetscCount ncoo;
305544be873SJed Brown   PetscInt  *rows, *cols;
306544be873SJed Brown 
307544be873SJed Brown   PetscFunctionBeginUser;
308544be873SJed Brown   if (pbdiagonal) {
309544be873SJed Brown     CeedSize l_size;
310544be873SJed Brown     CeedOperatorGetActiveVectorLengths(user->op_ijacobian, &l_size, NULL);
311544be873SJed Brown     ncoo = l_size * 5;
312544be873SJed Brown     rows = malloc(ncoo * sizeof(rows[0]));
313544be873SJed Brown     cols = malloc(ncoo * sizeof(cols[0]));
314544be873SJed Brown     for (PetscCount n = 0; n < l_size / 5; n++) {
315544be873SJed Brown       for (PetscInt i = 0; i < 5; i++) {
316544be873SJed Brown         for (PetscInt j = 0; j < 5; j++) {
317544be873SJed Brown           rows[(n * 5 + i) * 5 + j] = n * 5 + i;
318544be873SJed Brown           cols[(n * 5 + i) * 5 + j] = n * 5 + j;
319544be873SJed Brown         }
320544be873SJed Brown       }
321544be873SJed Brown     }
322544be873SJed Brown   } else {
3232b730f8bSJeremy L Thompson     PetscCall(CeedOperatorLinearAssembleSymbolic(user->op_ijacobian, &ncoo, &rows, &cols));
324544be873SJed Brown   }
325544be873SJed Brown   PetscCall(MatSetPreallocationCOOLocal(J, ncoo, rows, cols));
326544be873SJed Brown   free(rows);
327544be873SJed Brown   free(cols);
328544be873SJed Brown   CeedVectorCreate(user->ceed, ncoo, coo_values);
329544be873SJed Brown   PetscFunctionReturn(0);
330544be873SJed Brown }
331544be873SJed Brown 
3322b730f8bSJeremy L Thompson static PetscErrorCode FormSetValues(User user, PetscBool pbdiagonal, Mat J, CeedVector coo_values) {
333544be873SJed Brown   CeedMemType        mem_type = CEED_MEM_HOST;
334544be873SJed Brown   const PetscScalar *values;
335544be873SJed Brown   MatType            mat_type;
336544be873SJed Brown 
337544be873SJed Brown   PetscFunctionBeginUser;
338544be873SJed Brown   PetscCall(MatGetType(J, &mat_type));
3392b730f8bSJeremy L Thompson   if (strstr(mat_type, "kokkos") || strstr(mat_type, "cusparse")) mem_type = CEED_MEM_DEVICE;
340544be873SJed Brown   if (user->app_ctx->pmat_pbdiagonal) {
3412b730f8bSJeremy L Thompson     CeedOperatorLinearAssemblePointBlockDiagonal(user->op_ijacobian, coo_values, CEED_REQUEST_IMMEDIATE);
342544be873SJed Brown   } else {
343544be873SJed Brown     CeedOperatorLinearAssemble(user->op_ijacobian, coo_values);
344544be873SJed Brown   }
345544be873SJed Brown   CeedVectorGetArrayRead(coo_values, mem_type, &values);
346544be873SJed Brown   PetscCall(MatSetValuesCOO(J, values, INSERT_VALUES));
347544be873SJed Brown   CeedVectorRestoreArrayRead(coo_values, &values);
348544be873SJed Brown   PetscFunctionReturn(0);
349544be873SJed Brown }
350544be873SJed Brown 
3512b730f8bSJeremy L Thompson PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, PetscReal shift, Mat J, Mat J_pre, void *user_data) {
352e334ad8fSJed Brown   User      user = *(User *)user_data;
353d6bf345cSJed Brown   PetscBool J_is_shell, J_is_mffd, J_pre_is_shell;
354e334ad8fSJed Brown   PetscFunctionBeginUser;
35517b0d5c6SJeremy L Thompson   if (user->phys->ijacobian_time_shift_label) CeedOperatorSetContextDouble(user->op_ijacobian, user->phys->ijacobian_time_shift_label, &shift);
356d6bf345cSJed Brown   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATMFFD, &J_is_mffd));
357e334ad8fSJed Brown   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATSHELL, &J_is_shell));
3582b730f8bSJeremy L Thompson   PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATSHELL, &J_pre_is_shell));
359e334ad8fSJed Brown   if (!user->matrices_set_up) {
360e334ad8fSJed Brown     if (J_is_shell) {
361e334ad8fSJed Brown       PetscCall(MatShellSetContext(J, user));
3622b730f8bSJeremy L Thompson       PetscCall(MatShellSetOperation(J, MATOP_MULT, (void (*)(void))MatMult_NS_IJacobian));
3632b730f8bSJeremy L Thompson       PetscCall(MatShellSetOperation(J, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_NS_IJacobian));
364e334ad8fSJed Brown       PetscCall(MatSetUp(J));
365e334ad8fSJed Brown     }
366e334ad8fSJed Brown     if (!J_pre_is_shell) {
3672b730f8bSJeremy L Thompson       PetscCall(FormPreallocation(user, user->app_ctx->pmat_pbdiagonal, J_pre, &user->coo_values_pmat));
368544be873SJed Brown     }
369d6bf345cSJed Brown     if (J != J_pre && !J_is_shell && !J_is_mffd) {
370544be873SJed Brown       PetscCall(FormPreallocation(user, PETSC_FALSE, J, &user->coo_values_amat));
371544be873SJed Brown     }
372e334ad8fSJed Brown     user->matrices_set_up = true;
373e334ad8fSJed Brown   }
374e334ad8fSJed Brown   if (!J_pre_is_shell) {
3752b730f8bSJeremy L Thompson     PetscCall(FormSetValues(user, user->app_ctx->pmat_pbdiagonal, J_pre, user->coo_values_pmat));
376e334ad8fSJed Brown   }
377d6bf345cSJed Brown   if (user->coo_values_amat) {
378d6bf345cSJed Brown     PetscCall(FormSetValues(user, PETSC_FALSE, J, user->coo_values_amat));
379d6bf345cSJed Brown   } else if (J_is_mffd) {
380d6bf345cSJed Brown     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
381d6bf345cSJed Brown     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
382d6bf345cSJed Brown   }
383e334ad8fSJed Brown   PetscFunctionReturn(0);
384e334ad8fSJed Brown }
385e334ad8fSJed Brown 
3862b730f8bSJeremy L Thompson PetscErrorCode WriteOutput(User user, Vec Q, PetscInt step_no, PetscScalar time) {
38777841947SLeila Ghaffari   Vec         Q_loc;
38877841947SLeila Ghaffari   char        file_path[PETSC_MAX_PATH_LEN];
38977841947SLeila Ghaffari   PetscViewer viewer;
39077841947SLeila Ghaffari   PetscFunctionBeginUser;
39177841947SLeila Ghaffari 
39237cbb16aSJed Brown   if (user->app_ctx->checkpoint_vtk) {
39377841947SLeila Ghaffari     // Set up output
394f14660a4SJames Wright     PetscCall(DMGetLocalVector(user->dm, &Q_loc));
395f14660a4SJames Wright     PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec"));
396f14660a4SJames Wright     PetscCall(VecZeroEntries(Q_loc));
397f14660a4SJames Wright     PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
39877841947SLeila Ghaffari 
39977841947SLeila Ghaffari     // Output
40037cbb16aSJed Brown     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
401f14660a4SJames Wright 
4022b730f8bSJeremy L Thompson     PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path, FILE_MODE_WRITE, &viewer));
403f14660a4SJames Wright     PetscCall(VecView(Q_loc, viewer));
404f14660a4SJames Wright     PetscCall(PetscViewerDestroy(&viewer));
40577841947SLeila Ghaffari     if (user->dm_viz) {
40677841947SLeila Ghaffari       Vec         Q_refined, Q_refined_loc;
40777841947SLeila Ghaffari       char        file_path_refined[PETSC_MAX_PATH_LEN];
40877841947SLeila Ghaffari       PetscViewer viewer_refined;
40977841947SLeila Ghaffari 
410f14660a4SJames Wright       PetscCall(DMGetGlobalVector(user->dm_viz, &Q_refined));
411f14660a4SJames Wright       PetscCall(DMGetLocalVector(user->dm_viz, &Q_refined_loc));
412f14660a4SJames Wright       PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined"));
413f14660a4SJames Wright 
414f14660a4SJames Wright       PetscCall(MatInterpolate(user->interp_viz, Q, Q_refined));
415f14660a4SJames Wright       PetscCall(VecZeroEntries(Q_refined_loc));
4162b730f8bSJeremy L Thompson       PetscCall(DMGlobalToLocal(user->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc));
417f14660a4SJames Wright 
41837cbb16aSJed Brown       PetscCall(
41937cbb16aSJed Brown           PetscSNPrintf(file_path_refined, sizeof file_path_refined, "%s/nsrefined-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
420f14660a4SJames Wright 
4212b730f8bSJeremy L Thompson       PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined), file_path_refined, FILE_MODE_WRITE, &viewer_refined));
422f14660a4SJames Wright       PetscCall(VecView(Q_refined_loc, viewer_refined));
423f14660a4SJames Wright       PetscCall(DMRestoreLocalVector(user->dm_viz, &Q_refined_loc));
424f14660a4SJames Wright       PetscCall(DMRestoreGlobalVector(user->dm_viz, &Q_refined));
425f14660a4SJames Wright       PetscCall(PetscViewerDestroy(&viewer_refined));
42677841947SLeila Ghaffari     }
427f14660a4SJames Wright     PetscCall(DMRestoreLocalVector(user->dm, &Q_loc));
42837cbb16aSJed Brown   }
42977841947SLeila Ghaffari 
43077841947SLeila Ghaffari   // Save data in a binary file for continuation of simulations
43169293791SJames Wright   if (user->app_ctx->add_stepnum2bin) {
43237cbb16aSJed Brown     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", user->app_ctx->output_dir, step_no));
43369293791SJames Wright   } else {
4342b730f8bSJeremy L Thompson     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", user->app_ctx->output_dir));
43569293791SJames Wright   }
4362b730f8bSJeremy L Thompson   PetscCall(PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, &viewer));
437f14660a4SJames Wright 
4384de8550aSJed Brown   PetscInt token = FLUIDS_FILE_TOKEN;
4394de8550aSJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &token, 1, PETSC_INT));
4404de8550aSJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &step_no, 1, PETSC_INT));
4414de8550aSJed Brown   time /= user->units->second;  // Dimensionalize time back
4424de8550aSJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL));
443f14660a4SJames Wright   PetscCall(VecView(Q, viewer));
444f14660a4SJames Wright   PetscCall(PetscViewerDestroy(&viewer));
445f14660a4SJames Wright   PetscFunctionReturn(0);
446f14660a4SJames Wright }
447f14660a4SJames Wright 
448ca69d878SAdeleke O. Bankole // CSV Monitor
449ca69d878SAdeleke O. Bankole PetscErrorCode TSMonitor_WallForce(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
450ca69d878SAdeleke O. Bankole   User              user = ctx;
451ca69d878SAdeleke O. Bankole   Vec               G_loc;
452ca69d878SAdeleke O. Bankole   PetscInt          num_wall = user->app_ctx->wall_forces.num_wall, dim = 3;
453ca69d878SAdeleke O. Bankole   const PetscInt   *walls  = user->app_ctx->wall_forces.walls;
454ca69d878SAdeleke O. Bankole   PetscViewer       viewer = user->app_ctx->wall_forces.viewer;
455ca69d878SAdeleke O. Bankole   PetscViewerFormat format = user->app_ctx->wall_forces.viewer_format;
456ca69d878SAdeleke O. Bankole   PetscScalar      *reaction_force;
457ca69d878SAdeleke O. Bankole   PetscBool         iascii;
458ca69d878SAdeleke O. Bankole 
459ca69d878SAdeleke O. Bankole   PetscFunctionBeginUser;
460ca69d878SAdeleke O. Bankole   if (!viewer) PetscFunctionReturn(0);
461ca69d878SAdeleke O. Bankole   PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
462ca69d878SAdeleke O. Bankole   PetscCall(PetscMalloc1(num_wall * dim, &reaction_force));
463ca69d878SAdeleke O. Bankole   PetscCall(Surface_Forces_NS(user->dm, G_loc, num_wall, walls, reaction_force));
464ca69d878SAdeleke O. Bankole   PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
465ca69d878SAdeleke O. Bankole 
466ca69d878SAdeleke O. Bankole   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
467ca69d878SAdeleke O. Bankole 
468ca69d878SAdeleke O. Bankole   if (iascii) {
469ca69d878SAdeleke O. Bankole     if (format == PETSC_VIEWER_ASCII_CSV && !user->app_ctx->wall_forces.header_written) {
470ca69d878SAdeleke O. Bankole       PetscCall(PetscViewerASCIIPrintf(viewer, "Step,Time,Wall,ForceX,ForceY,ForceZ\n"));
471ca69d878SAdeleke O. Bankole       user->app_ctx->wall_forces.header_written = PETSC_TRUE;
472ca69d878SAdeleke O. Bankole     }
473ca69d878SAdeleke O. Bankole     for (PetscInt w = 0; w < num_wall; w++) {
474ca69d878SAdeleke O. Bankole       PetscInt wall = walls[w];
475ca69d878SAdeleke O. Bankole       if (format == PETSC_VIEWER_ASCII_CSV) {
476ca69d878SAdeleke O. Bankole         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT ",%g,%" PetscInt_FMT ",%g,%g,%g\n", step_no, time, wall,
477ca69d878SAdeleke O. Bankole                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
478ca69d878SAdeleke O. Bankole 
479ca69d878SAdeleke O. Bankole       } else {
480ca69d878SAdeleke O. Bankole         PetscCall(PetscViewerASCIIPrintf(viewer, "Wall %" PetscInt_FMT " Forces: Force_x = %12g, Force_y = %12g, Force_z = %12g\n", wall,
481ca69d878SAdeleke O. Bankole                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
482ca69d878SAdeleke O. Bankole       }
483ca69d878SAdeleke O. Bankole     }
484ca69d878SAdeleke O. Bankole   }
485ca69d878SAdeleke O. Bankole   PetscCall(PetscFree(reaction_force));
486ca69d878SAdeleke O. Bankole   PetscFunctionReturn(0);
487ca69d878SAdeleke O. Bankole }
488ca69d878SAdeleke O. Bankole 
489f14660a4SJames Wright // User provided TS Monitor
4902b730f8bSJeremy L Thompson PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
491f14660a4SJames Wright   User user = ctx;
492f14660a4SJames Wright   PetscFunctionBeginUser;
493f14660a4SJames Wright 
49437cbb16aSJed Brown   // Print every 'checkpoint_interval' steps
495894de27cSJames Wright   if (user->app_ctx->checkpoint_interval <= 0 || step_no % user->app_ctx->checkpoint_interval != 0 ||
496894de27cSJames Wright       (user->app_ctx->cont_steps == step_no && step_no != 0))
497894de27cSJames Wright     PetscFunctionReturn(0);
498f14660a4SJames Wright 
499f14660a4SJames Wright   PetscCall(WriteOutput(user, Q, step_no, time));
50077841947SLeila Ghaffari 
50177841947SLeila Ghaffari   PetscFunctionReturn(0);
50277841947SLeila Ghaffari }
50377841947SLeila Ghaffari 
50477841947SLeila Ghaffari // TS: Create, setup, and solve
5052b730f8bSJeremy L Thompson PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, Vec *Q, PetscScalar *f_time, TS *ts) {
50677841947SLeila Ghaffari   MPI_Comm    comm = user->comm;
50777841947SLeila Ghaffari   TSAdapt     adapt;
50877841947SLeila Ghaffari   PetscScalar final_time;
50977841947SLeila Ghaffari   PetscFunctionBeginUser;
51077841947SLeila Ghaffari 
5112b730f8bSJeremy L Thompson   PetscCall(TSCreate(comm, ts));
5122b730f8bSJeremy L Thompson   PetscCall(TSSetDM(*ts, dm));
51377841947SLeila Ghaffari   if (phys->implicit) {
5142b730f8bSJeremy L Thompson     PetscCall(TSSetType(*ts, TSBDF));
51577841947SLeila Ghaffari     if (user->op_ifunction) {
5162b730f8bSJeremy L Thompson       PetscCall(TSSetIFunction(*ts, NULL, IFunction_NS, &user));
51777841947SLeila Ghaffari     } else {  // Implicit integrators can fall back to using an RHSFunction
5182b730f8bSJeremy L Thompson       PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
51977841947SLeila Ghaffari     }
520e334ad8fSJed Brown     if (user->op_ijacobian) {
5212b730f8bSJeremy L Thompson       PetscCall(DMTSSetIJacobian(dm, FormIJacobian_NS, &user));
522544be873SJed Brown       if (app_ctx->amat_type) {
523544be873SJed Brown         Mat Pmat, Amat;
5242b730f8bSJeremy L Thompson         PetscCall(DMCreateMatrix(dm, &Pmat));
5252b730f8bSJeremy L Thompson         PetscCall(DMSetMatType(dm, app_ctx->amat_type));
5262b730f8bSJeremy L Thompson         PetscCall(DMCreateMatrix(dm, &Amat));
5272b730f8bSJeremy L Thompson         PetscCall(TSSetIJacobian(*ts, Amat, Pmat, NULL, NULL));
5282b730f8bSJeremy L Thompson         PetscCall(MatDestroy(&Amat));
5292b730f8bSJeremy L Thompson         PetscCall(MatDestroy(&Pmat));
530544be873SJed Brown       }
531e334ad8fSJed Brown     }
53277841947SLeila Ghaffari   } else {
5332b730f8bSJeremy L Thompson     if (!user->op_rhs) SETERRQ(comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction");
5342b730f8bSJeremy L Thompson     PetscCall(TSSetType(*ts, TSRK));
5352b730f8bSJeremy L Thompson     PetscCall(TSRKSetType(*ts, TSRK5F));
5362b730f8bSJeremy L Thompson     PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
53777841947SLeila Ghaffari   }
5382b730f8bSJeremy L Thompson   PetscCall(TSSetMaxTime(*ts, 500. * user->units->second));
5392b730f8bSJeremy L Thompson   PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER));
540cf7a0454SJed Brown   PetscCall(TSSetErrorIfStepFails(*ts, PETSC_FALSE));
5412b730f8bSJeremy L Thompson   PetscCall(TSSetTimeStep(*ts, 1.e-2 * user->units->second));
5428fb33541SJames Wright   if (app_ctx->test_type != TESTTYPE_NONE) {
5432b730f8bSJeremy L Thompson     PetscCall(TSSetMaxSteps(*ts, 10));
5442b730f8bSJeremy L Thompson   }
5452b730f8bSJeremy L Thompson   PetscCall(TSGetAdapt(*ts, &adapt));
5462b730f8bSJeremy L Thompson   PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * user->units->second, 1.e2 * user->units->second));
5472b730f8bSJeremy L Thompson   PetscCall(TSSetFromOptions(*ts));
548*15637395SJames Wright   user->time_bc_set = -1.0;    // require all BCs be updated
54977841947SLeila Ghaffari   if (!app_ctx->cont_steps) {  // print initial condition
5508fb33541SJames Wright     if (app_ctx->test_type == TESTTYPE_NONE) {
5512b730f8bSJeremy L Thompson       PetscCall(TSMonitor_NS(*ts, 0, 0., *Q, user));
55277841947SLeila Ghaffari     }
55377841947SLeila Ghaffari   } else {  // continue from time of last output
55477841947SLeila Ghaffari     PetscInt    count;
55577841947SLeila Ghaffari     PetscViewer viewer;
5562b730f8bSJeremy L Thompson 
5574de8550aSJed Brown     if (app_ctx->cont_time <= 0) {  // Legacy files did not include step number and time
5582b730f8bSJeremy L Thompson       PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_time_file, FILE_MODE_READ, &viewer));
5594de8550aSJed Brown       PetscCall(PetscViewerBinaryRead(viewer, &app_ctx->cont_time, 1, &count, PETSC_REAL));
5602b730f8bSJeremy L Thompson       PetscCall(PetscViewerDestroy(&viewer));
5614de8550aSJed Brown       PetscCheck(app_ctx->cont_steps != -1, comm, PETSC_ERR_ARG_INCOMP,
5624de8550aSJed Brown                  "-continue step number not specified, but checkpoint file does not contain a step number (likely written by older code version)");
5634de8550aSJed Brown     }
5644de8550aSJed Brown     PetscCall(TSSetTime(*ts, app_ctx->cont_time * user->units->second));
565d8e0aecdSJed Brown     PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps));
56677841947SLeila Ghaffari   }
5678fb33541SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) {
5682b730f8bSJeremy L Thompson     PetscCall(TSMonitorSet(*ts, TSMonitor_NS, user, NULL));
5698fb33541SJames Wright   }
570ca69d878SAdeleke O. Bankole   if (app_ctx->wall_forces.viewer) {
571ca69d878SAdeleke O. Bankole     PetscCall(TSMonitorSet(*ts, TSMonitor_WallForce, user, NULL));
572ca69d878SAdeleke O. Bankole   }
573b7d66439SJames Wright   if (app_ctx->turb_spanstats_enable) {
574a175e481SJames Wright     PetscCall(TSMonitorSet(*ts, TSMonitor_Statistics, user, NULL));
575495b9769SJames Wright     CeedScalar previous_time = app_ctx->cont_time * user->units->second;
57617b0d5c6SJeremy L Thompson     CeedOperatorSetContextDouble(user->spanstats.op_stats_collect, user->spanstats.previous_time_label, &previous_time);
577a175e481SJames Wright   }
57877841947SLeila Ghaffari 
57977841947SLeila Ghaffari   // Solve
580d8e0aecdSJed Brown   PetscReal start_time;
581d8e0aecdSJed Brown   PetscInt  start_step;
5822b730f8bSJeremy L Thompson   PetscCall(TSGetTime(*ts, &start_time));
583d8e0aecdSJed Brown   PetscCall(TSGetStepNumber(*ts, &start_step));
5847b39487dSJeremy L Thompson 
5857b39487dSJeremy L Thompson   PetscPreLoadBegin(PETSC_FALSE, "Fluids Solve");
5867b39487dSJeremy L Thompson   PetscCall(TSSetTime(*ts, start_time));
587d8e0aecdSJed Brown   PetscCall(TSSetStepNumber(*ts, start_step));
5887b39487dSJeremy L Thompson   if (PetscPreLoadingOn) {
5897b39487dSJeremy L Thompson     // LCOV_EXCL_START
5907b39487dSJeremy L Thompson     SNES      snes;
5917b39487dSJeremy L Thompson     Vec       Q_preload;
5927b39487dSJeremy L Thompson     PetscReal rtol;
5937b39487dSJeremy L Thompson     PetscCall(VecDuplicate(*Q, &Q_preload));
5947b39487dSJeremy L Thompson     PetscCall(VecCopy(*Q, Q_preload));
5957b39487dSJeremy L Thompson     PetscCall(TSGetSNES(*ts, &snes));
5967b39487dSJeremy L Thompson     PetscCall(SNESGetTolerances(snes, NULL, &rtol, NULL, NULL, NULL));
5972b730f8bSJeremy L Thompson     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, .99, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
5987b39487dSJeremy L Thompson     PetscCall(TSSetSolution(*ts, *Q));
5997b39487dSJeremy L Thompson     PetscCall(TSStep(*ts));
6002b730f8bSJeremy L Thompson     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, rtol, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
6017b39487dSJeremy L Thompson     PetscCall(VecDestroy(&Q_preload));
6027b39487dSJeremy L Thompson     // LCOV_EXCL_STOP
6037b39487dSJeremy L Thompson   } else {
6042b730f8bSJeremy L Thompson     PetscCall(PetscBarrier((PetscObject)*ts));
6052b730f8bSJeremy L Thompson     PetscCall(TSSolve(*ts, *Q));
6067b39487dSJeremy L Thompson   }
6077b39487dSJeremy L Thompson   PetscPreLoadEnd();
6087b39487dSJeremy L Thompson 
6097b39487dSJeremy L Thompson   PetscCall(TSGetSolveTime(*ts, &final_time));
61077841947SLeila Ghaffari   *f_time = final_time;
6117b39487dSJeremy L Thompson 
6128fb33541SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) {
613f14660a4SJames Wright     PetscInt step_no;
614f14660a4SJames Wright     PetscCall(TSGetStepNumber(*ts, &step_no));
615a175e481SJames Wright     if (user->app_ctx->checkpoint_interval > 0 || user->app_ctx->checkpoint_interval == -1) {
616f14660a4SJames Wright       PetscCall(WriteOutput(user, *Q, step_no, final_time));
617f14660a4SJames Wright     }
618f14660a4SJames Wright 
6197b39487dSJeremy L Thompson     PetscLogEvent stage_id;
6207b39487dSJeremy L Thompson     PetscStageLog stage_log;
6217b39487dSJeremy L Thompson 
6227b39487dSJeremy L Thompson     PetscCall(PetscLogStageGetId("Fluids Solve", &stage_id));
6237b39487dSJeremy L Thompson     PetscCall(PetscLogGetStageLog(&stage_log));
6242b730f8bSJeremy L Thompson     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time taken for solution (sec): %g\n", stage_log->stageInfo[stage_id].perfInfo.time));
62577841947SLeila Ghaffari   }
62677841947SLeila Ghaffari   PetscFunctionReturn(0);
62777841947SLeila Ghaffari }
628