xref: /honee/src/setupts.c (revision cb315d14a5f51fe9f0b5e1e4bb6d35e687e55c6e)
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 
11e419654dSJeremy L Thompson #include <ceed.h>
12e419654dSJeremy L Thompson #include <petscdmplex.h>
13e419654dSJeremy L Thompson #include <petscts.h>
14e419654dSJeremy L Thompson 
15a515125bSLeila Ghaffari #include "../navierstokes.h"
16c5e9980aSAdeleke O. Bankole #include "../qfunctions/newtonian_state.h"
17a515125bSLeila Ghaffari 
18a515125bSLeila Ghaffari // Compute mass matrix for explicit scheme
192b916ea7SJeremy L Thompson PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm, CeedData ceed_data, Vec M) {
20a515125bSLeila Ghaffari   CeedQFunction        qf_mass;
21a515125bSLeila Ghaffari   CeedOperator         op_mass;
220143e3daSJames Wright   OperatorApplyContext op_mass_ctx;
230143e3daSJames Wright   Vec                  Ones_loc;
24a515125bSLeila Ghaffari   CeedInt              num_comp_q, q_data_size;
25a515125bSLeila Ghaffari   PetscFunctionBeginUser;
26a515125bSLeila Ghaffari 
27a515125bSLeila Ghaffari   // CEED Restriction
28a515125bSLeila Ghaffari   CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q);
29a515125bSLeila Ghaffari   CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &q_data_size);
30a515125bSLeila Ghaffari 
31a515125bSLeila Ghaffari   // CEED QFunction
329f59f36eSJames Wright   PetscCall(CreateMassQFunction(ceed, num_comp_q, q_data_size, &qf_mass));
33a515125bSLeila Ghaffari 
34a515125bSLeila Ghaffari   // CEED Operator
35a515125bSLeila Ghaffari   CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass);
369f59f36eSJames Wright   CeedOperatorSetField(op_mass, "u", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
372b916ea7SJeremy L Thompson   CeedOperatorSetField(op_mass, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data);
382b916ea7SJeremy L Thompson   CeedOperatorSetField(op_mass, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
39a515125bSLeila Ghaffari 
400143e3daSJames Wright   PetscCall(OperatorApplyContextCreate(NULL, dm, ceed, op_mass, NULL, NULL, NULL, NULL, &op_mass_ctx));
41a515125bSLeila Ghaffari 
420143e3daSJames Wright   PetscCall(DMGetLocalVector(dm, &Ones_loc));
430143e3daSJames Wright   PetscCall(VecSet(Ones_loc, 1));
440143e3daSJames Wright   PetscCall(ApplyCeedOperatorLocalToGlobal(Ones_loc, M, op_mass_ctx));
45a515125bSLeila Ghaffari 
46a515125bSLeila Ghaffari   // Invert diagonally lumped mass vector for RHS function
472b916ea7SJeremy L Thompson   PetscCall(VecReciprocal(M));
48a515125bSLeila Ghaffari 
49a515125bSLeila Ghaffari   // Cleanup
500143e3daSJames Wright   PetscCall(OperatorApplyContextDestroy(op_mass_ctx));
510143e3daSJames Wright   PetscCall(DMRestoreLocalVector(dm, &Ones_loc));
52a515125bSLeila Ghaffari   CeedQFunctionDestroy(&qf_mass);
53a515125bSLeila Ghaffari   CeedOperatorDestroy(&op_mass);
54a515125bSLeila Ghaffari 
55a515125bSLeila Ghaffari   PetscFunctionReturn(0);
56a515125bSLeila Ghaffari }
57a515125bSLeila Ghaffari 
58c996854bSJames Wright // Insert Boundary values if it's a new time
59c996854bSJames Wright PetscErrorCode UpdateBoundaryValues(User user, Vec Q_loc, PetscReal t) {
60c996854bSJames Wright   PetscFunctionBeginUser;
61c996854bSJames Wright   if (user->time_bc_set != t) {
62c996854bSJames Wright     PetscCall(DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL));
63c996854bSJames Wright     user->time_bc_set = t;
64c996854bSJames Wright   }
65c996854bSJames Wright   PetscFunctionReturn(0);
66c996854bSJames Wright }
67c996854bSJames Wright 
6891f639d2SJames Wright // @brief Update the context label value to new value if necessary.
6991f639d2SJames Wright // @note This only supports labels with scalar label values (ie. not arrays)
7091f639d2SJames Wright PetscErrorCode UpdateContextLabel(MPI_Comm comm, PetscScalar update_value, CeedOperator op, CeedContextFieldLabel label) {
7191f639d2SJames Wright   PetscScalar label_value;
72c2cb7fc8SJames Wright 
7391f639d2SJames Wright   PetscFunctionBeginUser;
7491f639d2SJames Wright   PetscCheck(label, comm, PETSC_ERR_ARG_BADPTR, "Label should be non-NULL");
7591f639d2SJames Wright 
7691f639d2SJames Wright   {
7791f639d2SJames Wright     size_t             num_elements;
7891f639d2SJames Wright     const PetscScalar *label_values;
7991f639d2SJames Wright     CeedOperatorGetContextDoubleRead(op, label, &num_elements, &label_values);
8091f639d2SJames Wright     PetscCheck(num_elements == 1, comm, PETSC_ERR_SUP, "%s does not support labels with more than 1 value. Label has %zu values", __func__,
8191f639d2SJames Wright                num_elements);
8291f639d2SJames Wright     label_value = *label_values;
8391f639d2SJames Wright     CeedOperatorRestoreContextDoubleRead(op, label, &label_values);
84c2cb7fc8SJames Wright   }
8591f639d2SJames Wright 
8691f639d2SJames Wright   if (label_value != update_value) {
8791f639d2SJames Wright     CeedOperatorSetContextDouble(op, label, &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;
9791f639d2SJames Wright   MPI_Comm    comm = PetscObjectComm((PetscObject)ts);
98fd969b44SJames Wright   PetscScalar dt;
99da5fe0e4SJames Wright   Vec         Q_loc = user->Q_loc;
100a515125bSLeila Ghaffari   PetscFunctionBeginUser;
101a515125bSLeila Ghaffari 
102e2f84137SJeremy L Thompson   // Update time dependent data
103c996854bSJames Wright   PetscCall(UpdateBoundaryValues(user, Q_loc, t));
104da5fe0e4SJames Wright   if (user->phys->solution_time_label) PetscCall(UpdateContextLabel(comm, t, user->op_rhs_ctx->op, user->phys->solution_time_label));
1052b916ea7SJeremy L Thompson   PetscCall(TSGetTimeStep(ts, &dt));
106da5fe0e4SJames Wright   if (user->phys->timestep_size_label) PetscCall(UpdateContextLabel(comm, dt, user->op_rhs_ctx->op, user->phys->timestep_size_label));
107a515125bSLeila Ghaffari 
108da5fe0e4SJames Wright   PetscCall(ApplyCeedOperatorGlobalToGlobal(Q, G, user->op_rhs_ctx));
109a515125bSLeila Ghaffari 
110a515125bSLeila Ghaffari   // Inverse of the lumped mass matrix (M is Minv)
111cc9aa765SJames Wright   PetscCall(VecPointwiseMult(G, G, user->M_inv));
112a515125bSLeila Ghaffari 
113a515125bSLeila Ghaffari   PetscFunctionReturn(0);
114a515125bSLeila Ghaffari }
115a515125bSLeila Ghaffari 
116c5e9980aSAdeleke O. Bankole // Surface forces function setup
117c5e9980aSAdeleke O. Bankole static PetscErrorCode Surface_Forces_NS(DM dm, Vec G_loc, PetscInt num_walls, const PetscInt walls[], PetscScalar *reaction_force) {
118c5e9980aSAdeleke O. Bankole   DMLabel            face_label;
119c5e9980aSAdeleke O. Bankole   const PetscScalar *g;
1202004e3acSAdeleke O. Bankole   PetscInt           dof, dim = 3;
121c5e9980aSAdeleke O. Bankole   MPI_Comm           comm;
1222004e3acSAdeleke O. Bankole   PetscSection       s;
123c5e9980aSAdeleke O. Bankole 
124c5e9980aSAdeleke O. Bankole   PetscFunctionBeginUser;
125c5e9980aSAdeleke O. Bankole   PetscCall(PetscArrayzero(reaction_force, num_walls * dim));
126c5e9980aSAdeleke O. Bankole   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
127c5e9980aSAdeleke O. Bankole   PetscCall(DMGetLabel(dm, "Face Sets", &face_label));
128c5e9980aSAdeleke O. Bankole   PetscCall(VecGetArrayRead(G_loc, &g));
129c5e9980aSAdeleke O. Bankole   for (PetscInt w = 0; w < num_walls; w++) {
130c5e9980aSAdeleke O. Bankole     const PetscInt wall = walls[w];
131c5e9980aSAdeleke O. Bankole     IS             wall_is;
1322004e3acSAdeleke O. Bankole     PetscCall(DMGetLocalSection(dm, &s));
133c5e9980aSAdeleke O. Bankole     PetscCall(DMLabelGetStratumIS(face_label, wall, &wall_is));
134c5e9980aSAdeleke O. Bankole     if (wall_is) {  // There exist such points on this process
135c5e9980aSAdeleke O. Bankole       PetscInt        num_points;
1362004e3acSAdeleke O. Bankole       PetscInt        num_comp = 0;
137c5e9980aSAdeleke O. Bankole       const PetscInt *points;
1382004e3acSAdeleke O. Bankole       PetscCall(PetscSectionGetFieldComponents(s, 0, &num_comp));
139c5e9980aSAdeleke O. Bankole       PetscCall(ISGetSize(wall_is, &num_points));
140c5e9980aSAdeleke O. Bankole       PetscCall(ISGetIndices(wall_is, &points));
141c5e9980aSAdeleke O. Bankole       for (PetscInt i = 0; i < num_points; i++) {
142c5e9980aSAdeleke O. Bankole         const PetscInt           p = points[i];
143c5e9980aSAdeleke O. Bankole         const StateConservative *r;
144c5e9980aSAdeleke O. Bankole         PetscCall(DMPlexPointLocalRead(dm, p, g, &r));
1452004e3acSAdeleke O. Bankole         PetscCall(PetscSectionGetDof(s, p, &dof));
1462004e3acSAdeleke O. Bankole         for (PetscInt node = 0; node < dof / num_comp; node++) {
147c5e9980aSAdeleke O. Bankole           for (PetscInt j = 0; j < 3; j++) {
1482004e3acSAdeleke O. Bankole             reaction_force[w * dim + j] -= r[node].momentum[j];
1492004e3acSAdeleke O. Bankole           }
150c5e9980aSAdeleke O. Bankole         }
151c5e9980aSAdeleke O. Bankole       }
152c5e9980aSAdeleke O. Bankole       PetscCall(ISRestoreIndices(wall_is, &points));
153c5e9980aSAdeleke O. Bankole     }
154c5e9980aSAdeleke O. Bankole     PetscCall(ISDestroy(&wall_is));
155c5e9980aSAdeleke O. Bankole   }
156c5e9980aSAdeleke O. Bankole   PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, reaction_force, dim * num_walls, MPIU_SCALAR, MPI_SUM, comm));
157c5e9980aSAdeleke O. Bankole   //  Restore Vectors
158c5e9980aSAdeleke O. Bankole   PetscCall(VecRestoreArrayRead(G_loc, &g));
159c5e9980aSAdeleke O. Bankole 
160c5e9980aSAdeleke O. Bankole   PetscFunctionReturn(0);
161c5e9980aSAdeleke O. Bankole }
162c5e9980aSAdeleke O. Bankole 
163a515125bSLeila Ghaffari // Implicit time-stepper function setup
1642b916ea7SJeremy L Thompson PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) {
165a515125bSLeila Ghaffari   User         user = *(User *)user_data;
16691f639d2SJames Wright   MPI_Comm     comm = PetscObjectComm((PetscObject)ts);
167fd969b44SJames Wright   PetscScalar  dt;
168e2f84137SJeremy L Thompson   Vec          Q_loc = user->Q_loc, Q_dot_loc = user->Q_dot_loc, G_loc;
169a515125bSLeila Ghaffari   PetscMemType q_mem_type, q_dot_mem_type, g_mem_type;
170a515125bSLeila Ghaffari   PetscFunctionBeginUser;
171a515125bSLeila Ghaffari 
172e2f84137SJeremy L Thompson   // Get local vectors
173c5e9980aSAdeleke O. Bankole   PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
174e2f84137SJeremy L Thompson 
175e2f84137SJeremy L Thompson   // Update time dependent data
176c996854bSJames Wright   PetscCall(UpdateBoundaryValues(user, Q_loc, t));
17791f639d2SJames Wright   if (user->phys->solution_time_label) PetscCall(UpdateContextLabel(comm, t, user->op_ifunction, user->phys->solution_time_label));
1782b916ea7SJeremy L Thompson   PetscCall(TSGetTimeStep(ts, &dt));
17991f639d2SJames Wright   if (user->phys->timestep_size_label) PetscCall(UpdateContextLabel(comm, dt, user->op_ifunction, user->phys->timestep_size_label));
180a515125bSLeila Ghaffari 
181a515125bSLeila Ghaffari   // Global-to-local
18206108310SJames Wright   PetscCall(DMGlobalToLocalBegin(user->dm, Q, INSERT_VALUES, Q_loc));
18306108310SJames Wright   PetscCall(DMGlobalToLocalBegin(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
18406108310SJames Wright   PetscCall(DMGlobalToLocalEnd(user->dm, Q, INSERT_VALUES, Q_loc));
18506108310SJames Wright   PetscCall(DMGlobalToLocalEnd(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
186a515125bSLeila Ghaffari 
187a515125bSLeila Ghaffari   // Place PETSc vectors in CEED vectors
188fd969b44SJames Wright   PetscCall(VecReadP2C(Q_loc, &q_mem_type, user->q_ceed));
189fd969b44SJames Wright   PetscCall(VecReadP2C(Q_dot_loc, &q_dot_mem_type, user->q_dot_ceed));
190fd969b44SJames Wright   PetscCall(VecP2C(G_loc, &g_mem_type, user->g_ceed));
191a515125bSLeila Ghaffari 
192a515125bSLeila Ghaffari   // Apply CEED operator
1932b916ea7SJeremy L Thompson   CeedOperatorApply(user->op_ifunction, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE);
194a515125bSLeila Ghaffari 
195a515125bSLeila Ghaffari   // Restore vectors
196fd969b44SJames Wright   PetscCall(VecReadC2P(user->q_ceed, q_mem_type, Q_loc));
197fd969b44SJames Wright   PetscCall(VecReadC2P(user->q_dot_ceed, q_dot_mem_type, Q_dot_loc));
198fd969b44SJames Wright   PetscCall(VecC2P(user->g_ceed, g_mem_type, G_loc));
199a515125bSLeila Ghaffari 
20001ab89c1SJames Wright   if (user->app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) {
20101ab89c1SJames Wright     PetscCall(SGS_DD_ModelApplyIFunction(user, Q_loc, G_loc));
20201ab89c1SJames Wright   }
2039c678832SJames Wright 
204a515125bSLeila Ghaffari   // Local-to-Global
2052b916ea7SJeremy L Thompson   PetscCall(VecZeroEntries(G));
2062b916ea7SJeremy L Thompson   PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G));
207a515125bSLeila Ghaffari 
208a515125bSLeila Ghaffari   // Restore vectors
209c5e9980aSAdeleke O. Bankole   PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
210a515125bSLeila Ghaffari 
211a515125bSLeila Ghaffari   PetscFunctionReturn(0);
212a515125bSLeila Ghaffari }
213a515125bSLeila Ghaffari 
2142b916ea7SJeremy L Thompson static PetscErrorCode FormPreallocation(User user, PetscBool pbdiagonal, Mat J, CeedVector *coo_values) {
215b107fddaSJed Brown   PetscCount ncoo;
216b107fddaSJed Brown   PetscInt  *rows, *cols;
217b107fddaSJed Brown 
218b107fddaSJed Brown   PetscFunctionBeginUser;
219b107fddaSJed Brown   if (pbdiagonal) {
220b107fddaSJed Brown     CeedSize l_size;
221b107fddaSJed Brown     CeedOperatorGetActiveVectorLengths(user->op_ijacobian, &l_size, NULL);
222b107fddaSJed Brown     ncoo = l_size * 5;
223b107fddaSJed Brown     rows = malloc(ncoo * sizeof(rows[0]));
224b107fddaSJed Brown     cols = malloc(ncoo * sizeof(cols[0]));
225b107fddaSJed Brown     for (PetscCount n = 0; n < l_size / 5; n++) {
226b107fddaSJed Brown       for (PetscInt i = 0; i < 5; i++) {
227b107fddaSJed Brown         for (PetscInt j = 0; j < 5; j++) {
228b107fddaSJed Brown           rows[(n * 5 + i) * 5 + j] = n * 5 + i;
229b107fddaSJed Brown           cols[(n * 5 + i) * 5 + j] = n * 5 + j;
230b107fddaSJed Brown         }
231b107fddaSJed Brown       }
232b107fddaSJed Brown     }
233b107fddaSJed Brown   } else {
2342b916ea7SJeremy L Thompson     PetscCall(CeedOperatorLinearAssembleSymbolic(user->op_ijacobian, &ncoo, &rows, &cols));
235b107fddaSJed Brown   }
236b107fddaSJed Brown   PetscCall(MatSetPreallocationCOOLocal(J, ncoo, rows, cols));
237b107fddaSJed Brown   free(rows);
238b107fddaSJed Brown   free(cols);
239b107fddaSJed Brown   CeedVectorCreate(user->ceed, ncoo, coo_values);
240b107fddaSJed Brown   PetscFunctionReturn(0);
241b107fddaSJed Brown }
242b107fddaSJed Brown 
2432b916ea7SJeremy L Thompson static PetscErrorCode FormSetValues(User user, PetscBool pbdiagonal, Mat J, CeedVector coo_values) {
244b107fddaSJed Brown   CeedMemType        mem_type = CEED_MEM_HOST;
245b107fddaSJed Brown   const PetscScalar *values;
246b107fddaSJed Brown   MatType            mat_type;
247b107fddaSJed Brown 
248b107fddaSJed Brown   PetscFunctionBeginUser;
249b107fddaSJed Brown   PetscCall(MatGetType(J, &mat_type));
2502b916ea7SJeremy L Thompson   if (strstr(mat_type, "kokkos") || strstr(mat_type, "cusparse")) mem_type = CEED_MEM_DEVICE;
251*cb315d14SJames Wright   if (pbdiagonal) {
2522b916ea7SJeremy L Thompson     CeedOperatorLinearAssemblePointBlockDiagonal(user->op_ijacobian, coo_values, CEED_REQUEST_IMMEDIATE);
253b107fddaSJed Brown   } else {
254b107fddaSJed Brown     CeedOperatorLinearAssemble(user->op_ijacobian, coo_values);
255b107fddaSJed Brown   }
256b107fddaSJed Brown   CeedVectorGetArrayRead(coo_values, mem_type, &values);
257b107fddaSJed Brown   PetscCall(MatSetValuesCOO(J, values, INSERT_VALUES));
258b107fddaSJed Brown   CeedVectorRestoreArrayRead(coo_values, &values);
259b107fddaSJed Brown   PetscFunctionReturn(0);
260b107fddaSJed Brown }
261b107fddaSJed Brown 
2622b916ea7SJeremy L Thompson PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, PetscReal shift, Mat J, Mat J_pre, void *user_data) {
263f0b65372SJed Brown   User      user = *(User *)user_data;
26404855949SJed Brown   PetscBool J_is_shell, J_is_mffd, J_pre_is_shell;
265f0b65372SJed Brown   PetscFunctionBeginUser;
266a5f46a7bSJeremy L Thompson   if (user->phys->ijacobian_time_shift_label) CeedOperatorSetContextDouble(user->op_ijacobian, user->phys->ijacobian_time_shift_label, &shift);
26704855949SJed Brown   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATMFFD, &J_is_mffd));
268f0b65372SJed Brown   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATSHELL, &J_is_shell));
2692b916ea7SJeremy L Thompson   PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATSHELL, &J_pre_is_shell));
270f0b65372SJed Brown   if (!user->matrices_set_up) {
271f0b65372SJed Brown     if (J_is_shell) {
272f9028c3cSJames Wright       OperatorApplyContext op_ijacobian_ctx;
273f9028c3cSJames Wright       OperatorApplyContextCreate(user->dm, user->dm, user->ceed, user->op_ijacobian, user->q_ceed, user->g_ceed, user->Q_dot_loc, NULL,
274f9028c3cSJames Wright                                  &op_ijacobian_ctx);
275f9028c3cSJames Wright       PetscCall(MatShellSetContext(J, op_ijacobian_ctx));
276f9028c3cSJames Wright       PetscCall(MatShellSetContextDestroy(J, (PetscErrorCode(*)(void *))OperatorApplyContextDestroy));
277f9028c3cSJames Wright       PetscCall(MatShellSetOperation(J, MATOP_MULT, (void (*)(void))MatMult_Ceed));
278f9028c3cSJames Wright       PetscCall(MatShellSetOperation(J, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiag_Ceed));
279f0b65372SJed Brown       PetscCall(MatSetUp(J));
280f0b65372SJed Brown     }
281f0b65372SJed Brown     if (!J_pre_is_shell) {
2822b916ea7SJeremy L Thompson       PetscCall(FormPreallocation(user, user->app_ctx->pmat_pbdiagonal, J_pre, &user->coo_values_pmat));
283b107fddaSJed Brown     }
28404855949SJed Brown     if (J != J_pre && !J_is_shell && !J_is_mffd) {
285b107fddaSJed Brown       PetscCall(FormPreallocation(user, PETSC_FALSE, J, &user->coo_values_amat));
286b107fddaSJed Brown     }
287f0b65372SJed Brown     user->matrices_set_up = true;
288f0b65372SJed Brown   }
289f0b65372SJed Brown   if (!J_pre_is_shell) {
2902b916ea7SJeremy L Thompson     PetscCall(FormSetValues(user, user->app_ctx->pmat_pbdiagonal, J_pre, user->coo_values_pmat));
291f0b65372SJed Brown   }
29204855949SJed Brown   if (user->coo_values_amat) {
29304855949SJed Brown     PetscCall(FormSetValues(user, PETSC_FALSE, J, user->coo_values_amat));
29404855949SJed Brown   } else if (J_is_mffd) {
29504855949SJed Brown     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
29604855949SJed Brown     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
29704855949SJed Brown   }
298f0b65372SJed Brown   PetscFunctionReturn(0);
299f0b65372SJed Brown }
300f0b65372SJed Brown 
3012b916ea7SJeremy L Thompson PetscErrorCode WriteOutput(User user, Vec Q, PetscInt step_no, PetscScalar time) {
302a515125bSLeila Ghaffari   Vec         Q_loc;
303a515125bSLeila Ghaffari   char        file_path[PETSC_MAX_PATH_LEN];
304a515125bSLeila Ghaffari   PetscViewer viewer;
305a515125bSLeila Ghaffari   PetscFunctionBeginUser;
306a515125bSLeila Ghaffari 
307852e5969SJed Brown   if (user->app_ctx->checkpoint_vtk) {
308a515125bSLeila Ghaffari     // Set up output
3097538d537SJames Wright     PetscCall(DMGetLocalVector(user->dm, &Q_loc));
3107538d537SJames Wright     PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec"));
3117538d537SJames Wright     PetscCall(VecZeroEntries(Q_loc));
3127538d537SJames Wright     PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
313a515125bSLeila Ghaffari 
314a515125bSLeila Ghaffari     // Output
315852e5969SJed Brown     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
3167538d537SJames Wright 
3172b916ea7SJeremy L Thompson     PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path, FILE_MODE_WRITE, &viewer));
3187538d537SJames Wright     PetscCall(VecView(Q_loc, viewer));
3197538d537SJames Wright     PetscCall(PetscViewerDestroy(&viewer));
320a515125bSLeila Ghaffari     if (user->dm_viz) {
321a515125bSLeila Ghaffari       Vec         Q_refined, Q_refined_loc;
322a515125bSLeila Ghaffari       char        file_path_refined[PETSC_MAX_PATH_LEN];
323a515125bSLeila Ghaffari       PetscViewer viewer_refined;
324a515125bSLeila Ghaffari 
3257538d537SJames Wright       PetscCall(DMGetGlobalVector(user->dm_viz, &Q_refined));
3267538d537SJames Wright       PetscCall(DMGetLocalVector(user->dm_viz, &Q_refined_loc));
3277538d537SJames Wright       PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined"));
3287538d537SJames Wright 
3297538d537SJames Wright       PetscCall(MatInterpolate(user->interp_viz, Q, Q_refined));
3307538d537SJames Wright       PetscCall(VecZeroEntries(Q_refined_loc));
3312b916ea7SJeremy L Thompson       PetscCall(DMGlobalToLocal(user->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc));
3327538d537SJames Wright 
333852e5969SJed Brown       PetscCall(
334852e5969SJed Brown           PetscSNPrintf(file_path_refined, sizeof file_path_refined, "%s/nsrefined-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
3357538d537SJames Wright 
3362b916ea7SJeremy L Thompson       PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined), file_path_refined, FILE_MODE_WRITE, &viewer_refined));
3377538d537SJames Wright       PetscCall(VecView(Q_refined_loc, viewer_refined));
3387538d537SJames Wright       PetscCall(DMRestoreLocalVector(user->dm_viz, &Q_refined_loc));
3397538d537SJames Wright       PetscCall(DMRestoreGlobalVector(user->dm_viz, &Q_refined));
3407538d537SJames Wright       PetscCall(PetscViewerDestroy(&viewer_refined));
341a515125bSLeila Ghaffari     }
3427538d537SJames Wright     PetscCall(DMRestoreLocalVector(user->dm, &Q_loc));
343852e5969SJed Brown   }
344a515125bSLeila Ghaffari 
345a515125bSLeila Ghaffari   // Save data in a binary file for continuation of simulations
34691a36801SJames Wright   if (user->app_ctx->add_stepnum2bin) {
347852e5969SJed Brown     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", user->app_ctx->output_dir, step_no));
34891a36801SJames Wright   } else {
3492b916ea7SJeremy L Thompson     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", user->app_ctx->output_dir));
35091a36801SJames Wright   }
3512b916ea7SJeremy L Thompson   PetscCall(PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, &viewer));
3527538d537SJames Wright 
3539293eaa1SJed Brown   PetscInt token = FLUIDS_FILE_TOKEN;
3549293eaa1SJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &token, 1, PETSC_INT));
3559293eaa1SJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &step_no, 1, PETSC_INT));
3569293eaa1SJed Brown   time /= user->units->second;  // Dimensionalize time back
3579293eaa1SJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL));
3587538d537SJames Wright   PetscCall(VecView(Q, viewer));
3597538d537SJames Wright   PetscCall(PetscViewerDestroy(&viewer));
3607538d537SJames Wright   PetscFunctionReturn(0);
3617538d537SJames Wright }
3627538d537SJames Wright 
363c5e9980aSAdeleke O. Bankole // CSV Monitor
364c5e9980aSAdeleke O. Bankole PetscErrorCode TSMonitor_WallForce(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
365c5e9980aSAdeleke O. Bankole   User              user = ctx;
366c5e9980aSAdeleke O. Bankole   Vec               G_loc;
367c5e9980aSAdeleke O. Bankole   PetscInt          num_wall = user->app_ctx->wall_forces.num_wall, dim = 3;
368c5e9980aSAdeleke O. Bankole   const PetscInt   *walls  = user->app_ctx->wall_forces.walls;
369c5e9980aSAdeleke O. Bankole   PetscViewer       viewer = user->app_ctx->wall_forces.viewer;
370c5e9980aSAdeleke O. Bankole   PetscViewerFormat format = user->app_ctx->wall_forces.viewer_format;
371c5e9980aSAdeleke O. Bankole   PetscScalar      *reaction_force;
372c5e9980aSAdeleke O. Bankole   PetscBool         iascii;
373c5e9980aSAdeleke O. Bankole 
374c5e9980aSAdeleke O. Bankole   PetscFunctionBeginUser;
375c5e9980aSAdeleke O. Bankole   if (!viewer) PetscFunctionReturn(0);
376c5e9980aSAdeleke O. Bankole   PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
377c5e9980aSAdeleke O. Bankole   PetscCall(PetscMalloc1(num_wall * dim, &reaction_force));
378c5e9980aSAdeleke O. Bankole   PetscCall(Surface_Forces_NS(user->dm, G_loc, num_wall, walls, reaction_force));
379c5e9980aSAdeleke O. Bankole   PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
380c5e9980aSAdeleke O. Bankole 
381c5e9980aSAdeleke O. Bankole   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
382c5e9980aSAdeleke O. Bankole 
383c5e9980aSAdeleke O. Bankole   if (iascii) {
384c5e9980aSAdeleke O. Bankole     if (format == PETSC_VIEWER_ASCII_CSV && !user->app_ctx->wall_forces.header_written) {
385c5e9980aSAdeleke O. Bankole       PetscCall(PetscViewerASCIIPrintf(viewer, "Step,Time,Wall,ForceX,ForceY,ForceZ\n"));
386c5e9980aSAdeleke O. Bankole       user->app_ctx->wall_forces.header_written = PETSC_TRUE;
387c5e9980aSAdeleke O. Bankole     }
388c5e9980aSAdeleke O. Bankole     for (PetscInt w = 0; w < num_wall; w++) {
389c5e9980aSAdeleke O. Bankole       PetscInt wall = walls[w];
390c5e9980aSAdeleke O. Bankole       if (format == PETSC_VIEWER_ASCII_CSV) {
391c5e9980aSAdeleke O. Bankole         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT ",%g,%" PetscInt_FMT ",%g,%g,%g\n", step_no, time, wall,
392c5e9980aSAdeleke O. Bankole                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
393c5e9980aSAdeleke O. Bankole 
394c5e9980aSAdeleke O. Bankole       } else {
395c5e9980aSAdeleke O. Bankole         PetscCall(PetscViewerASCIIPrintf(viewer, "Wall %" PetscInt_FMT " Forces: Force_x = %12g, Force_y = %12g, Force_z = %12g\n", wall,
396c5e9980aSAdeleke O. Bankole                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
397c5e9980aSAdeleke O. Bankole       }
398c5e9980aSAdeleke O. Bankole     }
399c5e9980aSAdeleke O. Bankole   }
400c5e9980aSAdeleke O. Bankole   PetscCall(PetscFree(reaction_force));
401c5e9980aSAdeleke O. Bankole   PetscFunctionReturn(0);
402c5e9980aSAdeleke O. Bankole }
403c5e9980aSAdeleke O. Bankole 
4047538d537SJames Wright // User provided TS Monitor
4052b916ea7SJeremy L Thompson PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
4067538d537SJames Wright   User user = ctx;
4077538d537SJames Wright   PetscFunctionBeginUser;
4087538d537SJames Wright 
409852e5969SJed Brown   // Print every 'checkpoint_interval' steps
410c539088bSJames Wright   if (user->app_ctx->checkpoint_interval <= 0 || step_no % user->app_ctx->checkpoint_interval != 0 ||
411e419654dSJeremy L Thompson       (user->app_ctx->cont_steps == step_no && step_no != 0)) {
412c539088bSJames Wright     PetscFunctionReturn(0);
413e419654dSJeremy L Thompson   }
4147538d537SJames Wright 
4157538d537SJames Wright   PetscCall(WriteOutput(user, Q, step_no, time));
416a515125bSLeila Ghaffari 
417a515125bSLeila Ghaffari   PetscFunctionReturn(0);
418a515125bSLeila Ghaffari }
419a515125bSLeila Ghaffari 
420a515125bSLeila Ghaffari // TS: Create, setup, and solve
4212b916ea7SJeremy L Thompson PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, Vec *Q, PetscScalar *f_time, TS *ts) {
422a515125bSLeila Ghaffari   MPI_Comm    comm = user->comm;
423a515125bSLeila Ghaffari   TSAdapt     adapt;
424a515125bSLeila Ghaffari   PetscScalar final_time;
425a515125bSLeila Ghaffari   PetscFunctionBeginUser;
426a515125bSLeila Ghaffari 
4272b916ea7SJeremy L Thompson   PetscCall(TSCreate(comm, ts));
4282b916ea7SJeremy L Thompson   PetscCall(TSSetDM(*ts, dm));
429a515125bSLeila Ghaffari   if (phys->implicit) {
4302b916ea7SJeremy L Thompson     PetscCall(TSSetType(*ts, TSBDF));
431a515125bSLeila Ghaffari     if (user->op_ifunction) {
4322b916ea7SJeremy L Thompson       PetscCall(TSSetIFunction(*ts, NULL, IFunction_NS, &user));
433a515125bSLeila Ghaffari     } else {  // Implicit integrators can fall back to using an RHSFunction
4342b916ea7SJeremy L Thompson       PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
435a515125bSLeila Ghaffari     }
436f0b65372SJed Brown     if (user->op_ijacobian) {
4372b916ea7SJeremy L Thompson       PetscCall(DMTSSetIJacobian(dm, FormIJacobian_NS, &user));
438b107fddaSJed Brown       if (app_ctx->amat_type) {
439b107fddaSJed Brown         Mat Pmat, Amat;
4402b916ea7SJeremy L Thompson         PetscCall(DMCreateMatrix(dm, &Pmat));
4412b916ea7SJeremy L Thompson         PetscCall(DMSetMatType(dm, app_ctx->amat_type));
4422b916ea7SJeremy L Thompson         PetscCall(DMCreateMatrix(dm, &Amat));
4432b916ea7SJeremy L Thompson         PetscCall(TSSetIJacobian(*ts, Amat, Pmat, NULL, NULL));
4442b916ea7SJeremy L Thompson         PetscCall(MatDestroy(&Amat));
4452b916ea7SJeremy L Thompson         PetscCall(MatDestroy(&Pmat));
446b107fddaSJed Brown       }
447f0b65372SJed Brown     }
448a515125bSLeila Ghaffari   } else {
449da5fe0e4SJames Wright     PetscCheck(user->op_rhs_ctx, comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction");
4502b916ea7SJeremy L Thompson     PetscCall(TSSetType(*ts, TSRK));
4512b916ea7SJeremy L Thompson     PetscCall(TSRKSetType(*ts, TSRK5F));
4522b916ea7SJeremy L Thompson     PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
453a515125bSLeila Ghaffari   }
4542b916ea7SJeremy L Thompson   PetscCall(TSSetMaxTime(*ts, 500. * user->units->second));
4552b916ea7SJeremy L Thompson   PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER));
45622387d3aSJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(TSSetErrorIfStepFails(*ts, PETSC_FALSE));
4572b916ea7SJeremy L Thompson   PetscCall(TSSetTimeStep(*ts, 1.e-2 * user->units->second));
4580e1e9333SJames Wright   if (app_ctx->test_type != TESTTYPE_NONE) {
4592b916ea7SJeremy L Thompson     PetscCall(TSSetMaxSteps(*ts, 10));
4602b916ea7SJeremy L Thompson   }
4612b916ea7SJeremy L Thompson   PetscCall(TSGetAdapt(*ts, &adapt));
4622b916ea7SJeremy L Thompson   PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * user->units->second, 1.e2 * user->units->second));
4632b916ea7SJeremy L Thompson   PetscCall(TSSetFromOptions(*ts));
46491f639d2SJames Wright   user->time_bc_set = -1.0;    // require all BCs be updated
465a515125bSLeila Ghaffari   if (!app_ctx->cont_steps) {  // print initial condition
4660e1e9333SJames Wright     if (app_ctx->test_type == TESTTYPE_NONE) {
4672b916ea7SJeremy L Thompson       PetscCall(TSMonitor_NS(*ts, 0, 0., *Q, user));
468a515125bSLeila Ghaffari     }
469a515125bSLeila Ghaffari   } else {  // continue from time of last output
470a515125bSLeila Ghaffari     PetscInt    count;
471a515125bSLeila Ghaffari     PetscViewer viewer;
4722b916ea7SJeremy L Thompson 
4739293eaa1SJed Brown     if (app_ctx->cont_time <= 0) {  // Legacy files did not include step number and time
4742b916ea7SJeremy L Thompson       PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_time_file, FILE_MODE_READ, &viewer));
4759293eaa1SJed Brown       PetscCall(PetscViewerBinaryRead(viewer, &app_ctx->cont_time, 1, &count, PETSC_REAL));
4762b916ea7SJeremy L Thompson       PetscCall(PetscViewerDestroy(&viewer));
4779293eaa1SJed Brown       PetscCheck(app_ctx->cont_steps != -1, comm, PETSC_ERR_ARG_INCOMP,
4789293eaa1SJed Brown                  "-continue step number not specified, but checkpoint file does not contain a step number (likely written by older code version)");
4799293eaa1SJed Brown     }
4809293eaa1SJed Brown     PetscCall(TSSetTime(*ts, app_ctx->cont_time * user->units->second));
48174a6f4ddSJed Brown     PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps));
482a515125bSLeila Ghaffari   }
4830e1e9333SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) {
4842b916ea7SJeremy L Thompson     PetscCall(TSMonitorSet(*ts, TSMonitor_NS, user, NULL));
4850e1e9333SJames Wright   }
486c5e9980aSAdeleke O. Bankole   if (app_ctx->wall_forces.viewer) {
487c5e9980aSAdeleke O. Bankole     PetscCall(TSMonitorSet(*ts, TSMonitor_WallForce, user, NULL));
488c5e9980aSAdeleke O. Bankole   }
489c931fa59SJames Wright   if (app_ctx->turb_spanstats_enable) {
49091933550SJames Wright     PetscCall(TSMonitorSet(*ts, TSMonitor_TurbulenceStatistics, user, NULL));
491b8daee98SJames Wright     CeedScalar previous_time = app_ctx->cont_time * user->units->second;
4925cfe26e6SJames Wright     CeedOperatorSetContextDouble(user->spanstats.op_stats_collect_ctx->op, user->spanstats.previous_time_label, &previous_time);
493b0488d1fSJames Wright   }
494a515125bSLeila Ghaffari 
495a515125bSLeila Ghaffari   // Solve
49674a6f4ddSJed Brown   PetscReal start_time;
49774a6f4ddSJed Brown   PetscInt  start_step;
4982b916ea7SJeremy L Thompson   PetscCall(TSGetTime(*ts, &start_time));
49974a6f4ddSJed Brown   PetscCall(TSGetStepNumber(*ts, &start_step));
50091982731SJeremy L Thompson 
50191982731SJeremy L Thompson   PetscPreLoadBegin(PETSC_FALSE, "Fluids Solve");
50291982731SJeremy L Thompson   PetscCall(TSSetTime(*ts, start_time));
50374a6f4ddSJed Brown   PetscCall(TSSetStepNumber(*ts, start_step));
50491982731SJeremy L Thompson   if (PetscPreLoadingOn) {
50591982731SJeremy L Thompson     // LCOV_EXCL_START
50691982731SJeremy L Thompson     SNES      snes;
50791982731SJeremy L Thompson     Vec       Q_preload;
50891982731SJeremy L Thompson     PetscReal rtol;
50991982731SJeremy L Thompson     PetscCall(VecDuplicate(*Q, &Q_preload));
51091982731SJeremy L Thompson     PetscCall(VecCopy(*Q, Q_preload));
51191982731SJeremy L Thompson     PetscCall(TSGetSNES(*ts, &snes));
51291982731SJeremy L Thompson     PetscCall(SNESGetTolerances(snes, NULL, &rtol, NULL, NULL, NULL));
5132b916ea7SJeremy L Thompson     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, .99, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
51491982731SJeremy L Thompson     PetscCall(TSSetSolution(*ts, *Q));
51591982731SJeremy L Thompson     PetscCall(TSStep(*ts));
5162b916ea7SJeremy L Thompson     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, rtol, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
51791982731SJeremy L Thompson     PetscCall(VecDestroy(&Q_preload));
51891982731SJeremy L Thompson     // LCOV_EXCL_STOP
51991982731SJeremy L Thompson   } else {
5202b916ea7SJeremy L Thompson     PetscCall(PetscBarrier((PetscObject)*ts));
5212b916ea7SJeremy L Thompson     PetscCall(TSSolve(*ts, *Q));
52291982731SJeremy L Thompson   }
52391982731SJeremy L Thompson   PetscPreLoadEnd();
52491982731SJeremy L Thompson 
52591982731SJeremy L Thompson   PetscCall(TSGetSolveTime(*ts, &final_time));
526a515125bSLeila Ghaffari   *f_time = final_time;
52791982731SJeremy L Thompson 
5280e1e9333SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) {
5297538d537SJames Wright     PetscInt step_no;
5307538d537SJames Wright     PetscCall(TSGetStepNumber(*ts, &step_no));
531b0488d1fSJames Wright     if (user->app_ctx->checkpoint_interval > 0 || user->app_ctx->checkpoint_interval == -1) {
5327538d537SJames Wright       PetscCall(WriteOutput(user, *Q, step_no, final_time));
5337538d537SJames Wright     }
5347538d537SJames Wright 
53591982731SJeremy L Thompson     PetscLogEvent stage_id;
53691982731SJeremy L Thompson     PetscStageLog stage_log;
53791982731SJeremy L Thompson 
53891982731SJeremy L Thompson     PetscCall(PetscLogStageGetId("Fluids Solve", &stage_id));
53991982731SJeremy L Thompson     PetscCall(PetscLogGetStageLog(&stage_log));
5402b916ea7SJeremy L Thompson     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time taken for solution (sec): %g\n", stage_log->stageInfo[stage_id].perfInfo.time));
541a515125bSLeila Ghaffari   }
542a515125bSLeila Ghaffari   PetscFunctionReturn(0);
543a515125bSLeila Ghaffari }
544