xref: /honee/src/setupts.c (revision ebfabadf5ae225a12a8c93971166571a57c343d2)
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 
18f8e2d240SJames Wright // @brief Create KSP to solve the inverse mass operator for explicit time stepping schemes
19d6ca8aaaSJames Wright PetscErrorCode CreateKSPMassOperator(User user, CeedData ceed_data) {
20f8e2d240SJames Wright   Ceed                 ceed = user->ceed;
21f8e2d240SJames Wright   DM                   dm   = user->dm;
22a515125bSLeila Ghaffari   CeedQFunction        qf_mass;
23a515125bSLeila Ghaffari   CeedOperator         op_mass;
24f8e2d240SJames Wright   OperatorApplyContext mass_matop_ctx;
25a515125bSLeila Ghaffari   CeedInt              num_comp_q, q_data_size;
26a515125bSLeila Ghaffari 
2706f41313SJames Wright   PetscFunctionBeginUser;
28b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q));
29b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &q_data_size));
30a515125bSLeila Ghaffari 
319f59f36eSJames Wright   PetscCall(CreateMassQFunction(ceed, num_comp_q, q_data_size, &qf_mass));
32b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass));
33b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
3458e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
35b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
36a515125bSLeila Ghaffari 
37f8e2d240SJames Wright   {  // -- Setup KSP for mass operator
3845101827SJames Wright     Mat      mat_mass = NULL;
39614823f4SJames Wright     Vec      Zeros_loc;
40f8e2d240SJames Wright     MPI_Comm comm = PetscObjectComm((PetscObject)dm);
41a515125bSLeila Ghaffari 
42614823f4SJames Wright     PetscCall(DMCreateLocalVector(dm, &Zeros_loc));
43614823f4SJames Wright     PetscCall(VecZeroEntries(Zeros_loc));
44614823f4SJames Wright     PetscCall(OperatorApplyContextCreate(dm, dm, ceed, op_mass, NULL, NULL, Zeros_loc, NULL, &mass_matop_ctx));
45f8e2d240SJames Wright     PetscCall(CreateMatShell_Ceed(mass_matop_ctx, &mat_mass));
46a515125bSLeila Ghaffari 
47f8e2d240SJames Wright     PetscCall(KSPCreate(comm, &user->mass_ksp));
48f8e2d240SJames Wright     PetscCall(KSPSetOptionsPrefix(user->mass_ksp, "mass_"));
49f8e2d240SJames Wright     {  // lumped by default
50f8e2d240SJames Wright       PC pc;
51f8e2d240SJames Wright       PetscCall(KSPGetPC(user->mass_ksp, &pc));
52f8e2d240SJames Wright       PetscCall(PCSetType(pc, PCJACOBI));
53f8e2d240SJames Wright       PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
54f8e2d240SJames Wright       PetscCall(KSPSetType(user->mass_ksp, KSPPREONLY));
55f8e2d240SJames Wright     }
56f8e2d240SJames Wright     PetscCall(KSPSetOperators(user->mass_ksp, mat_mass, mat_mass));
57f8e2d240SJames Wright     PetscCall(KSPSetFromOptions(user->mass_ksp));
58614823f4SJames Wright     PetscCall(VecDestroy(&Zeros_loc));
59614823f4SJames Wright     PetscCall(MatDestroy(&mat_mass));
60f8e2d240SJames Wright   }
61a515125bSLeila Ghaffari 
62a515125bSLeila Ghaffari   // Cleanup
63b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
64b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
65d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
66a515125bSLeila Ghaffari }
67a515125bSLeila Ghaffari 
68c996854bSJames Wright // Insert Boundary values if it's a new time
69c996854bSJames Wright PetscErrorCode UpdateBoundaryValues(User user, Vec Q_loc, PetscReal t) {
70c996854bSJames Wright   PetscFunctionBeginUser;
71c996854bSJames Wright   if (user->time_bc_set != t) {
72c996854bSJames Wright     PetscCall(DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL));
73c996854bSJames Wright     user->time_bc_set = t;
74c996854bSJames Wright   }
75d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
76c996854bSJames Wright }
77c996854bSJames Wright 
78a515125bSLeila Ghaffari // RHS (Explicit time-stepper) function setup
79a515125bSLeila Ghaffari //   This is the RHS of the ODE, given as u_t = G(t,u)
80a515125bSLeila Ghaffari //   This function takes in a state vector Q and writes into G
81a515125bSLeila Ghaffari PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) {
82a515125bSLeila Ghaffari   User        user = *(User *)user_data;
83701e5830SJames Wright   Ceed        ceed = user->ceed;
84fd969b44SJames Wright   PetscScalar dt;
85da5fe0e4SJames Wright   Vec         Q_loc = user->Q_loc;
86a515125bSLeila Ghaffari 
8706f41313SJames Wright   PetscFunctionBeginUser;
88e2f84137SJeremy L Thompson   // Update time dependent data
89c996854bSJames Wright   PetscCall(UpdateBoundaryValues(user, Q_loc, t));
90701e5830SJames Wright   if (user->phys->solution_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_rhs_ctx->op, user->phys->solution_time_label, &t));
912b916ea7SJeremy L Thompson   PetscCall(TSGetTimeStep(ts, &dt));
92701e5830SJames Wright   if (user->phys->timestep_size_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_rhs_ctx->op, user->phys->timestep_size_label, &dt));
93a515125bSLeila Ghaffari 
94da5fe0e4SJames Wright   PetscCall(ApplyCeedOperatorGlobalToGlobal(Q, G, user->op_rhs_ctx));
95a515125bSLeila Ghaffari 
96f3fcf8f4SJames Wright   // Inverse of the lumped mass matrix
97f8e2d240SJames Wright   PetscCall(KSPSolve(user->mass_ksp, G, G));
98d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
99a515125bSLeila Ghaffari }
100a515125bSLeila Ghaffari 
101c5e9980aSAdeleke O. Bankole // Surface forces function setup
102c5e9980aSAdeleke O. Bankole static PetscErrorCode Surface_Forces_NS(DM dm, Vec G_loc, PetscInt num_walls, const PetscInt walls[], PetscScalar *reaction_force) {
103c5e9980aSAdeleke O. Bankole   DMLabel            face_label;
104c5e9980aSAdeleke O. Bankole   const PetscScalar *g;
1052004e3acSAdeleke O. Bankole   PetscInt           dof, dim = 3;
106c5e9980aSAdeleke O. Bankole   MPI_Comm           comm;
1072004e3acSAdeleke O. Bankole   PetscSection       s;
108c5e9980aSAdeleke O. Bankole 
109c5e9980aSAdeleke O. Bankole   PetscFunctionBeginUser;
110c5e9980aSAdeleke O. Bankole   PetscCall(PetscArrayzero(reaction_force, num_walls * dim));
111c5e9980aSAdeleke O. Bankole   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
112c5e9980aSAdeleke O. Bankole   PetscCall(DMGetLabel(dm, "Face Sets", &face_label));
113c5e9980aSAdeleke O. Bankole   PetscCall(VecGetArrayRead(G_loc, &g));
114c5e9980aSAdeleke O. Bankole   for (PetscInt w = 0; w < num_walls; w++) {
115c5e9980aSAdeleke O. Bankole     const PetscInt wall = walls[w];
116c5e9980aSAdeleke O. Bankole     IS             wall_is;
1172004e3acSAdeleke O. Bankole     PetscCall(DMGetLocalSection(dm, &s));
118c5e9980aSAdeleke O. Bankole     PetscCall(DMLabelGetStratumIS(face_label, wall, &wall_is));
119c5e9980aSAdeleke O. Bankole     if (wall_is) {  // There exist such points on this process
120c5e9980aSAdeleke O. Bankole       PetscInt        num_points;
1212004e3acSAdeleke O. Bankole       PetscInt        num_comp = 0;
122c5e9980aSAdeleke O. Bankole       const PetscInt *points;
1232004e3acSAdeleke O. Bankole       PetscCall(PetscSectionGetFieldComponents(s, 0, &num_comp));
124c5e9980aSAdeleke O. Bankole       PetscCall(ISGetSize(wall_is, &num_points));
125c5e9980aSAdeleke O. Bankole       PetscCall(ISGetIndices(wall_is, &points));
126c5e9980aSAdeleke O. Bankole       for (PetscInt i = 0; i < num_points; i++) {
127c5e9980aSAdeleke O. Bankole         const PetscInt           p = points[i];
128c5e9980aSAdeleke O. Bankole         const StateConservative *r;
129c5e9980aSAdeleke O. Bankole         PetscCall(DMPlexPointLocalRead(dm, p, g, &r));
1302004e3acSAdeleke O. Bankole         PetscCall(PetscSectionGetDof(s, p, &dof));
1312004e3acSAdeleke O. Bankole         for (PetscInt node = 0; node < dof / num_comp; node++) {
132c5e9980aSAdeleke O. Bankole           for (PetscInt j = 0; j < 3; j++) {
1332004e3acSAdeleke O. Bankole             reaction_force[w * dim + j] -= r[node].momentum[j];
1342004e3acSAdeleke O. Bankole           }
135c5e9980aSAdeleke O. Bankole         }
136c5e9980aSAdeleke O. Bankole       }
137c5e9980aSAdeleke O. Bankole       PetscCall(ISRestoreIndices(wall_is, &points));
138c5e9980aSAdeleke O. Bankole     }
139c5e9980aSAdeleke O. Bankole     PetscCall(ISDestroy(&wall_is));
140c5e9980aSAdeleke O. Bankole   }
141c5e9980aSAdeleke O. Bankole   PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, reaction_force, dim * num_walls, MPIU_SCALAR, MPI_SUM, comm));
142c5e9980aSAdeleke O. Bankole   //  Restore Vectors
143c5e9980aSAdeleke O. Bankole   PetscCall(VecRestoreArrayRead(G_loc, &g));
144d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
145c5e9980aSAdeleke O. Bankole }
146c5e9980aSAdeleke O. Bankole 
147a515125bSLeila Ghaffari // Implicit time-stepper function setup
1482b916ea7SJeremy L Thompson PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) {
149a515125bSLeila Ghaffari   User         user = *(User *)user_data;
150701e5830SJames Wright   Ceed         ceed = user->ceed;
151fd969b44SJames Wright   PetscScalar  dt;
152e2f84137SJeremy L Thompson   Vec          Q_loc = user->Q_loc, Q_dot_loc = user->Q_dot_loc, G_loc;
153a515125bSLeila Ghaffari   PetscMemType q_mem_type, q_dot_mem_type, g_mem_type;
154a515125bSLeila Ghaffari 
15506f41313SJames Wright   PetscFunctionBeginUser;
156e2f84137SJeremy L Thompson   // Get local vectors
157c5e9980aSAdeleke O. Bankole   PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
158e2f84137SJeremy L Thompson 
159e2f84137SJeremy L Thompson   // Update time dependent data
160c996854bSJames Wright   PetscCall(UpdateBoundaryValues(user, Q_loc, t));
161701e5830SJames Wright   if (user->phys->solution_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_ifunction, user->phys->solution_time_label, &t));
1622b916ea7SJeremy L Thompson   PetscCall(TSGetTimeStep(ts, &dt));
163701e5830SJames Wright   if (user->phys->timestep_size_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_ifunction, user->phys->timestep_size_label, &dt));
164a515125bSLeila Ghaffari 
165a515125bSLeila Ghaffari   // Global-to-local
16606108310SJames Wright   PetscCall(DMGlobalToLocalBegin(user->dm, Q, INSERT_VALUES, Q_loc));
16706108310SJames Wright   PetscCall(DMGlobalToLocalBegin(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
16806108310SJames Wright   PetscCall(DMGlobalToLocalEnd(user->dm, Q, INSERT_VALUES, Q_loc));
16906108310SJames Wright   PetscCall(DMGlobalToLocalEnd(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
170a515125bSLeila Ghaffari 
171a515125bSLeila Ghaffari   // Place PETSc vectors in CEED vectors
172fd969b44SJames Wright   PetscCall(VecReadP2C(Q_loc, &q_mem_type, user->q_ceed));
173fd969b44SJames Wright   PetscCall(VecReadP2C(Q_dot_loc, &q_dot_mem_type, user->q_dot_ceed));
174fd969b44SJames Wright   PetscCall(VecP2C(G_loc, &g_mem_type, user->g_ceed));
175a515125bSLeila Ghaffari 
176a515125bSLeila Ghaffari   // Apply CEED operator
1777eedc94cSJames Wright   PetscCall(PetscLogEventBegin(FLUIDS_CeedOperatorApply, Q, G, 0, 0));
1787eedc94cSJames Wright   PetscCall(PetscLogGpuTimeBegin());
179b4c37c5cSJames Wright   PetscCallCeed(user->ceed, CeedOperatorApply(user->op_ifunction, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE));
1807eedc94cSJames Wright   PetscCall(PetscLogGpuTimeEnd());
1817eedc94cSJames Wright   PetscCall(PetscLogEventEnd(FLUIDS_CeedOperatorApply, Q, G, 0, 0));
182a515125bSLeila Ghaffari 
183a515125bSLeila Ghaffari   // Restore vectors
184fd969b44SJames Wright   PetscCall(VecReadC2P(user->q_ceed, q_mem_type, Q_loc));
185fd969b44SJames Wright   PetscCall(VecReadC2P(user->q_dot_ceed, q_dot_mem_type, Q_dot_loc));
186fd969b44SJames Wright   PetscCall(VecC2P(user->g_ceed, g_mem_type, G_loc));
187a515125bSLeila Ghaffari 
18801ab89c1SJames Wright   if (user->app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) {
189ad494f68SJames Wright     PetscCall(SgsDDApplyIFunction(user, Q_loc, G_loc));
19001ab89c1SJames Wright   }
1919c678832SJames Wright 
192a515125bSLeila Ghaffari   // Local-to-Global
1932b916ea7SJeremy L Thompson   PetscCall(VecZeroEntries(G));
1942b916ea7SJeremy L Thompson   PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G));
195a515125bSLeila Ghaffari 
196a515125bSLeila Ghaffari   // Restore vectors
197c5e9980aSAdeleke O. Bankole   PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
198d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
199a515125bSLeila Ghaffari }
200a515125bSLeila Ghaffari 
2012b916ea7SJeremy L Thompson PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, PetscReal shift, Mat J, Mat J_pre, void *user_data) {
202f0b65372SJed Brown   User      user = *(User *)user_data;
203b4c37c5cSJames Wright   Ceed      ceed = user->ceed;
204*ebfabadfSJames Wright   PetscBool J_is_matceed, J_is_mffd, J_pre_is_matceed, J_pre_is_mffd;
20506f41313SJames Wright 
206f0b65372SJed Brown   PetscFunctionBeginUser;
20704855949SJed Brown   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATMFFD, &J_is_mffd));
208*ebfabadfSJames Wright   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATCEED, &J_is_matceed));
209*ebfabadfSJames Wright   PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATMFFD, &J_pre_is_mffd));
210*ebfabadfSJames Wright   PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATCEED, &J_pre_is_matceed));
211*ebfabadfSJames Wright   if (user->phys->ijacobian_time_shift_label) {
212*ebfabadfSJames Wright     CeedOperator op_ijacobian;
213*ebfabadfSJames Wright 
214*ebfabadfSJames Wright     PetscCall(MatCeedGetCeedOperators(user->mat_ijacobian, &op_ijacobian, NULL));
215*ebfabadfSJames Wright     PetscCallCeed(ceed, CeedOperatorSetContextDouble(op_ijacobian, user->phys->ijacobian_time_shift_label, &shift));
216f0b65372SJed Brown   }
217*ebfabadfSJames Wright 
218*ebfabadfSJames Wright   if (J_is_matceed || J_is_mffd) {
21904855949SJed Brown     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
22004855949SJed Brown     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
221*ebfabadfSJames Wright   } else PetscCall(MatCeedAssembleCOO(user->mat_ijacobian, J));
222*ebfabadfSJames Wright 
223*ebfabadfSJames Wright   if (J_pre_is_matceed && J != J_pre) {
224*ebfabadfSJames Wright     PetscCall(MatAssemblyBegin(J_pre, MAT_FINAL_ASSEMBLY));
225*ebfabadfSJames Wright     PetscCall(MatAssemblyEnd(J_pre, MAT_FINAL_ASSEMBLY));
226*ebfabadfSJames Wright   } else if (!J_pre_is_matceed && !J_pre_is_mffd && J != J_pre) {
227*ebfabadfSJames Wright     PetscCall(MatCeedAssembleCOO(user->mat_ijacobian, J_pre));
22804855949SJed Brown   }
229d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
230f0b65372SJed Brown }
231f0b65372SJed Brown 
2322b916ea7SJeremy L Thompson PetscErrorCode WriteOutput(User user, Vec Q, PetscInt step_no, PetscScalar time) {
233a515125bSLeila Ghaffari   Vec         Q_loc;
234a515125bSLeila Ghaffari   char        file_path[PETSC_MAX_PATH_LEN];
235a515125bSLeila Ghaffari   PetscViewer viewer;
236a515125bSLeila Ghaffari 
23706f41313SJames Wright   PetscFunctionBeginUser;
238852e5969SJed Brown   if (user->app_ctx->checkpoint_vtk) {
239a515125bSLeila Ghaffari     // Set up output
2407538d537SJames Wright     PetscCall(DMGetLocalVector(user->dm, &Q_loc));
2417538d537SJames Wright     PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec"));
2427538d537SJames Wright     PetscCall(VecZeroEntries(Q_loc));
2437538d537SJames Wright     PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
244a515125bSLeila Ghaffari 
245a515125bSLeila Ghaffari     // Output
246852e5969SJed Brown     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
2477538d537SJames Wright 
2482b916ea7SJeremy L Thompson     PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path, FILE_MODE_WRITE, &viewer));
2497538d537SJames Wright     PetscCall(VecView(Q_loc, viewer));
2507538d537SJames Wright     PetscCall(PetscViewerDestroy(&viewer));
251a515125bSLeila Ghaffari     if (user->dm_viz) {
252a515125bSLeila Ghaffari       Vec         Q_refined, Q_refined_loc;
253a515125bSLeila Ghaffari       char        file_path_refined[PETSC_MAX_PATH_LEN];
254a515125bSLeila Ghaffari       PetscViewer viewer_refined;
255a515125bSLeila Ghaffari 
2567538d537SJames Wright       PetscCall(DMGetGlobalVector(user->dm_viz, &Q_refined));
2577538d537SJames Wright       PetscCall(DMGetLocalVector(user->dm_viz, &Q_refined_loc));
2587538d537SJames Wright       PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined"));
2597538d537SJames Wright 
2607538d537SJames Wright       PetscCall(MatInterpolate(user->interp_viz, Q, Q_refined));
2617538d537SJames Wright       PetscCall(VecZeroEntries(Q_refined_loc));
2622b916ea7SJeremy L Thompson       PetscCall(DMGlobalToLocal(user->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc));
2637538d537SJames Wright 
264852e5969SJed Brown       PetscCall(
265852e5969SJed Brown           PetscSNPrintf(file_path_refined, sizeof file_path_refined, "%s/nsrefined-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
2667538d537SJames Wright 
2672b916ea7SJeremy L Thompson       PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined), file_path_refined, FILE_MODE_WRITE, &viewer_refined));
2687538d537SJames Wright       PetscCall(VecView(Q_refined_loc, viewer_refined));
2697538d537SJames Wright       PetscCall(DMRestoreLocalVector(user->dm_viz, &Q_refined_loc));
2707538d537SJames Wright       PetscCall(DMRestoreGlobalVector(user->dm_viz, &Q_refined));
2717538d537SJames Wright       PetscCall(PetscViewerDestroy(&viewer_refined));
272a515125bSLeila Ghaffari     }
2737538d537SJames Wright     PetscCall(DMRestoreLocalVector(user->dm, &Q_loc));
274852e5969SJed Brown   }
275a515125bSLeila Ghaffari 
276a515125bSLeila Ghaffari   // Save data in a binary file for continuation of simulations
27791a36801SJames Wright   if (user->app_ctx->add_stepnum2bin) {
278852e5969SJed Brown     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", user->app_ctx->output_dir, step_no));
27991a36801SJames Wright   } else {
2802b916ea7SJeremy L Thompson     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", user->app_ctx->output_dir));
28191a36801SJames Wright   }
2822b916ea7SJeremy L Thompson   PetscCall(PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, &viewer));
2837538d537SJames Wright 
284e1233009SJames Wright   PetscInt32 token = PetscDefined(USE_64BIT_INDICES) ? FLUIDS_FILE_TOKEN_64 : FLUIDS_FILE_TOKEN_32;
285e1233009SJames Wright   PetscCall(PetscViewerBinaryWrite(viewer, &token, 1, PETSC_INT32));
2869293eaa1SJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &step_no, 1, PETSC_INT));
2879293eaa1SJed Brown   time /= user->units->second;  // Dimensionalize time back
2889293eaa1SJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL));
2897538d537SJames Wright   PetscCall(VecView(Q, viewer));
2907538d537SJames Wright   PetscCall(PetscViewerDestroy(&viewer));
291d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2927538d537SJames Wright }
2937538d537SJames Wright 
294c5e9980aSAdeleke O. Bankole // CSV Monitor
295c5e9980aSAdeleke O. Bankole PetscErrorCode TSMonitor_WallForce(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
296c5e9980aSAdeleke O. Bankole   User              user = ctx;
297c5e9980aSAdeleke O. Bankole   Vec               G_loc;
298c5e9980aSAdeleke O. Bankole   PetscInt          num_wall = user->app_ctx->wall_forces.num_wall, dim = 3;
299c5e9980aSAdeleke O. Bankole   const PetscInt   *walls  = user->app_ctx->wall_forces.walls;
300c5e9980aSAdeleke O. Bankole   PetscViewer       viewer = user->app_ctx->wall_forces.viewer;
301c5e9980aSAdeleke O. Bankole   PetscViewerFormat format = user->app_ctx->wall_forces.viewer_format;
302c5e9980aSAdeleke O. Bankole   PetscScalar      *reaction_force;
303c5e9980aSAdeleke O. Bankole   PetscBool         iascii;
304c5e9980aSAdeleke O. Bankole 
305c5e9980aSAdeleke O. Bankole   PetscFunctionBeginUser;
306d949ddfcSJames Wright   if (!viewer) PetscFunctionReturn(PETSC_SUCCESS);
307c5e9980aSAdeleke O. Bankole   PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
308c5e9980aSAdeleke O. Bankole   PetscCall(PetscMalloc1(num_wall * dim, &reaction_force));
309c5e9980aSAdeleke O. Bankole   PetscCall(Surface_Forces_NS(user->dm, G_loc, num_wall, walls, reaction_force));
310c5e9980aSAdeleke O. Bankole   PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
311c5e9980aSAdeleke O. Bankole 
312c5e9980aSAdeleke O. Bankole   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
313c5e9980aSAdeleke O. Bankole 
314c5e9980aSAdeleke O. Bankole   if (iascii) {
315c5e9980aSAdeleke O. Bankole     if (format == PETSC_VIEWER_ASCII_CSV && !user->app_ctx->wall_forces.header_written) {
316c5e9980aSAdeleke O. Bankole       PetscCall(PetscViewerASCIIPrintf(viewer, "Step,Time,Wall,ForceX,ForceY,ForceZ\n"));
317c5e9980aSAdeleke O. Bankole       user->app_ctx->wall_forces.header_written = PETSC_TRUE;
318c5e9980aSAdeleke O. Bankole     }
319c5e9980aSAdeleke O. Bankole     for (PetscInt w = 0; w < num_wall; w++) {
320c5e9980aSAdeleke O. Bankole       PetscInt wall = walls[w];
321c5e9980aSAdeleke O. Bankole       if (format == PETSC_VIEWER_ASCII_CSV) {
322c5e9980aSAdeleke O. Bankole         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT ",%g,%" PetscInt_FMT ",%g,%g,%g\n", step_no, time, wall,
323c5e9980aSAdeleke O. Bankole                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
324c5e9980aSAdeleke O. Bankole 
325c5e9980aSAdeleke O. Bankole       } else {
326c5e9980aSAdeleke O. Bankole         PetscCall(PetscViewerASCIIPrintf(viewer, "Wall %" PetscInt_FMT " Forces: Force_x = %12g, Force_y = %12g, Force_z = %12g\n", wall,
327c5e9980aSAdeleke O. Bankole                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
328c5e9980aSAdeleke O. Bankole       }
329c5e9980aSAdeleke O. Bankole     }
330c5e9980aSAdeleke O. Bankole   }
331c5e9980aSAdeleke O. Bankole   PetscCall(PetscFree(reaction_force));
332d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
333c5e9980aSAdeleke O. Bankole }
334c5e9980aSAdeleke O. Bankole 
3357538d537SJames Wright // User provided TS Monitor
3362b916ea7SJeremy L Thompson PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
3377538d537SJames Wright   User user = ctx;
3387538d537SJames Wright 
33906f41313SJames Wright   PetscFunctionBeginUser;
340852e5969SJed Brown   // Print every 'checkpoint_interval' steps
341c539088bSJames Wright   if (user->app_ctx->checkpoint_interval <= 0 || step_no % user->app_ctx->checkpoint_interval != 0 ||
342e419654dSJeremy L Thompson       (user->app_ctx->cont_steps == step_no && step_no != 0)) {
343d949ddfcSJames Wright     PetscFunctionReturn(PETSC_SUCCESS);
344e419654dSJeremy L Thompson   }
3457538d537SJames Wright 
3467538d537SJames Wright   PetscCall(WriteOutput(user, Q, step_no, time));
347d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
348a515125bSLeila Ghaffari }
349a515125bSLeila Ghaffari 
350a515125bSLeila Ghaffari // TS: Create, setup, and solve
3512b916ea7SJeremy L Thompson PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, Vec *Q, PetscScalar *f_time, TS *ts) {
352a515125bSLeila Ghaffari   MPI_Comm    comm = user->comm;
353a515125bSLeila Ghaffari   TSAdapt     adapt;
354a515125bSLeila Ghaffari   PetscScalar final_time;
355a515125bSLeila Ghaffari 
35606f41313SJames Wright   PetscFunctionBeginUser;
3572b916ea7SJeremy L Thompson   PetscCall(TSCreate(comm, ts));
3582b916ea7SJeremy L Thompson   PetscCall(TSSetDM(*ts, dm));
359632a41e1SJames Wright   PetscCall(TSSetApplicationContext(*ts, user));
360a515125bSLeila Ghaffari   if (phys->implicit) {
3612b916ea7SJeremy L Thompson     PetscCall(TSSetType(*ts, TSBDF));
362a515125bSLeila Ghaffari     if (user->op_ifunction) {
3632b916ea7SJeremy L Thompson       PetscCall(TSSetIFunction(*ts, NULL, IFunction_NS, &user));
364a515125bSLeila Ghaffari     } else {  // Implicit integrators can fall back to using an RHSFunction
3652b916ea7SJeremy L Thompson       PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
366a515125bSLeila Ghaffari     }
367*ebfabadfSJames Wright     if (user->mat_ijacobian) {
3682b916ea7SJeremy L Thompson       PetscCall(DMTSSetIJacobian(dm, FormIJacobian_NS, &user));
369f0b65372SJed Brown     }
370a515125bSLeila Ghaffari   } else {
371da5fe0e4SJames Wright     PetscCheck(user->op_rhs_ctx, comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction");
3722b916ea7SJeremy L Thompson     PetscCall(TSSetType(*ts, TSRK));
3732b916ea7SJeremy L Thompson     PetscCall(TSRKSetType(*ts, TSRK5F));
3742b916ea7SJeremy L Thompson     PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
375a515125bSLeila Ghaffari   }
3762b916ea7SJeremy L Thompson   PetscCall(TSSetMaxTime(*ts, 500. * user->units->second));
3772b916ea7SJeremy L Thompson   PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER));
37822387d3aSJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(TSSetErrorIfStepFails(*ts, PETSC_FALSE));
3792b916ea7SJeremy L Thompson   PetscCall(TSSetTimeStep(*ts, 1.e-2 * user->units->second));
3802b916ea7SJeremy L Thompson   PetscCall(TSGetAdapt(*ts, &adapt));
3812b916ea7SJeremy L Thompson   PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * user->units->second, 1.e2 * user->units->second));
3822b916ea7SJeremy L Thompson   PetscCall(TSSetFromOptions(*ts));
383*ebfabadfSJames Wright   if (user->mat_ijacobian) {
384*ebfabadfSJames Wright     if (app_ctx->amat_type && !strcmp(app_ctx->amat_type, MATSHELL)) {
385*ebfabadfSJames Wright       PetscBool use_matceed_pmat;
386*ebfabadfSJames Wright       SNES      snes;
387*ebfabadfSJames Wright       KSP       ksp;
388*ebfabadfSJames Wright       PC        pc;
389*ebfabadfSJames Wright       PCType    pc_type;
390*ebfabadfSJames Wright       Mat       Pmat;
391*ebfabadfSJames Wright 
392*ebfabadfSJames Wright       PetscCall(TSGetSNES(*ts, &snes));
393*ebfabadfSJames Wright       PetscCall(SNESGetKSP(snes, &ksp));
394*ebfabadfSJames Wright       PetscCall(KSPGetPC(ksp, &pc));
395*ebfabadfSJames Wright       PetscCall(PCGetType(pc, &pc_type));
396*ebfabadfSJames Wright       PetscCall(PetscStrcmpAny(pc_type, &use_matceed_pmat, PCJACOBI, PCVPBJACOBI, PCPBJACOBI, ""));
397*ebfabadfSJames Wright       Pmat = use_matceed_pmat ? user->mat_ijacobian : NULL;
398*ebfabadfSJames Wright       PetscCall(TSSetIJacobian(*ts, user->mat_ijacobian, Pmat, NULL, NULL));
399*ebfabadfSJames Wright     }
400*ebfabadfSJames Wright   }
40191f639d2SJames Wright   user->time_bc_set = -1.0;   // require all BCs be updated
402c26b555cSJames Wright   if (app_ctx->cont_steps) {  // continue from previous timestep data
403a515125bSLeila Ghaffari     PetscInt    count;
404a515125bSLeila Ghaffari     PetscViewer viewer;
4052b916ea7SJeremy L Thompson 
4069293eaa1SJed Brown     if (app_ctx->cont_time <= 0) {  // Legacy files did not include step number and time
4072b916ea7SJeremy L Thompson       PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_time_file, FILE_MODE_READ, &viewer));
4089293eaa1SJed Brown       PetscCall(PetscViewerBinaryRead(viewer, &app_ctx->cont_time, 1, &count, PETSC_REAL));
4092b916ea7SJeremy L Thompson       PetscCall(PetscViewerDestroy(&viewer));
4109293eaa1SJed Brown       PetscCheck(app_ctx->cont_steps != -1, comm, PETSC_ERR_ARG_INCOMP,
4119293eaa1SJed Brown                  "-continue step number not specified, but checkpoint file does not contain a step number (likely written by older code version)");
4129293eaa1SJed Brown     }
4139293eaa1SJed Brown     PetscCall(TSSetTime(*ts, app_ctx->cont_time * user->units->second));
41474a6f4ddSJed Brown     PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps));
415a515125bSLeila Ghaffari   }
4160e1e9333SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) {
4172b916ea7SJeremy L Thompson     PetscCall(TSMonitorSet(*ts, TSMonitor_NS, user, NULL));
4180e1e9333SJames Wright   }
419c5e9980aSAdeleke O. Bankole   if (app_ctx->wall_forces.viewer) {
420c5e9980aSAdeleke O. Bankole     PetscCall(TSMonitorSet(*ts, TSMonitor_WallForce, user, NULL));
421c5e9980aSAdeleke O. Bankole   }
422c931fa59SJames Wright   if (app_ctx->turb_spanstats_enable) {
42391933550SJames Wright     PetscCall(TSMonitorSet(*ts, TSMonitor_TurbulenceStatistics, user, NULL));
424b8daee98SJames Wright     CeedScalar previous_time = app_ctx->cont_time * user->units->second;
425b4c37c5cSJames Wright     PetscCallCeed(user->ceed,
426b4c37c5cSJames Wright                   CeedOperatorSetContextDouble(user->spanstats.op_stats_collect_ctx->op, user->spanstats.previous_time_label, &previous_time));
427b0488d1fSJames Wright   }
42888b07121SJames Wright   if (app_ctx->diff_filter_monitor) PetscCall(TSMonitorSet(*ts, TSMonitor_DifferentialFilter, user, NULL));
429a515125bSLeila Ghaffari 
4301c17f66aSJames Wright   if (app_ctx->sgs_train_enable) {
4311c17f66aSJames Wright     PetscCall(TSMonitorSet(*ts, TSMonitor_SGS_DD_Training, user, NULL));
4321c17f66aSJames Wright     PetscCall(TSSetPostStep(*ts, TSPostStep_SGS_DD_Training));
4331c17f66aSJames Wright   }
434a515125bSLeila Ghaffari   // Solve
43574a6f4ddSJed Brown   PetscReal start_time;
43674a6f4ddSJed Brown   PetscInt  start_step;
4372b916ea7SJeremy L Thompson   PetscCall(TSGetTime(*ts, &start_time));
43874a6f4ddSJed Brown   PetscCall(TSGetStepNumber(*ts, &start_step));
43991982731SJeremy L Thompson 
440df4304b5SJed Brown   PetscCall(PetscLogDefaultBegin());  // So we can use PetscLogStageGetPerfInfo without -log_view
44191982731SJeremy L Thompson   PetscPreLoadBegin(PETSC_FALSE, "Fluids Solve");
44291982731SJeremy L Thompson   PetscCall(TSSetTime(*ts, start_time));
44374a6f4ddSJed Brown   PetscCall(TSSetStepNumber(*ts, start_step));
44491982731SJeremy L Thompson   if (PetscPreLoadingOn) {
44591982731SJeremy L Thompson     // LCOV_EXCL_START
44691982731SJeremy L Thompson     SNES      snes;
44791982731SJeremy L Thompson     Vec       Q_preload;
44891982731SJeremy L Thompson     PetscReal rtol;
44991982731SJeremy L Thompson     PetscCall(VecDuplicate(*Q, &Q_preload));
45091982731SJeremy L Thompson     PetscCall(VecCopy(*Q, Q_preload));
45191982731SJeremy L Thompson     PetscCall(TSGetSNES(*ts, &snes));
45291982731SJeremy L Thompson     PetscCall(SNESGetTolerances(snes, NULL, &rtol, NULL, NULL, NULL));
4532b916ea7SJeremy L Thompson     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, .99, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
45422c1b34eSJames Wright     PetscCall(TSSetSolution(*ts, Q_preload));
45591982731SJeremy L Thompson     PetscCall(TSStep(*ts));
4562b916ea7SJeremy L Thompson     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, rtol, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
45791982731SJeremy L Thompson     PetscCall(VecDestroy(&Q_preload));
45891982731SJeremy L Thompson     // LCOV_EXCL_STOP
45991982731SJeremy L Thompson   } else {
4602b916ea7SJeremy L Thompson     PetscCall(PetscBarrier((PetscObject)*ts));
4612b916ea7SJeremy L Thompson     PetscCall(TSSolve(*ts, *Q));
46291982731SJeremy L Thompson   }
46391982731SJeremy L Thompson   PetscPreLoadEnd();
46491982731SJeremy L Thompson 
46591982731SJeremy L Thompson   PetscCall(TSGetSolveTime(*ts, &final_time));
466a515125bSLeila Ghaffari   *f_time = final_time;
46791982731SJeremy L Thompson 
4680e1e9333SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) {
4697538d537SJames Wright     PetscInt step_no;
4707538d537SJames Wright     PetscCall(TSGetStepNumber(*ts, &step_no));
471b0488d1fSJames Wright     if (user->app_ctx->checkpoint_interval > 0 || user->app_ctx->checkpoint_interval == -1) {
4727538d537SJames Wright       PetscCall(WriteOutput(user, *Q, step_no, final_time));
4737538d537SJames Wright     }
4747538d537SJames Wright 
4757eedc94cSJames Wright     PetscLogStage      stage_id;
476df4304b5SJed Brown     PetscEventPerfInfo stage_perf;
47791982731SJeremy L Thompson 
47891982731SJeremy L Thompson     PetscCall(PetscLogStageGetId("Fluids Solve", &stage_id));
479df4304b5SJed Brown     PetscCall(PetscLogStageGetPerfInfo(stage_id, &stage_perf));
480df4304b5SJed Brown     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time taken for solution (sec): %g\n", stage_perf.time));
481a515125bSLeila Ghaffari   }
482d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
483a515125bSLeila Ghaffari }
484