xref: /libCEED/examples/fluids/src/setupts.c (revision d0593705e733b5bdd5e4c173fe0008b11db2ed29)
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 
1149aac155SJeremy L Thompson #include <ceed.h>
1249aac155SJeremy L Thompson #include <petscdmplex.h>
1349aac155SJeremy L Thompson #include <petscts.h>
1449aac155SJeremy L Thompson 
1577841947SLeila Ghaffari #include "../navierstokes.h"
16ca69d878SAdeleke O. Bankole #include "../qfunctions/newtonian_state.h"
1777841947SLeila Ghaffari 
1804292f7dSJames Wright // @brief Create KSP to solve the inverse mass operator for explicit time stepping schemes
19dfcf44b0SJames Wright PetscErrorCode CreateKSPMassOperator(User user, CeedData ceed_data) {
2004292f7dSJames Wright   Ceed          ceed = user->ceed;
2104292f7dSJames Wright   DM            dm   = user->dm;
2277841947SLeila Ghaffari   CeedQFunction qf_mass;
2377841947SLeila Ghaffari   CeedOperator  op_mass;
2477841947SLeila Ghaffari   CeedInt       num_comp_q, q_data_size;
2577841947SLeila Ghaffari 
26f17d818dSJames Wright   PetscFunctionBeginUser;
27a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q));
28a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &q_data_size));
2977841947SLeila Ghaffari 
30ef080ff9SJames Wright   PetscCall(CreateMassQFunction(ceed, num_comp_q, q_data_size, &qf_mass));
31a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass));
32a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
33356036faSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
34a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
3577841947SLeila Ghaffari 
3604292f7dSJames Wright   {  // -- Setup KSP for mass operator
377f2a9303SJames Wright     Mat      mat_mass;
380f2fa9b4SJames Wright     Vec      Zeros_loc;
3904292f7dSJames Wright     MPI_Comm comm = PetscObjectComm((PetscObject)dm);
4077841947SLeila Ghaffari 
410f2fa9b4SJames Wright     PetscCall(DMCreateLocalVector(dm, &Zeros_loc));
420f2fa9b4SJames Wright     PetscCall(VecZeroEntries(Zeros_loc));
437f2a9303SJames Wright     PetscCall(MatCeedCreate(dm, dm, op_mass, NULL, &mat_mass));
447f2a9303SJames Wright     PetscCall(MatCeedSetLocalVectors(mat_mass, Zeros_loc, NULL));
4577841947SLeila Ghaffari 
4604292f7dSJames Wright     PetscCall(KSPCreate(comm, &user->mass_ksp));
4704292f7dSJames Wright     PetscCall(KSPSetOptionsPrefix(user->mass_ksp, "mass_"));
4804292f7dSJames Wright     {  // lumped by default
4904292f7dSJames Wright       PC pc;
5004292f7dSJames Wright       PetscCall(KSPGetPC(user->mass_ksp, &pc));
5104292f7dSJames Wright       PetscCall(PCSetType(pc, PCJACOBI));
5204292f7dSJames Wright       PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
5304292f7dSJames Wright       PetscCall(KSPSetType(user->mass_ksp, KSPPREONLY));
5404292f7dSJames Wright     }
557f2a9303SJames Wright     PetscCall(KSPSetFromOptions_WithMatCeed(user->mass_ksp, mat_mass));
560f2fa9b4SJames Wright     PetscCall(VecDestroy(&Zeros_loc));
570f2fa9b4SJames Wright     PetscCall(MatDestroy(&mat_mass));
5804292f7dSJames Wright   }
5977841947SLeila Ghaffari 
60a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
61a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
62ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
6377841947SLeila Ghaffari }
6477841947SLeila Ghaffari 
65a0b9a424SJames Wright // Insert Boundary values if it's a new time
66a0b9a424SJames Wright PetscErrorCode UpdateBoundaryValues(User user, Vec Q_loc, PetscReal t) {
67a0b9a424SJames Wright   PetscFunctionBeginUser;
68a0b9a424SJames Wright   if (user->time_bc_set != t) {
69a0b9a424SJames Wright     PetscCall(DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL));
70a0b9a424SJames Wright     user->time_bc_set = t;
71a0b9a424SJames Wright   }
72ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
73a0b9a424SJames Wright }
74a0b9a424SJames Wright 
7577841947SLeila Ghaffari // RHS (Explicit time-stepper) function setup
7677841947SLeila Ghaffari //   This is the RHS of the ODE, given as u_t = G(t,u)
7777841947SLeila Ghaffari //   This function takes in a state vector Q and writes into G
7877841947SLeila Ghaffari PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) {
7977841947SLeila Ghaffari   User        user = *(User *)user_data;
807d4c6defSJames Wright   Ceed        ceed = user->ceed;
81c798d55aSJames Wright   PetscScalar dt;
829ad5e8e4SJames Wright   Vec         Q_loc = user->Q_loc;
8377841947SLeila Ghaffari 
84f17d818dSJames Wright   PetscFunctionBeginUser;
855e82a6e1SJeremy L Thompson   // Update time dependent data
86a0b9a424SJames Wright   PetscCall(UpdateBoundaryValues(user, Q_loc, t));
877d4c6defSJames Wright   if (user->phys->solution_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_rhs_ctx->op, user->phys->solution_time_label, &t));
882b730f8bSJeremy L Thompson   PetscCall(TSGetTimeStep(ts, &dt));
897d4c6defSJames Wright   if (user->phys->timestep_size_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_rhs_ctx->op, user->phys->timestep_size_label, &dt));
9077841947SLeila Ghaffari 
919ad5e8e4SJames Wright   PetscCall(ApplyCeedOperatorGlobalToGlobal(Q, G, user->op_rhs_ctx));
9277841947SLeila Ghaffari 
93d6e67e47SJames Wright   // Inverse of the lumped mass matrix
9404292f7dSJames Wright   PetscCall(KSPSolve(user->mass_ksp, G, G));
95ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
9677841947SLeila Ghaffari }
9777841947SLeila Ghaffari 
98ca69d878SAdeleke O. Bankole // Surface forces function setup
99ca69d878SAdeleke O. Bankole static PetscErrorCode Surface_Forces_NS(DM dm, Vec G_loc, PetscInt num_walls, const PetscInt walls[], PetscScalar *reaction_force) {
100ca69d878SAdeleke O. Bankole   DMLabel            face_label;
101ca69d878SAdeleke O. Bankole   const PetscScalar *g;
102d6734f85SAdeleke O. Bankole   PetscInt           dof, dim = 3;
103ca69d878SAdeleke O. Bankole   MPI_Comm           comm;
104d6734f85SAdeleke O. Bankole   PetscSection       s;
105ca69d878SAdeleke O. Bankole 
106ca69d878SAdeleke O. Bankole   PetscFunctionBeginUser;
107ca69d878SAdeleke O. Bankole   PetscCall(PetscArrayzero(reaction_force, num_walls * dim));
108ca69d878SAdeleke O. Bankole   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
109ca69d878SAdeleke O. Bankole   PetscCall(DMGetLabel(dm, "Face Sets", &face_label));
110ca69d878SAdeleke O. Bankole   PetscCall(VecGetArrayRead(G_loc, &g));
111ca69d878SAdeleke O. Bankole   for (PetscInt w = 0; w < num_walls; w++) {
112ca69d878SAdeleke O. Bankole     const PetscInt wall = walls[w];
113ca69d878SAdeleke O. Bankole     IS             wall_is;
114d6734f85SAdeleke O. Bankole     PetscCall(DMGetLocalSection(dm, &s));
115ca69d878SAdeleke O. Bankole     PetscCall(DMLabelGetStratumIS(face_label, wall, &wall_is));
116ca69d878SAdeleke O. Bankole     if (wall_is) {  // There exist such points on this process
117ca69d878SAdeleke O. Bankole       PetscInt        num_points;
118d6734f85SAdeleke O. Bankole       PetscInt        num_comp = 0;
119ca69d878SAdeleke O. Bankole       const PetscInt *points;
120d6734f85SAdeleke O. Bankole       PetscCall(PetscSectionGetFieldComponents(s, 0, &num_comp));
121ca69d878SAdeleke O. Bankole       PetscCall(ISGetSize(wall_is, &num_points));
122ca69d878SAdeleke O. Bankole       PetscCall(ISGetIndices(wall_is, &points));
123ca69d878SAdeleke O. Bankole       for (PetscInt i = 0; i < num_points; i++) {
124ca69d878SAdeleke O. Bankole         const PetscInt           p = points[i];
125ca69d878SAdeleke O. Bankole         const StateConservative *r;
126ca69d878SAdeleke O. Bankole         PetscCall(DMPlexPointLocalRead(dm, p, g, &r));
127d6734f85SAdeleke O. Bankole         PetscCall(PetscSectionGetDof(s, p, &dof));
128d6734f85SAdeleke O. Bankole         for (PetscInt node = 0; node < dof / num_comp; node++) {
129ca69d878SAdeleke O. Bankole           for (PetscInt j = 0; j < 3; j++) {
130d6734f85SAdeleke O. Bankole             reaction_force[w * dim + j] -= r[node].momentum[j];
131d6734f85SAdeleke O. Bankole           }
132ca69d878SAdeleke O. Bankole         }
133ca69d878SAdeleke O. Bankole       }
134ca69d878SAdeleke O. Bankole       PetscCall(ISRestoreIndices(wall_is, &points));
135ca69d878SAdeleke O. Bankole     }
136ca69d878SAdeleke O. Bankole     PetscCall(ISDestroy(&wall_is));
137ca69d878SAdeleke O. Bankole   }
138ca69d878SAdeleke O. Bankole   PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, reaction_force, dim * num_walls, MPIU_SCALAR, MPI_SUM, comm));
139ca69d878SAdeleke O. Bankole   //  Restore Vectors
140ca69d878SAdeleke O. Bankole   PetscCall(VecRestoreArrayRead(G_loc, &g));
141ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
142ca69d878SAdeleke O. Bankole }
143ca69d878SAdeleke O. Bankole 
14477841947SLeila Ghaffari // Implicit time-stepper function setup
1452b730f8bSJeremy L Thompson PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) {
14677841947SLeila Ghaffari   User         user = *(User *)user_data;
1477d4c6defSJames Wright   Ceed         ceed = user->ceed;
148c798d55aSJames Wright   PetscScalar  dt;
1495e82a6e1SJeremy L Thompson   Vec          Q_loc = user->Q_loc, Q_dot_loc = user->Q_dot_loc, G_loc;
15077841947SLeila Ghaffari   PetscMemType q_mem_type, q_dot_mem_type, g_mem_type;
15177841947SLeila Ghaffari 
152f17d818dSJames Wright   PetscFunctionBeginUser;
1535e82a6e1SJeremy L Thompson   // Get local vectors
154ca69d878SAdeleke O. Bankole   PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
1555e82a6e1SJeremy L Thompson 
1565e82a6e1SJeremy L Thompson   // Update time dependent data
157a0b9a424SJames Wright   PetscCall(UpdateBoundaryValues(user, Q_loc, t));
1587d4c6defSJames Wright   if (user->phys->solution_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_ifunction, user->phys->solution_time_label, &t));
1592b730f8bSJeremy L Thompson   PetscCall(TSGetTimeStep(ts, &dt));
1607d4c6defSJames Wright   if (user->phys->timestep_size_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_ifunction, user->phys->timestep_size_label, &dt));
16177841947SLeila Ghaffari 
16277841947SLeila Ghaffari   // Global-to-local
163878eb0e7SJames Wright   PetscCall(DMGlobalToLocalBegin(user->dm, Q, INSERT_VALUES, Q_loc));
164878eb0e7SJames Wright   PetscCall(DMGlobalToLocalBegin(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
165878eb0e7SJames Wright   PetscCall(DMGlobalToLocalEnd(user->dm, Q, INSERT_VALUES, Q_loc));
166878eb0e7SJames Wright   PetscCall(DMGlobalToLocalEnd(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
16777841947SLeila Ghaffari 
16877841947SLeila Ghaffari   // Place PETSc vectors in CEED vectors
169*d0593705SJames Wright   PetscCall(VecReadPetscToCeed(Q_loc, &q_mem_type, user->q_ceed));
170*d0593705SJames Wright   PetscCall(VecReadPetscToCeed(Q_dot_loc, &q_dot_mem_type, user->q_dot_ceed));
171*d0593705SJames Wright   PetscCall(VecPetscToCeed(G_loc, &g_mem_type, user->g_ceed));
17277841947SLeila Ghaffari 
17377841947SLeila Ghaffari   // Apply CEED operator
17475d1798cSJames Wright   PetscCall(PetscLogEventBegin(FLUIDS_CeedOperatorApply, Q, G, 0, 0));
17575d1798cSJames Wright   PetscCall(PetscLogGpuTimeBegin());
176a424bcd0SJames Wright   PetscCallCeed(user->ceed, CeedOperatorApply(user->op_ifunction, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE));
17775d1798cSJames Wright   PetscCall(PetscLogGpuTimeEnd());
17875d1798cSJames Wright   PetscCall(PetscLogEventEnd(FLUIDS_CeedOperatorApply, Q, G, 0, 0));
17977841947SLeila Ghaffari 
18077841947SLeila Ghaffari   // Restore vectors
181*d0593705SJames Wright   PetscCall(VecReadCeedToPetsc(user->q_ceed, q_mem_type, Q_loc));
182*d0593705SJames Wright   PetscCall(VecReadCeedToPetsc(user->q_dot_ceed, q_dot_mem_type, Q_dot_loc));
183*d0593705SJames Wright   PetscCall(VecCeedToPetsc(user->g_ceed, g_mem_type, G_loc));
18477841947SLeila Ghaffari 
18519ffbc25SJames Wright   if (user->app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) {
1868f5ab23bSJames Wright     PetscCall(SgsDDApplyIFunction(user, Q_loc, G_loc));
18719ffbc25SJames Wright   }
188be34b3b7SJames Wright 
18977841947SLeila Ghaffari   // Local-to-Global
1902b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(G));
1912b730f8bSJeremy L Thompson   PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G));
19277841947SLeila Ghaffari 
19377841947SLeila Ghaffari   // Restore vectors
194ca69d878SAdeleke O. Bankole   PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
195ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
19677841947SLeila Ghaffari }
19777841947SLeila Ghaffari 
1982b730f8bSJeremy L Thompson PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, PetscReal shift, Mat J, Mat J_pre, void *user_data) {
199e334ad8fSJed Brown   User      user = *(User *)user_data;
200a424bcd0SJames Wright   Ceed      ceed = user->ceed;
20191c97f41SJames Wright   PetscBool J_is_matceed, J_is_mffd, J_pre_is_matceed, J_pre_is_mffd;
202f17d818dSJames Wright 
203e334ad8fSJed Brown   PetscFunctionBeginUser;
204d6bf345cSJed Brown   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATMFFD, &J_is_mffd));
20591c97f41SJames Wright   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATCEED, &J_is_matceed));
20691c97f41SJames Wright   PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATMFFD, &J_pre_is_mffd));
20791c97f41SJames Wright   PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATCEED, &J_pre_is_matceed));
20891c97f41SJames Wright   if (user->phys->ijacobian_time_shift_label) {
20991c97f41SJames Wright     CeedOperator op_ijacobian;
21091c97f41SJames Wright 
21191c97f41SJames Wright     PetscCall(MatCeedGetCeedOperators(user->mat_ijacobian, &op_ijacobian, NULL));
21291c97f41SJames Wright     PetscCallCeed(ceed, CeedOperatorSetContextDouble(op_ijacobian, user->phys->ijacobian_time_shift_label, &shift));
213e334ad8fSJed Brown   }
21491c97f41SJames Wright 
21591c97f41SJames Wright   if (J_is_matceed || J_is_mffd) {
216d6bf345cSJed Brown     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
217d6bf345cSJed Brown     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
21891c97f41SJames Wright   } else PetscCall(MatCeedAssembleCOO(user->mat_ijacobian, J));
21991c97f41SJames Wright 
22091c97f41SJames Wright   if (J_pre_is_matceed && J != J_pre) {
22191c97f41SJames Wright     PetscCall(MatAssemblyBegin(J_pre, MAT_FINAL_ASSEMBLY));
22291c97f41SJames Wright     PetscCall(MatAssemblyEnd(J_pre, MAT_FINAL_ASSEMBLY));
22391c97f41SJames Wright   } else if (!J_pre_is_matceed && !J_pre_is_mffd && J != J_pre) {
22491c97f41SJames Wright     PetscCall(MatCeedAssembleCOO(user->mat_ijacobian, J_pre));
225d6bf345cSJed Brown   }
226ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
227e334ad8fSJed Brown }
228e334ad8fSJed Brown 
2292b730f8bSJeremy L Thompson PetscErrorCode WriteOutput(User user, Vec Q, PetscInt step_no, PetscScalar time) {
23077841947SLeila Ghaffari   Vec         Q_loc;
23177841947SLeila Ghaffari   char        file_path[PETSC_MAX_PATH_LEN];
23277841947SLeila Ghaffari   PetscViewer viewer;
23377841947SLeila Ghaffari 
234f17d818dSJames Wright   PetscFunctionBeginUser;
23537cbb16aSJed Brown   if (user->app_ctx->checkpoint_vtk) {
23677841947SLeila Ghaffari     // Set up output
237f14660a4SJames Wright     PetscCall(DMGetLocalVector(user->dm, &Q_loc));
238f14660a4SJames Wright     PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec"));
239f14660a4SJames Wright     PetscCall(VecZeroEntries(Q_loc));
240f14660a4SJames Wright     PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
24177841947SLeila Ghaffari 
24277841947SLeila Ghaffari     // Output
24337cbb16aSJed Brown     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
244f14660a4SJames Wright 
2452b730f8bSJeremy L Thompson     PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path, FILE_MODE_WRITE, &viewer));
246f14660a4SJames Wright     PetscCall(VecView(Q_loc, viewer));
247f14660a4SJames Wright     PetscCall(PetscViewerDestroy(&viewer));
24877841947SLeila Ghaffari     if (user->dm_viz) {
24977841947SLeila Ghaffari       Vec         Q_refined, Q_refined_loc;
25077841947SLeila Ghaffari       char        file_path_refined[PETSC_MAX_PATH_LEN];
25177841947SLeila Ghaffari       PetscViewer viewer_refined;
25277841947SLeila Ghaffari 
253f14660a4SJames Wright       PetscCall(DMGetGlobalVector(user->dm_viz, &Q_refined));
254f14660a4SJames Wright       PetscCall(DMGetLocalVector(user->dm_viz, &Q_refined_loc));
255f14660a4SJames Wright       PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined"));
256f14660a4SJames Wright 
257f14660a4SJames Wright       PetscCall(MatInterpolate(user->interp_viz, Q, Q_refined));
258f14660a4SJames Wright       PetscCall(VecZeroEntries(Q_refined_loc));
2592b730f8bSJeremy L Thompson       PetscCall(DMGlobalToLocal(user->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc));
260f14660a4SJames Wright 
26137cbb16aSJed Brown       PetscCall(
26237cbb16aSJed Brown           PetscSNPrintf(file_path_refined, sizeof file_path_refined, "%s/nsrefined-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
263f14660a4SJames Wright 
2642b730f8bSJeremy L Thompson       PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined), file_path_refined, FILE_MODE_WRITE, &viewer_refined));
265f14660a4SJames Wright       PetscCall(VecView(Q_refined_loc, viewer_refined));
266f14660a4SJames Wright       PetscCall(DMRestoreLocalVector(user->dm_viz, &Q_refined_loc));
267f14660a4SJames Wright       PetscCall(DMRestoreGlobalVector(user->dm_viz, &Q_refined));
268f14660a4SJames Wright       PetscCall(PetscViewerDestroy(&viewer_refined));
26977841947SLeila Ghaffari     }
270f14660a4SJames Wright     PetscCall(DMRestoreLocalVector(user->dm, &Q_loc));
27137cbb16aSJed Brown   }
27277841947SLeila Ghaffari 
27377841947SLeila Ghaffari   // Save data in a binary file for continuation of simulations
27469293791SJames Wright   if (user->app_ctx->add_stepnum2bin) {
27537cbb16aSJed Brown     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", user->app_ctx->output_dir, step_no));
27669293791SJames Wright   } else {
2772b730f8bSJeremy L Thompson     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", user->app_ctx->output_dir));
27869293791SJames Wright   }
2792b730f8bSJeremy L Thompson   PetscCall(PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, &viewer));
280f14660a4SJames Wright 
2810de6a49fSJames Wright   PetscInt32 token = PetscDefined(USE_64BIT_INDICES) ? FLUIDS_FILE_TOKEN_64 : FLUIDS_FILE_TOKEN_32;
2820de6a49fSJames Wright   PetscCall(PetscViewerBinaryWrite(viewer, &token, 1, PETSC_INT32));
2834de8550aSJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &step_no, 1, PETSC_INT));
2844de8550aSJed Brown   time /= user->units->second;  // Dimensionalize time back
2854de8550aSJed Brown   PetscCall(PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL));
286f14660a4SJames Wright   PetscCall(VecView(Q, viewer));
287f14660a4SJames Wright   PetscCall(PetscViewerDestroy(&viewer));
288ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
289f14660a4SJames Wright }
290f14660a4SJames Wright 
291ca69d878SAdeleke O. Bankole // CSV Monitor
292ca69d878SAdeleke O. Bankole PetscErrorCode TSMonitor_WallForce(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
293ca69d878SAdeleke O. Bankole   User              user = ctx;
294ca69d878SAdeleke O. Bankole   Vec               G_loc;
295ca69d878SAdeleke O. Bankole   PetscInt          num_wall = user->app_ctx->wall_forces.num_wall, dim = 3;
296ca69d878SAdeleke O. Bankole   const PetscInt   *walls  = user->app_ctx->wall_forces.walls;
297ca69d878SAdeleke O. Bankole   PetscViewer       viewer = user->app_ctx->wall_forces.viewer;
298ca69d878SAdeleke O. Bankole   PetscViewerFormat format = user->app_ctx->wall_forces.viewer_format;
299ca69d878SAdeleke O. Bankole   PetscScalar      *reaction_force;
300ca69d878SAdeleke O. Bankole   PetscBool         iascii;
301ca69d878SAdeleke O. Bankole 
302ca69d878SAdeleke O. Bankole   PetscFunctionBeginUser;
303ee4ca9cbSJames Wright   if (!viewer) PetscFunctionReturn(PETSC_SUCCESS);
304ca69d878SAdeleke O. Bankole   PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
305ca69d878SAdeleke O. Bankole   PetscCall(PetscMalloc1(num_wall * dim, &reaction_force));
306ca69d878SAdeleke O. Bankole   PetscCall(Surface_Forces_NS(user->dm, G_loc, num_wall, walls, reaction_force));
307ca69d878SAdeleke O. Bankole   PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
308ca69d878SAdeleke O. Bankole 
309ca69d878SAdeleke O. Bankole   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
310ca69d878SAdeleke O. Bankole 
311ca69d878SAdeleke O. Bankole   if (iascii) {
312ca69d878SAdeleke O. Bankole     if (format == PETSC_VIEWER_ASCII_CSV && !user->app_ctx->wall_forces.header_written) {
313ca69d878SAdeleke O. Bankole       PetscCall(PetscViewerASCIIPrintf(viewer, "Step,Time,Wall,ForceX,ForceY,ForceZ\n"));
314ca69d878SAdeleke O. Bankole       user->app_ctx->wall_forces.header_written = PETSC_TRUE;
315ca69d878SAdeleke O. Bankole     }
316ca69d878SAdeleke O. Bankole     for (PetscInt w = 0; w < num_wall; w++) {
317ca69d878SAdeleke O. Bankole       PetscInt wall = walls[w];
318ca69d878SAdeleke O. Bankole       if (format == PETSC_VIEWER_ASCII_CSV) {
319ca69d878SAdeleke O. Bankole         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT ",%g,%" PetscInt_FMT ",%g,%g,%g\n", step_no, time, wall,
320ca69d878SAdeleke O. Bankole                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
321ca69d878SAdeleke O. Bankole 
322ca69d878SAdeleke O. Bankole       } else {
323ca69d878SAdeleke O. Bankole         PetscCall(PetscViewerASCIIPrintf(viewer, "Wall %" PetscInt_FMT " Forces: Force_x = %12g, Force_y = %12g, Force_z = %12g\n", wall,
324ca69d878SAdeleke O. Bankole                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
325ca69d878SAdeleke O. Bankole       }
326ca69d878SAdeleke O. Bankole     }
327ca69d878SAdeleke O. Bankole   }
328ca69d878SAdeleke O. Bankole   PetscCall(PetscFree(reaction_force));
329ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
330ca69d878SAdeleke O. Bankole }
331ca69d878SAdeleke O. Bankole 
332f14660a4SJames Wright // User provided TS Monitor
3332b730f8bSJeremy L Thompson PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
334f14660a4SJames Wright   User user = ctx;
335f14660a4SJames Wright 
336f17d818dSJames Wright   PetscFunctionBeginUser;
33737cbb16aSJed Brown   // Print every 'checkpoint_interval' steps
338894de27cSJames Wright   if (user->app_ctx->checkpoint_interval <= 0 || step_no % user->app_ctx->checkpoint_interval != 0 ||
33949aac155SJeremy L Thompson       (user->app_ctx->cont_steps == step_no && step_no != 0)) {
340ee4ca9cbSJames Wright     PetscFunctionReturn(PETSC_SUCCESS);
34149aac155SJeremy L Thompson   }
342f14660a4SJames Wright 
343f14660a4SJames Wright   PetscCall(WriteOutput(user, Q, step_no, time));
344ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
34577841947SLeila Ghaffari }
34677841947SLeila Ghaffari 
34777841947SLeila Ghaffari // TS: Create, setup, and solve
3482b730f8bSJeremy L Thompson PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, Vec *Q, PetscScalar *f_time, TS *ts) {
34977841947SLeila Ghaffari   MPI_Comm    comm = user->comm;
35077841947SLeila Ghaffari   TSAdapt     adapt;
35177841947SLeila Ghaffari   PetscScalar final_time;
35277841947SLeila Ghaffari 
353f17d818dSJames Wright   PetscFunctionBeginUser;
3542b730f8bSJeremy L Thompson   PetscCall(TSCreate(comm, ts));
3552b730f8bSJeremy L Thompson   PetscCall(TSSetDM(*ts, dm));
3561a7db67cSJames Wright   PetscCall(TSSetApplicationContext(*ts, user));
35777841947SLeila Ghaffari   if (phys->implicit) {
3582b730f8bSJeremy L Thompson     PetscCall(TSSetType(*ts, TSBDF));
35977841947SLeila Ghaffari     if (user->op_ifunction) {
3602b730f8bSJeremy L Thompson       PetscCall(TSSetIFunction(*ts, NULL, IFunction_NS, &user));
36177841947SLeila Ghaffari     } else {  // Implicit integrators can fall back to using an RHSFunction
3622b730f8bSJeremy L Thompson       PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
36377841947SLeila Ghaffari     }
36491c97f41SJames Wright     if (user->mat_ijacobian) {
3652b730f8bSJeremy L Thompson       PetscCall(DMTSSetIJacobian(dm, FormIJacobian_NS, &user));
366e334ad8fSJed Brown     }
36777841947SLeila Ghaffari   } else {
3689ad5e8e4SJames Wright     PetscCheck(user->op_rhs_ctx, comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction");
3692b730f8bSJeremy L Thompson     PetscCall(TSSetType(*ts, TSRK));
3702b730f8bSJeremy L Thompson     PetscCall(TSRKSetType(*ts, TSRK5F));
3712b730f8bSJeremy L Thompson     PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
37277841947SLeila Ghaffari   }
3732b730f8bSJeremy L Thompson   PetscCall(TSSetMaxTime(*ts, 500. * user->units->second));
3742b730f8bSJeremy L Thompson   PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER));
37551b00d91SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(TSSetErrorIfStepFails(*ts, PETSC_FALSE));
3762b730f8bSJeremy L Thompson   PetscCall(TSSetTimeStep(*ts, 1.e-2 * user->units->second));
3772b730f8bSJeremy L Thompson   PetscCall(TSGetAdapt(*ts, &adapt));
3782b730f8bSJeremy L Thompson   PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * user->units->second, 1.e2 * user->units->second));
3792b730f8bSJeremy L Thompson   PetscCall(TSSetFromOptions(*ts));
38091c97f41SJames Wright   if (user->mat_ijacobian) {
38191c97f41SJames Wright     if (app_ctx->amat_type && !strcmp(app_ctx->amat_type, MATSHELL)) {
38291c97f41SJames Wright       SNES snes;
38391c97f41SJames Wright       KSP  ksp;
3847f2a9303SJames Wright       Mat  Pmat, Amat;
38591c97f41SJames Wright 
38691c97f41SJames Wright       PetscCall(TSGetSNES(*ts, &snes));
38791c97f41SJames Wright       PetscCall(SNESGetKSP(snes, &ksp));
3887f2a9303SJames Wright       PetscCall(CreateSolveOperatorsFromMatCeed(ksp, user->mat_ijacobian, PETSC_FALSE, &Amat, &Pmat));
38991c97f41SJames Wright       PetscCall(TSSetIJacobian(*ts, user->mat_ijacobian, Pmat, NULL, NULL));
3907f2a9303SJames Wright       PetscCall(MatDestroy(&Amat));
3917f2a9303SJames Wright       PetscCall(MatDestroy(&Pmat));
39291c97f41SJames Wright     }
39391c97f41SJames Wright   }
39415637395SJames Wright   user->time_bc_set = -1.0;   // require all BCs be updated
395d55646a4SJames Wright   if (app_ctx->cont_steps) {  // continue from previous timestep data
39677841947SLeila Ghaffari     PetscInt    count;
39777841947SLeila Ghaffari     PetscViewer viewer;
3982b730f8bSJeremy L Thompson 
3994de8550aSJed Brown     if (app_ctx->cont_time <= 0) {  // Legacy files did not include step number and time
4002b730f8bSJeremy L Thompson       PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_time_file, FILE_MODE_READ, &viewer));
4014de8550aSJed Brown       PetscCall(PetscViewerBinaryRead(viewer, &app_ctx->cont_time, 1, &count, PETSC_REAL));
4022b730f8bSJeremy L Thompson       PetscCall(PetscViewerDestroy(&viewer));
4034de8550aSJed Brown       PetscCheck(app_ctx->cont_steps != -1, comm, PETSC_ERR_ARG_INCOMP,
4044de8550aSJed Brown                  "-continue step number not specified, but checkpoint file does not contain a step number (likely written by older code version)");
4054de8550aSJed Brown     }
4064de8550aSJed Brown     PetscCall(TSSetTime(*ts, app_ctx->cont_time * user->units->second));
407d8e0aecdSJed Brown     PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps));
40877841947SLeila Ghaffari   }
4098fb33541SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) {
4102b730f8bSJeremy L Thompson     PetscCall(TSMonitorSet(*ts, TSMonitor_NS, user, NULL));
4118fb33541SJames Wright   }
412ca69d878SAdeleke O. Bankole   if (app_ctx->wall_forces.viewer) {
413ca69d878SAdeleke O. Bankole     PetscCall(TSMonitorSet(*ts, TSMonitor_WallForce, user, NULL));
414ca69d878SAdeleke O. Bankole   }
415b7d66439SJames Wright   if (app_ctx->turb_spanstats_enable) {
416f5452247SJames Wright     PetscCall(TSMonitorSet(*ts, TSMonitor_TurbulenceStatistics, user, NULL));
417495b9769SJames Wright     CeedScalar previous_time = app_ctx->cont_time * user->units->second;
418a424bcd0SJames Wright     PetscCallCeed(user->ceed,
419a424bcd0SJames Wright                   CeedOperatorSetContextDouble(user->spanstats.op_stats_collect_ctx->op, user->spanstats.previous_time_label, &previous_time));
420a175e481SJames Wright   }
4214e9802d1SJames Wright   if (app_ctx->diff_filter_monitor) PetscCall(TSMonitorSet(*ts, TSMonitor_DifferentialFilter, user, NULL));
42277841947SLeila Ghaffari 
42328134cfaSJames Wright   if (app_ctx->sgs_train_enable) {
42428134cfaSJames Wright     PetscCall(TSMonitorSet(*ts, TSMonitor_SGS_DD_Training, user, NULL));
42528134cfaSJames Wright     PetscCall(TSSetPostStep(*ts, TSPostStep_SGS_DD_Training));
42628134cfaSJames Wright   }
42777841947SLeila Ghaffari   // Solve
428d8e0aecdSJed Brown   PetscReal start_time;
429d8e0aecdSJed Brown   PetscInt  start_step;
4302b730f8bSJeremy L Thompson   PetscCall(TSGetTime(*ts, &start_time));
431d8e0aecdSJed Brown   PetscCall(TSGetStepNumber(*ts, &start_step));
4327b39487dSJeremy L Thompson 
433711a423aSJed Brown   PetscCall(PetscLogDefaultBegin());  // So we can use PetscLogStageGetPerfInfo without -log_view
4347b39487dSJeremy L Thompson   PetscPreLoadBegin(PETSC_FALSE, "Fluids Solve");
4357b39487dSJeremy L Thompson   PetscCall(TSSetTime(*ts, start_time));
436d8e0aecdSJed Brown   PetscCall(TSSetStepNumber(*ts, start_step));
4377b39487dSJeremy L Thompson   if (PetscPreLoadingOn) {
4387b39487dSJeremy L Thompson     // LCOV_EXCL_START
4397b39487dSJeremy L Thompson     SNES      snes;
4407b39487dSJeremy L Thompson     Vec       Q_preload;
4417b39487dSJeremy L Thompson     PetscReal rtol;
4427b39487dSJeremy L Thompson     PetscCall(VecDuplicate(*Q, &Q_preload));
4437b39487dSJeremy L Thompson     PetscCall(VecCopy(*Q, Q_preload));
4447b39487dSJeremy L Thompson     PetscCall(TSGetSNES(*ts, &snes));
4457b39487dSJeremy L Thompson     PetscCall(SNESGetTolerances(snes, NULL, &rtol, NULL, NULL, NULL));
4462b730f8bSJeremy L Thompson     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, .99, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
447fc4c3d69SJames Wright     PetscCall(TSSetSolution(*ts, Q_preload));
4487b39487dSJeremy L Thompson     PetscCall(TSStep(*ts));
4492b730f8bSJeremy L Thompson     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, rtol, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
4507b39487dSJeremy L Thompson     PetscCall(VecDestroy(&Q_preload));
4517b39487dSJeremy L Thompson     // LCOV_EXCL_STOP
4527b39487dSJeremy L Thompson   } else {
4532b730f8bSJeremy L Thompson     PetscCall(PetscBarrier((PetscObject)*ts));
4542b730f8bSJeremy L Thompson     PetscCall(TSSolve(*ts, *Q));
4557b39487dSJeremy L Thompson   }
4567b39487dSJeremy L Thompson   PetscPreLoadEnd();
4577b39487dSJeremy L Thompson 
4587b39487dSJeremy L Thompson   PetscCall(TSGetSolveTime(*ts, &final_time));
45977841947SLeila Ghaffari   *f_time = final_time;
4607b39487dSJeremy L Thompson 
4618fb33541SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) {
462f14660a4SJames Wright     PetscInt step_no;
463f14660a4SJames Wright     PetscCall(TSGetStepNumber(*ts, &step_no));
464a175e481SJames Wright     if (user->app_ctx->checkpoint_interval > 0 || user->app_ctx->checkpoint_interval == -1) {
465f14660a4SJames Wright       PetscCall(WriteOutput(user, *Q, step_no, final_time));
466f14660a4SJames Wright     }
467f14660a4SJames Wright 
46875d1798cSJames Wright     PetscLogStage      stage_id;
469711a423aSJed Brown     PetscEventPerfInfo stage_perf;
4707b39487dSJeremy L Thompson 
4717b39487dSJeremy L Thompson     PetscCall(PetscLogStageGetId("Fluids Solve", &stage_id));
472711a423aSJed Brown     PetscCall(PetscLogStageGetPerfInfo(stage_id, &stage_perf));
473711a423aSJed Brown     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time taken for solution (sec): %g\n", stage_perf.time));
47477841947SLeila Ghaffari   }
475ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
47677841947SLeila Ghaffari }
477