xref: /honee/src/setupts.c (revision 62ab3c0b513036ca23e682b988b54f2e08122db7)
1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3a515125bSLeila Ghaffari 
4a515125bSLeila Ghaffari /// @file
5a515125bSLeila Ghaffari /// Time-stepping functions for Navier-Stokes example using PETSc
6a515125bSLeila Ghaffari 
7e419654dSJeremy L Thompson #include <ceed.h>
8e419654dSJeremy L Thompson #include <petscdmplex.h>
9e419654dSJeremy L Thompson #include <petscts.h>
10e419654dSJeremy L Thompson 
11149fb536SJames Wright #include <navierstokes.h>
12c5e9980aSAdeleke O. Bankole #include "../qfunctions/newtonian_state.h"
13a515125bSLeila Ghaffari 
144cb09fddSJames Wright // @brief Insert Boundary values if it's a new time
15c996854bSJames Wright PetscErrorCode UpdateBoundaryValues(User user, Vec Q_loc, PetscReal t) {
16c996854bSJames Wright   PetscFunctionBeginUser;
17c996854bSJames Wright   if (user->time_bc_set != t) {
18c996854bSJames Wright     PetscCall(DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL));
19c996854bSJames Wright     user->time_bc_set = t;
20c996854bSJames Wright   }
21d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
22c996854bSJames Wright }
23c996854bSJames Wright 
24a515125bSLeila Ghaffari // RHS (Explicit time-stepper) function setup
25a515125bSLeila Ghaffari //   This is the RHS of the ODE, given as u_t = G(t,u)
26a515125bSLeila Ghaffari //   This function takes in a state vector Q and writes into G
27a515125bSLeila Ghaffari PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) {
28a515125bSLeila Ghaffari   User         user = *(User *)user_data;
29701e5830SJames Wright   Ceed         ceed = user->ceed;
30fd969b44SJames Wright   PetscScalar  dt;
31da5fe0e4SJames Wright   Vec          Q_loc = user->Q_loc;
32a78efa86SJames Wright   PetscMemType q_mem_type;
33a515125bSLeila Ghaffari 
3406f41313SJames Wright   PetscFunctionBeginUser;
35e2f84137SJeremy L Thompson   // Update time dependent data
36c996854bSJames Wright   PetscCall(UpdateBoundaryValues(user, Q_loc, t));
37701e5830SJames Wright   if (user->phys->solution_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_rhs_ctx->op, user->phys->solution_time_label, &t));
382b916ea7SJeremy L Thompson   PetscCall(TSGetTimeStep(ts, &dt));
39701e5830SJames Wright   if (user->phys->timestep_size_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_rhs_ctx->op, user->phys->timestep_size_label, &dt));
40a515125bSLeila Ghaffari 
41da5fe0e4SJames Wright   PetscCall(ApplyCeedOperatorGlobalToGlobal(Q, G, user->op_rhs_ctx));
42a515125bSLeila Ghaffari 
43a78efa86SJames Wright   PetscCall(VecReadPetscToCeed(Q_loc, &q_mem_type, user->q_ceed));
44a78efa86SJames Wright 
45a78efa86SJames Wright   // Inverse of the mass matrix
46f8e2d240SJames Wright   PetscCall(KSPSolve(user->mass_ksp, G, G));
47a78efa86SJames Wright 
48a78efa86SJames Wright   PetscCall(VecReadCeedToPetsc(user->q_ceed, q_mem_type, Q_loc));
49d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
50a515125bSLeila Ghaffari }
51a515125bSLeila Ghaffari 
52c5e9980aSAdeleke O. Bankole // Surface forces function setup
53c5e9980aSAdeleke O. Bankole static PetscErrorCode Surface_Forces_NS(DM dm, Vec G_loc, PetscInt num_walls, const PetscInt walls[], PetscScalar *reaction_force) {
54c5e9980aSAdeleke O. Bankole   DMLabel            face_label;
55*62ab3c0bSJames Wright   const PetscScalar *g_array;
56*62ab3c0bSJames Wright   PetscInt           dim  = 3;
57*62ab3c0bSJames Wright   MPI_Comm           comm = PetscObjectComm((PetscObject)dm);
58*62ab3c0bSJames Wright   PetscSection       section;
59c5e9980aSAdeleke O. Bankole 
60c5e9980aSAdeleke O. Bankole   PetscFunctionBeginUser;
61c5e9980aSAdeleke O. Bankole   PetscCall(DMGetLabel(dm, "Face Sets", &face_label));
62*62ab3c0bSJames Wright   PetscCall(VecGetArrayRead(G_loc, &g_array));
63c5e9980aSAdeleke O. Bankole   for (PetscInt w = 0; w < num_walls; w++) {
64*62ab3c0bSJames Wright     const PetscInt wall = walls[w], *points;
65c5e9980aSAdeleke O. Bankole     IS             wall_is;
66*62ab3c0bSJames Wright     PetscInt       num_points, num_comp = 0;
67*62ab3c0bSJames Wright 
68c5e9980aSAdeleke O. Bankole     PetscCall(DMLabelGetStratumIS(face_label, wall, &wall_is));
69*62ab3c0bSJames Wright     if (!wall_is) continue;  // No wall points on this process, skip
70*62ab3c0bSJames Wright 
71*62ab3c0bSJames Wright     PetscCall(DMGetLocalSection(dm, &section));
72*62ab3c0bSJames Wright     PetscCall(PetscSectionGetFieldComponents(section, 0, &num_comp));
73c5e9980aSAdeleke O. Bankole     PetscCall(ISGetSize(wall_is, &num_points));
74c5e9980aSAdeleke O. Bankole     PetscCall(ISGetIndices(wall_is, &points));
75c5e9980aSAdeleke O. Bankole     for (PetscInt i = 0; i < num_points; i++) {
76c5e9980aSAdeleke O. Bankole       const PetscInt           p = points[i];
77c5e9980aSAdeleke O. Bankole       const StateConservative *r;
78*62ab3c0bSJames Wright       PetscInt                 dof;
79*62ab3c0bSJames Wright 
80*62ab3c0bSJames Wright       PetscCall(DMPlexPointLocalRead(dm, p, g_array, &r));
81*62ab3c0bSJames Wright       PetscCall(PetscSectionGetDof(section, p, &dof));
822004e3acSAdeleke O. Bankole       for (PetscInt node = 0; node < dof / num_comp; node++) {
83*62ab3c0bSJames Wright         for (PetscInt j = 0; j < dim; j++) {
842004e3acSAdeleke O. Bankole           reaction_force[w * dim + j] -= r[node].momentum[j];
852004e3acSAdeleke O. Bankole         }
86c5e9980aSAdeleke O. Bankole       }
87c5e9980aSAdeleke O. Bankole     }
88c5e9980aSAdeleke O. Bankole     PetscCall(ISRestoreIndices(wall_is, &points));
89c5e9980aSAdeleke O. Bankole     PetscCall(ISDestroy(&wall_is));
90c5e9980aSAdeleke O. Bankole   }
91c5e9980aSAdeleke O. Bankole   PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, reaction_force, dim * num_walls, MPIU_SCALAR, MPI_SUM, comm));
92*62ab3c0bSJames Wright   PetscCall(VecRestoreArrayRead(G_loc, &g_array));
93d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
94c5e9980aSAdeleke O. Bankole }
95c5e9980aSAdeleke O. Bankole 
96a515125bSLeila Ghaffari // Implicit time-stepper function setup
972b916ea7SJeremy L Thompson PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) {
98a515125bSLeila Ghaffari   User         user = *(User *)user_data;
99701e5830SJames Wright   Ceed         ceed = user->ceed;
100fd969b44SJames Wright   PetscScalar  dt;
101e2f84137SJeremy L Thompson   Vec          Q_loc = user->Q_loc, Q_dot_loc = user->Q_dot_loc, G_loc;
102a515125bSLeila Ghaffari   PetscMemType q_mem_type, q_dot_mem_type, g_mem_type;
103a515125bSLeila Ghaffari 
10406f41313SJames Wright   PetscFunctionBeginUser;
105e2f84137SJeremy L Thompson   // Get local vectors
106c5e9980aSAdeleke O. Bankole   PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
107e2f84137SJeremy L Thompson 
108e2f84137SJeremy L Thompson   // Update time dependent data
109c996854bSJames Wright   PetscCall(UpdateBoundaryValues(user, Q_loc, t));
110701e5830SJames Wright   if (user->phys->solution_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_ifunction, user->phys->solution_time_label, &t));
1112b916ea7SJeremy L Thompson   PetscCall(TSGetTimeStep(ts, &dt));
112701e5830SJames Wright   if (user->phys->timestep_size_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_ifunction, user->phys->timestep_size_label, &dt));
113a515125bSLeila Ghaffari 
114a515125bSLeila Ghaffari   // Global-to-local
11506108310SJames Wright   PetscCall(DMGlobalToLocalBegin(user->dm, Q, INSERT_VALUES, Q_loc));
11606108310SJames Wright   PetscCall(DMGlobalToLocalBegin(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
11706108310SJames Wright   PetscCall(DMGlobalToLocalEnd(user->dm, Q, INSERT_VALUES, Q_loc));
11806108310SJames Wright   PetscCall(DMGlobalToLocalEnd(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
119a515125bSLeila Ghaffari 
120a515125bSLeila Ghaffari   // Place PETSc vectors in CEED vectors
121a7dac1d5SJames Wright   PetscCall(VecReadPetscToCeed(Q_loc, &q_mem_type, user->q_ceed));
122a7dac1d5SJames Wright   PetscCall(VecReadPetscToCeed(Q_dot_loc, &q_dot_mem_type, user->q_dot_ceed));
123a7dac1d5SJames Wright   PetscCall(VecPetscToCeed(G_loc, &g_mem_type, user->g_ceed));
124a515125bSLeila Ghaffari 
125a515125bSLeila Ghaffari   // Apply CEED operator
1267eedc94cSJames Wright   PetscCall(PetscLogEventBegin(FLUIDS_CeedOperatorApply, Q, G, 0, 0));
1277eedc94cSJames Wright   PetscCall(PetscLogGpuTimeBegin());
128b4c37c5cSJames Wright   PetscCallCeed(user->ceed, CeedOperatorApply(user->op_ifunction, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE));
1297eedc94cSJames Wright   PetscCall(PetscLogGpuTimeEnd());
1307eedc94cSJames Wright   PetscCall(PetscLogEventEnd(FLUIDS_CeedOperatorApply, Q, G, 0, 0));
131a515125bSLeila Ghaffari 
132a515125bSLeila Ghaffari   // Restore vectors
133a7dac1d5SJames Wright   PetscCall(VecReadCeedToPetsc(user->q_ceed, q_mem_type, Q_loc));
134a7dac1d5SJames Wright   PetscCall(VecReadCeedToPetsc(user->q_dot_ceed, q_dot_mem_type, Q_dot_loc));
135a7dac1d5SJames Wright   PetscCall(VecCeedToPetsc(user->g_ceed, g_mem_type, G_loc));
136a515125bSLeila Ghaffari 
13701ab89c1SJames Wright   if (user->app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) {
138ad494f68SJames Wright     PetscCall(SgsDDApplyIFunction(user, Q_loc, G_loc));
13901ab89c1SJames Wright   }
1409c678832SJames Wright 
141a515125bSLeila Ghaffari   // Local-to-Global
1422b916ea7SJeremy L Thompson   PetscCall(VecZeroEntries(G));
1432b916ea7SJeremy L Thompson   PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G));
144a515125bSLeila Ghaffari 
145a515125bSLeila Ghaffari   // Restore vectors
146c5e9980aSAdeleke O. Bankole   PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
147d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
148a515125bSLeila Ghaffari }
149a515125bSLeila Ghaffari 
1502b916ea7SJeremy L Thompson PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, PetscReal shift, Mat J, Mat J_pre, void *user_data) {
151f0b65372SJed Brown   User      user = *(User *)user_data;
152ebfabadfSJames Wright   PetscBool J_is_matceed, J_is_mffd, J_pre_is_matceed, J_pre_is_mffd;
15306f41313SJames Wright 
154f0b65372SJed Brown   PetscFunctionBeginUser;
15504855949SJed Brown   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATMFFD, &J_is_mffd));
156ebfabadfSJames Wright   PetscCall(PetscObjectTypeCompare((PetscObject)J, MATCEED, &J_is_matceed));
157ebfabadfSJames Wright   PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATMFFD, &J_pre_is_mffd));
158ebfabadfSJames Wright   PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATCEED, &J_pre_is_matceed));
159ebfabadfSJames Wright 
16051bb547fSJames Wright   PetscCall(MatCeedSetContextReal(user->mat_ijacobian, "ijacobian time shift", shift));
161ebfabadfSJames Wright 
162ebfabadfSJames Wright   if (J_is_matceed || J_is_mffd) {
16304855949SJed Brown     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
16404855949SJed Brown     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
165ebfabadfSJames Wright   } else PetscCall(MatCeedAssembleCOO(user->mat_ijacobian, J));
166ebfabadfSJames Wright 
167ebfabadfSJames Wright   if (J_pre_is_matceed && J != J_pre) {
168ebfabadfSJames Wright     PetscCall(MatAssemblyBegin(J_pre, MAT_FINAL_ASSEMBLY));
169ebfabadfSJames Wright     PetscCall(MatAssemblyEnd(J_pre, MAT_FINAL_ASSEMBLY));
170ebfabadfSJames Wright   } else if (!J_pre_is_matceed && !J_pre_is_mffd && J != J_pre) {
171ebfabadfSJames Wright     PetscCall(MatCeedAssembleCOO(user->mat_ijacobian, J_pre));
17204855949SJed Brown   }
173d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
174f0b65372SJed Brown }
175f0b65372SJed Brown 
1762b916ea7SJeremy L Thompson PetscErrorCode WriteOutput(User user, Vec Q, PetscInt step_no, PetscScalar time) {
177a515125bSLeila Ghaffari   Vec         Q_loc;
178a515125bSLeila Ghaffari   char        file_path[PETSC_MAX_PATH_LEN];
179a515125bSLeila Ghaffari   PetscViewer viewer;
180a515125bSLeila Ghaffari 
18106f41313SJames Wright   PetscFunctionBeginUser;
182852e5969SJed Brown   if (user->app_ctx->checkpoint_vtk) {
183a515125bSLeila Ghaffari     // Set up output
1847538d537SJames Wright     PetscCall(DMGetLocalVector(user->dm, &Q_loc));
1857538d537SJames Wright     PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec"));
1867538d537SJames Wright     PetscCall(VecZeroEntries(Q_loc));
1877538d537SJames Wright     PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
188a515125bSLeila Ghaffari 
189a515125bSLeila Ghaffari     // Output
190852e5969SJed Brown     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
1917538d537SJames Wright 
1922b916ea7SJeremy L Thompson     PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path, FILE_MODE_WRITE, &viewer));
1937538d537SJames Wright     PetscCall(VecView(Q_loc, viewer));
1947538d537SJames Wright     PetscCall(PetscViewerDestroy(&viewer));
195a515125bSLeila Ghaffari     if (user->dm_viz) {
196a515125bSLeila Ghaffari       Vec         Q_refined, Q_refined_loc;
197a515125bSLeila Ghaffari       char        file_path_refined[PETSC_MAX_PATH_LEN];
198a515125bSLeila Ghaffari       PetscViewer viewer_refined;
199a515125bSLeila Ghaffari 
2007538d537SJames Wright       PetscCall(DMGetGlobalVector(user->dm_viz, &Q_refined));
2017538d537SJames Wright       PetscCall(DMGetLocalVector(user->dm_viz, &Q_refined_loc));
2027538d537SJames Wright       PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined"));
2037538d537SJames Wright 
2047538d537SJames Wright       PetscCall(MatInterpolate(user->interp_viz, Q, Q_refined));
2057538d537SJames Wright       PetscCall(VecZeroEntries(Q_refined_loc));
2062b916ea7SJeremy L Thompson       PetscCall(DMGlobalToLocal(user->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc));
2077538d537SJames Wright 
208852e5969SJed Brown       PetscCall(
209852e5969SJed Brown           PetscSNPrintf(file_path_refined, sizeof file_path_refined, "%s/nsrefined-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
2107538d537SJames Wright 
2112b916ea7SJeremy L Thompson       PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined), file_path_refined, FILE_MODE_WRITE, &viewer_refined));
2127538d537SJames Wright       PetscCall(VecView(Q_refined_loc, viewer_refined));
2137538d537SJames Wright       PetscCall(DMRestoreLocalVector(user->dm_viz, &Q_refined_loc));
2147538d537SJames Wright       PetscCall(DMRestoreGlobalVector(user->dm_viz, &Q_refined));
2157538d537SJames Wright       PetscCall(PetscViewerDestroy(&viewer_refined));
216a515125bSLeila Ghaffari     }
2177538d537SJames Wright     PetscCall(DMRestoreLocalVector(user->dm, &Q_loc));
218852e5969SJed Brown   }
219a515125bSLeila Ghaffari 
220a515125bSLeila Ghaffari   // Save data in a binary file for continuation of simulations
22191a36801SJames Wright   if (user->app_ctx->add_stepnum2bin) {
222852e5969SJed Brown     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", user->app_ctx->output_dir, step_no));
22391a36801SJames Wright   } else {
2242b916ea7SJeremy L Thompson     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", user->app_ctx->output_dir));
22591a36801SJames Wright   }
2262b916ea7SJeremy L Thompson   PetscCall(PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, &viewer));
2277538d537SJames Wright 
2289293eaa1SJed Brown   time /= user->units->second;  // Dimensionalize time back
229b237916aSJames Wright   PetscCall(HoneeWriteBinaryVec(viewer, Q, time, step_no));
2307538d537SJames Wright   PetscCall(PetscViewerDestroy(&viewer));
231d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2327538d537SJames Wright }
2337538d537SJames Wright 
234c5e9980aSAdeleke O. Bankole // CSV Monitor
235c5e9980aSAdeleke O. Bankole PetscErrorCode TSMonitor_WallForce(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
236c5e9980aSAdeleke O. Bankole   User              user = ctx;
237c5e9980aSAdeleke O. Bankole   Vec               G_loc;
238c5e9980aSAdeleke O. Bankole   PetscInt          num_wall = user->app_ctx->wall_forces.num_wall, dim = 3;
239c5e9980aSAdeleke O. Bankole   const PetscInt   *walls  = user->app_ctx->wall_forces.walls;
240c5e9980aSAdeleke O. Bankole   PetscViewer       viewer = user->app_ctx->wall_forces.viewer;
241c5e9980aSAdeleke O. Bankole   PetscViewerFormat format = user->app_ctx->wall_forces.viewer_format;
242c5e9980aSAdeleke O. Bankole   PetscScalar      *reaction_force;
243*62ab3c0bSJames Wright   PetscBool         is_ascii;
244c5e9980aSAdeleke O. Bankole 
245c5e9980aSAdeleke O. Bankole   PetscFunctionBeginUser;
246d949ddfcSJames Wright   if (!viewer) PetscFunctionReturn(PETSC_SUCCESS);
247c5e9980aSAdeleke O. Bankole   PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
248*62ab3c0bSJames Wright   PetscCall(PetscCalloc1(num_wall * dim, &reaction_force));
249c5e9980aSAdeleke O. Bankole   PetscCall(Surface_Forces_NS(user->dm, G_loc, num_wall, walls, reaction_force));
250c5e9980aSAdeleke O. Bankole   PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
251c5e9980aSAdeleke O. Bankole 
252*62ab3c0bSJames Wright   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &is_ascii));
253c5e9980aSAdeleke O. Bankole 
254*62ab3c0bSJames Wright   if (is_ascii) {
255c5e9980aSAdeleke O. Bankole     if (format == PETSC_VIEWER_ASCII_CSV && !user->app_ctx->wall_forces.header_written) {
256c5e9980aSAdeleke O. Bankole       PetscCall(PetscViewerASCIIPrintf(viewer, "Step,Time,Wall,ForceX,ForceY,ForceZ\n"));
257c5e9980aSAdeleke O. Bankole       user->app_ctx->wall_forces.header_written = PETSC_TRUE;
258c5e9980aSAdeleke O. Bankole     }
259c5e9980aSAdeleke O. Bankole     for (PetscInt w = 0; w < num_wall; w++) {
260c5e9980aSAdeleke O. Bankole       PetscInt wall = walls[w];
261c5e9980aSAdeleke O. Bankole       if (format == PETSC_VIEWER_ASCII_CSV) {
262c5e9980aSAdeleke O. Bankole         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT ",%g,%" PetscInt_FMT ",%g,%g,%g\n", step_no, time, wall,
263c5e9980aSAdeleke O. Bankole                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
264c5e9980aSAdeleke O. Bankole 
265c5e9980aSAdeleke O. Bankole       } else {
266c5e9980aSAdeleke O. Bankole         PetscCall(PetscViewerASCIIPrintf(viewer, "Wall %" PetscInt_FMT " Forces: Force_x = %12g, Force_y = %12g, Force_z = %12g\n", wall,
267c5e9980aSAdeleke O. Bankole                                          reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
268c5e9980aSAdeleke O. Bankole       }
269c5e9980aSAdeleke O. Bankole     }
270c5e9980aSAdeleke O. Bankole   }
271c5e9980aSAdeleke O. Bankole   PetscCall(PetscFree(reaction_force));
272d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
273c5e9980aSAdeleke O. Bankole }
274c5e9980aSAdeleke O. Bankole 
2757538d537SJames Wright // User provided TS Monitor
2762b916ea7SJeremy L Thompson PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
2777538d537SJames Wright   User user = ctx;
2787538d537SJames Wright 
27906f41313SJames Wright   PetscFunctionBeginUser;
280852e5969SJed Brown   // Print every 'checkpoint_interval' steps
281c539088bSJames Wright   if (user->app_ctx->checkpoint_interval <= 0 || step_no % user->app_ctx->checkpoint_interval != 0 ||
282e419654dSJeremy L Thompson       (user->app_ctx->cont_steps == step_no && step_no != 0)) {
283d949ddfcSJames Wright     PetscFunctionReturn(PETSC_SUCCESS);
284e419654dSJeremy L Thompson   }
2857538d537SJames Wright 
2867538d537SJames Wright   PetscCall(WriteOutput(user, Q, step_no, time));
287d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
288a515125bSLeila Ghaffari }
289a515125bSLeila Ghaffari 
290a515125bSLeila Ghaffari // TS: Create, setup, and solve
291e539bbbaSJames Wright PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, ProblemData problem, Vec *Q, PetscScalar *f_time, TS *ts) {
292a515125bSLeila Ghaffari   MPI_Comm    comm = user->comm;
293a515125bSLeila Ghaffari   TSAdapt     adapt;
294a515125bSLeila Ghaffari   PetscScalar final_time;
295a515125bSLeila Ghaffari 
29606f41313SJames Wright   PetscFunctionBeginUser;
2972b916ea7SJeremy L Thompson   PetscCall(TSCreate(comm, ts));
2982b916ea7SJeremy L Thompson   PetscCall(TSSetDM(*ts, dm));
299632a41e1SJames Wright   PetscCall(TSSetApplicationContext(*ts, user));
300a515125bSLeila Ghaffari   if (phys->implicit) {
3012b916ea7SJeremy L Thompson     PetscCall(TSSetType(*ts, TSBDF));
302a515125bSLeila Ghaffari     if (user->op_ifunction) {
3032b916ea7SJeremy L Thompson       PetscCall(TSSetIFunction(*ts, NULL, IFunction_NS, &user));
304a515125bSLeila Ghaffari     } else {  // Implicit integrators can fall back to using an RHSFunction
3052b916ea7SJeremy L Thompson       PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
306a515125bSLeila Ghaffari     }
307ebfabadfSJames Wright     if (user->mat_ijacobian) {
3082b916ea7SJeremy L Thompson       PetscCall(DMTSSetIJacobian(dm, FormIJacobian_NS, &user));
309f0b65372SJed Brown     }
310a515125bSLeila Ghaffari   } else {
311da5fe0e4SJames Wright     PetscCheck(user->op_rhs_ctx, comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction");
3122b916ea7SJeremy L Thompson     PetscCall(TSSetType(*ts, TSRK));
3132b916ea7SJeremy L Thompson     PetscCall(TSRKSetType(*ts, TSRK5F));
3142b916ea7SJeremy L Thompson     PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
315a515125bSLeila Ghaffari   }
3162b916ea7SJeremy L Thompson   PetscCall(TSSetMaxTime(*ts, 500. * user->units->second));
3172b916ea7SJeremy L Thompson   PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER));
31822387d3aSJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(TSSetErrorIfStepFails(*ts, PETSC_FALSE));
3192b916ea7SJeremy L Thompson   PetscCall(TSSetTimeStep(*ts, 1.e-2 * user->units->second));
3202b916ea7SJeremy L Thompson   PetscCall(TSGetAdapt(*ts, &adapt));
3212b916ea7SJeremy L Thompson   PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * user->units->second, 1.e2 * user->units->second));
3222b916ea7SJeremy L Thompson   PetscCall(TSSetFromOptions(*ts));
323ebfabadfSJames Wright   if (user->mat_ijacobian) {
324ebfabadfSJames Wright     if (app_ctx->amat_type && !strcmp(app_ctx->amat_type, MATSHELL)) {
325ebfabadfSJames Wright       SNES snes;
326ebfabadfSJames Wright       KSP  ksp;
3273170c09fSJames Wright       Mat  Pmat, Amat;
328ebfabadfSJames Wright 
329ebfabadfSJames Wright       PetscCall(TSGetSNES(*ts, &snes));
330ebfabadfSJames Wright       PetscCall(SNESGetKSP(snes, &ksp));
3313170c09fSJames Wright       PetscCall(CreateSolveOperatorsFromMatCeed(ksp, user->mat_ijacobian, PETSC_FALSE, &Amat, &Pmat));
332ebfabadfSJames Wright       PetscCall(TSSetIJacobian(*ts, user->mat_ijacobian, Pmat, NULL, NULL));
3333170c09fSJames Wright       PetscCall(MatDestroy(&Amat));
3343170c09fSJames Wright       PetscCall(MatDestroy(&Pmat));
335ebfabadfSJames Wright     }
336ebfabadfSJames Wright   }
33791f639d2SJames Wright   user->time_bc_set = -1.0;   // require all BCs be updated
338c26b555cSJames Wright   if (app_ctx->cont_steps) {  // continue from previous timestep data
339a515125bSLeila Ghaffari     PetscInt    count;
340a515125bSLeila Ghaffari     PetscViewer viewer;
3412b916ea7SJeremy L Thompson 
3429293eaa1SJed Brown     if (app_ctx->cont_time <= 0) {  // Legacy files did not include step number and time
3432b916ea7SJeremy L Thompson       PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_time_file, FILE_MODE_READ, &viewer));
3449293eaa1SJed Brown       PetscCall(PetscViewerBinaryRead(viewer, &app_ctx->cont_time, 1, &count, PETSC_REAL));
3452b916ea7SJeremy L Thompson       PetscCall(PetscViewerDestroy(&viewer));
3469293eaa1SJed Brown       PetscCheck(app_ctx->cont_steps != -1, comm, PETSC_ERR_ARG_INCOMP,
3479293eaa1SJed Brown                  "-continue step number not specified, but checkpoint file does not contain a step number (likely written by older code version)");
3489293eaa1SJed Brown     }
3499293eaa1SJed Brown     PetscCall(TSSetTime(*ts, app_ctx->cont_time * user->units->second));
35074a6f4ddSJed Brown     PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps));
351a515125bSLeila Ghaffari   }
3520e1e9333SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) {
3532b916ea7SJeremy L Thompson     PetscCall(TSMonitorSet(*ts, TSMonitor_NS, user, NULL));
3540e1e9333SJames Wright   }
355c5e9980aSAdeleke O. Bankole   if (app_ctx->wall_forces.viewer) {
356c5e9980aSAdeleke O. Bankole     PetscCall(TSMonitorSet(*ts, TSMonitor_WallForce, user, NULL));
357c5e9980aSAdeleke O. Bankole   }
358c931fa59SJames Wright   if (app_ctx->turb_spanstats_enable) {
35991933550SJames Wright     PetscCall(TSMonitorSet(*ts, TSMonitor_TurbulenceStatistics, user, NULL));
360b8daee98SJames Wright     CeedScalar previous_time = app_ctx->cont_time * user->units->second;
361b4c37c5cSJames Wright     PetscCallCeed(user->ceed,
362b4c37c5cSJames Wright                   CeedOperatorSetContextDouble(user->spanstats.op_stats_collect_ctx->op, user->spanstats.previous_time_label, &previous_time));
363b0488d1fSJames Wright   }
36488b07121SJames Wright   if (app_ctx->diff_filter_monitor) PetscCall(TSMonitorSet(*ts, TSMonitor_DifferentialFilter, user, NULL));
365a515125bSLeila Ghaffari 
3661c17f66aSJames Wright   if (app_ctx->sgs_train_enable) {
3671c17f66aSJames Wright     PetscCall(TSMonitorSet(*ts, TSMonitor_SGS_DD_Training, user, NULL));
3681c17f66aSJames Wright     PetscCall(TSSetPostStep(*ts, TSPostStep_SGS_DD_Training));
3691c17f66aSJames Wright   }
370e539bbbaSJames Wright 
3713678fae3SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(PrintRunInfo(user, user->phys, problem, *ts));
372a515125bSLeila Ghaffari   // Solve
37374a6f4ddSJed Brown   PetscReal start_time;
37474a6f4ddSJed Brown   PetscInt  start_step;
3752b916ea7SJeremy L Thompson   PetscCall(TSGetTime(*ts, &start_time));
37674a6f4ddSJed Brown   PetscCall(TSGetStepNumber(*ts, &start_step));
37791982731SJeremy L Thompson 
378df4304b5SJed Brown   PetscCall(PetscLogDefaultBegin());  // So we can use PetscLogStageGetPerfInfo without -log_view
37991982731SJeremy L Thompson   PetscPreLoadBegin(PETSC_FALSE, "Fluids Solve");
38091982731SJeremy L Thompson   PetscCall(TSSetTime(*ts, start_time));
38174a6f4ddSJed Brown   PetscCall(TSSetStepNumber(*ts, start_step));
38291982731SJeremy L Thompson   if (PetscPreLoadingOn) {
38391982731SJeremy L Thompson     // LCOV_EXCL_START
38491982731SJeremy L Thompson     SNES      snes;
38591982731SJeremy L Thompson     Vec       Q_preload;
38691982731SJeremy L Thompson     PetscReal rtol;
38791982731SJeremy L Thompson     PetscCall(VecDuplicate(*Q, &Q_preload));
38891982731SJeremy L Thompson     PetscCall(VecCopy(*Q, Q_preload));
38991982731SJeremy L Thompson     PetscCall(TSGetSNES(*ts, &snes));
39091982731SJeremy L Thompson     PetscCall(SNESGetTolerances(snes, NULL, &rtol, NULL, NULL, NULL));
3912b916ea7SJeremy L Thompson     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, .99, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
39222c1b34eSJames Wright     PetscCall(TSSetSolution(*ts, Q_preload));
39391982731SJeremy L Thompson     PetscCall(TSStep(*ts));
3942b916ea7SJeremy L Thompson     PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, rtol, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
39591982731SJeremy L Thompson     PetscCall(VecDestroy(&Q_preload));
39691982731SJeremy L Thompson     // LCOV_EXCL_STOP
39791982731SJeremy L Thompson   } else {
3982b916ea7SJeremy L Thompson     PetscCall(PetscBarrier((PetscObject)*ts));
3992b916ea7SJeremy L Thompson     PetscCall(TSSolve(*ts, *Q));
40091982731SJeremy L Thompson   }
40191982731SJeremy L Thompson   PetscPreLoadEnd();
40291982731SJeremy L Thompson 
40391982731SJeremy L Thompson   PetscCall(TSGetSolveTime(*ts, &final_time));
404a515125bSLeila Ghaffari   *f_time = final_time;
40591982731SJeremy L Thompson 
4060e1e9333SJames Wright   if (app_ctx->test_type == TESTTYPE_NONE) {
4077538d537SJames Wright     PetscInt step_no;
4087538d537SJames Wright     PetscCall(TSGetStepNumber(*ts, &step_no));
409b0488d1fSJames Wright     if (user->app_ctx->checkpoint_interval > 0 || user->app_ctx->checkpoint_interval == -1) {
4107538d537SJames Wright       PetscCall(WriteOutput(user, *Q, step_no, final_time));
4117538d537SJames Wright     }
4127538d537SJames Wright 
4137eedc94cSJames Wright     PetscLogStage      stage_id;
414df4304b5SJed Brown     PetscEventPerfInfo stage_perf;
41591982731SJeremy L Thompson 
41691982731SJeremy L Thompson     PetscCall(PetscLogStageGetId("Fluids Solve", &stage_id));
417df4304b5SJed Brown     PetscCall(PetscLogStageGetPerfInfo(stage_id, &stage_perf));
418df4304b5SJed Brown     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time taken for solution (sec): %g\n", stage_perf.time));
419a515125bSLeila Ghaffari   }
420d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
421a515125bSLeila Ghaffari }
422