1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, 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
18a0b9a424SJames Wright // Insert Boundary values if it's a new time
UpdateBoundaryValues(User user,Vec Q_loc,PetscReal t)19a0b9a424SJames Wright PetscErrorCode UpdateBoundaryValues(User user, Vec Q_loc, PetscReal t) {
20a0b9a424SJames Wright PetscFunctionBeginUser;
21a0b9a424SJames Wright if (user->time_bc_set != t) {
22a0b9a424SJames Wright PetscCall(DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL));
23a0b9a424SJames Wright user->time_bc_set = t;
24a0b9a424SJames Wright }
25ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
26a0b9a424SJames Wright }
27a0b9a424SJames Wright
2877841947SLeila Ghaffari // RHS (Explicit time-stepper) function setup
2977841947SLeila Ghaffari // This is the RHS of the ODE, given as u_t = G(t,u)
3077841947SLeila Ghaffari // This function takes in a state vector Q and writes into G
RHS_NS(TS ts,PetscReal t,Vec Q,Vec G,void * user_data)3177841947SLeila Ghaffari PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) {
3277841947SLeila Ghaffari User user = *(User *)user_data;
337d4c6defSJames Wright Ceed ceed = user->ceed;
34c798d55aSJames Wright PetscScalar dt;
359ad5e8e4SJames Wright Vec Q_loc = user->Q_loc;
36b4e9a8f8SJames Wright PetscMemType q_mem_type;
3777841947SLeila Ghaffari
38f17d818dSJames Wright PetscFunctionBeginUser;
395e82a6e1SJeremy L Thompson // Update time dependent data
40a0b9a424SJames Wright PetscCall(UpdateBoundaryValues(user, Q_loc, t));
417d4c6defSJames Wright if (user->phys->solution_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_rhs_ctx->op, user->phys->solution_time_label, &t));
422b730f8bSJeremy L Thompson PetscCall(TSGetTimeStep(ts, &dt));
437d4c6defSJames Wright if (user->phys->timestep_size_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_rhs_ctx->op, user->phys->timestep_size_label, &dt));
4477841947SLeila Ghaffari
459ad5e8e4SJames Wright PetscCall(ApplyCeedOperatorGlobalToGlobal(Q, G, user->op_rhs_ctx));
4677841947SLeila Ghaffari
47b4e9a8f8SJames Wright PetscCall(VecReadPetscToCeed(Q_loc, &q_mem_type, user->q_ceed));
48b4e9a8f8SJames Wright
49b4e9a8f8SJames Wright // Inverse of the mass matrix
5004292f7dSJames Wright PetscCall(KSPSolve(user->mass_ksp, G, G));
51b4e9a8f8SJames Wright
52b4e9a8f8SJames Wright PetscCall(VecReadCeedToPetsc(user->q_ceed, q_mem_type, Q_loc));
53ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
5477841947SLeila Ghaffari }
5577841947SLeila Ghaffari
56ca69d878SAdeleke O. Bankole // Surface forces function setup
Surface_Forces_NS(DM dm,Vec G_loc,PetscInt num_walls,const PetscInt walls[],PetscScalar * reaction_force)57ca69d878SAdeleke O. Bankole static PetscErrorCode Surface_Forces_NS(DM dm, Vec G_loc, PetscInt num_walls, const PetscInt walls[], PetscScalar *reaction_force) {
58ca69d878SAdeleke O. Bankole DMLabel face_label;
59ca69d878SAdeleke O. Bankole const PetscScalar *g;
60d6734f85SAdeleke O. Bankole PetscInt dof, dim = 3;
61ca69d878SAdeleke O. Bankole MPI_Comm comm;
62d6734f85SAdeleke O. Bankole PetscSection s;
63ca69d878SAdeleke O. Bankole
64ca69d878SAdeleke O. Bankole PetscFunctionBeginUser;
65ca69d878SAdeleke O. Bankole PetscCall(PetscArrayzero(reaction_force, num_walls * dim));
66ca69d878SAdeleke O. Bankole PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
67ca69d878SAdeleke O. Bankole PetscCall(DMGetLabel(dm, "Face Sets", &face_label));
68ca69d878SAdeleke O. Bankole PetscCall(VecGetArrayRead(G_loc, &g));
69ca69d878SAdeleke O. Bankole for (PetscInt w = 0; w < num_walls; w++) {
70ca69d878SAdeleke O. Bankole const PetscInt wall = walls[w];
71ca69d878SAdeleke O. Bankole IS wall_is;
72d6734f85SAdeleke O. Bankole PetscCall(DMGetLocalSection(dm, &s));
73ca69d878SAdeleke O. Bankole PetscCall(DMLabelGetStratumIS(face_label, wall, &wall_is));
74ca69d878SAdeleke O. Bankole if (wall_is) { // There exist such points on this process
75ca69d878SAdeleke O. Bankole PetscInt num_points;
76d6734f85SAdeleke O. Bankole PetscInt num_comp = 0;
77ca69d878SAdeleke O. Bankole const PetscInt *points;
78d6734f85SAdeleke O. Bankole PetscCall(PetscSectionGetFieldComponents(s, 0, &num_comp));
79ca69d878SAdeleke O. Bankole PetscCall(ISGetSize(wall_is, &num_points));
80ca69d878SAdeleke O. Bankole PetscCall(ISGetIndices(wall_is, &points));
81ca69d878SAdeleke O. Bankole for (PetscInt i = 0; i < num_points; i++) {
82ca69d878SAdeleke O. Bankole const PetscInt p = points[i];
83ca69d878SAdeleke O. Bankole const StateConservative *r;
84ca69d878SAdeleke O. Bankole PetscCall(DMPlexPointLocalRead(dm, p, g, &r));
85d6734f85SAdeleke O. Bankole PetscCall(PetscSectionGetDof(s, p, &dof));
86d6734f85SAdeleke O. Bankole for (PetscInt node = 0; node < dof / num_comp; node++) {
87ca69d878SAdeleke O. Bankole for (PetscInt j = 0; j < 3; j++) {
88d6734f85SAdeleke O. Bankole reaction_force[w * dim + j] -= r[node].momentum[j];
89d6734f85SAdeleke O. Bankole }
90ca69d878SAdeleke O. Bankole }
91ca69d878SAdeleke O. Bankole }
92ca69d878SAdeleke O. Bankole PetscCall(ISRestoreIndices(wall_is, &points));
93ca69d878SAdeleke O. Bankole }
94ca69d878SAdeleke O. Bankole PetscCall(ISDestroy(&wall_is));
95ca69d878SAdeleke O. Bankole }
96ca69d878SAdeleke O. Bankole PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, reaction_force, dim * num_walls, MPIU_SCALAR, MPI_SUM, comm));
97ca69d878SAdeleke O. Bankole // Restore Vectors
98ca69d878SAdeleke O. Bankole PetscCall(VecRestoreArrayRead(G_loc, &g));
99ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
100ca69d878SAdeleke O. Bankole }
101ca69d878SAdeleke O. Bankole
10277841947SLeila Ghaffari // Implicit time-stepper function setup
IFunction_NS(TS ts,PetscReal t,Vec Q,Vec Q_dot,Vec G,void * user_data)1032b730f8bSJeremy L Thompson PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) {
10477841947SLeila Ghaffari User user = *(User *)user_data;
1057d4c6defSJames Wright Ceed ceed = user->ceed;
106c798d55aSJames Wright PetscScalar dt;
1075e82a6e1SJeremy L Thompson Vec Q_loc = user->Q_loc, Q_dot_loc = user->Q_dot_loc, G_loc;
10877841947SLeila Ghaffari PetscMemType q_mem_type, q_dot_mem_type, g_mem_type;
10977841947SLeila Ghaffari
110f17d818dSJames Wright PetscFunctionBeginUser;
1115e82a6e1SJeremy L Thompson // Get local vectors
112ca69d878SAdeleke O. Bankole PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
1135e82a6e1SJeremy L Thompson
1145e82a6e1SJeremy L Thompson // Update time dependent data
115a0b9a424SJames Wright PetscCall(UpdateBoundaryValues(user, Q_loc, t));
1167d4c6defSJames Wright if (user->phys->solution_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_ifunction, user->phys->solution_time_label, &t));
1172b730f8bSJeremy L Thompson PetscCall(TSGetTimeStep(ts, &dt));
1187d4c6defSJames Wright if (user->phys->timestep_size_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_ifunction, user->phys->timestep_size_label, &dt));
11977841947SLeila Ghaffari
12077841947SLeila Ghaffari // Global-to-local
121878eb0e7SJames Wright PetscCall(DMGlobalToLocalBegin(user->dm, Q, INSERT_VALUES, Q_loc));
122878eb0e7SJames Wright PetscCall(DMGlobalToLocalBegin(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
123878eb0e7SJames Wright PetscCall(DMGlobalToLocalEnd(user->dm, Q, INSERT_VALUES, Q_loc));
124878eb0e7SJames Wright PetscCall(DMGlobalToLocalEnd(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc));
12577841947SLeila Ghaffari
12677841947SLeila Ghaffari // Place PETSc vectors in CEED vectors
127d0593705SJames Wright PetscCall(VecReadPetscToCeed(Q_loc, &q_mem_type, user->q_ceed));
128d0593705SJames Wright PetscCall(VecReadPetscToCeed(Q_dot_loc, &q_dot_mem_type, user->q_dot_ceed));
129d0593705SJames Wright PetscCall(VecPetscToCeed(G_loc, &g_mem_type, user->g_ceed));
13077841947SLeila Ghaffari
13177841947SLeila Ghaffari // Apply CEED operator
13275d1798cSJames Wright PetscCall(PetscLogEventBegin(FLUIDS_CeedOperatorApply, Q, G, 0, 0));
13375d1798cSJames Wright PetscCall(PetscLogGpuTimeBegin());
134a424bcd0SJames Wright PetscCallCeed(user->ceed, CeedOperatorApply(user->op_ifunction, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE));
13575d1798cSJames Wright PetscCall(PetscLogGpuTimeEnd());
13675d1798cSJames Wright PetscCall(PetscLogEventEnd(FLUIDS_CeedOperatorApply, Q, G, 0, 0));
13777841947SLeila Ghaffari
13877841947SLeila Ghaffari // Restore vectors
139d0593705SJames Wright PetscCall(VecReadCeedToPetsc(user->q_ceed, q_mem_type, Q_loc));
140d0593705SJames Wright PetscCall(VecReadCeedToPetsc(user->q_dot_ceed, q_dot_mem_type, Q_dot_loc));
141d0593705SJames Wright PetscCall(VecCeedToPetsc(user->g_ceed, g_mem_type, G_loc));
14277841947SLeila Ghaffari
14377841947SLeila Ghaffari // Local-to-Global
1442b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(G));
1452b730f8bSJeremy L Thompson PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G));
14677841947SLeila Ghaffari
14777841947SLeila Ghaffari // Restore vectors
148ca69d878SAdeleke O. Bankole PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
149ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
15077841947SLeila Ghaffari }
15177841947SLeila Ghaffari
FormIJacobian_NS(TS ts,PetscReal t,Vec Q,Vec Q_dot,PetscReal shift,Mat J,Mat J_pre,void * user_data)1522b730f8bSJeremy L Thompson PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, PetscReal shift, Mat J, Mat J_pre, void *user_data) {
153e334ad8fSJed Brown User user = *(User *)user_data;
15491c97f41SJames Wright PetscBool J_is_matceed, J_is_mffd, J_pre_is_matceed, J_pre_is_mffd;
155f17d818dSJames Wright
156e334ad8fSJed Brown PetscFunctionBeginUser;
157d6bf345cSJed Brown PetscCall(PetscObjectTypeCompare((PetscObject)J, MATMFFD, &J_is_mffd));
15891c97f41SJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)J, MATCEED, &J_is_matceed));
15991c97f41SJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATMFFD, &J_pre_is_mffd));
16091c97f41SJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATCEED, &J_pre_is_matceed));
16191c97f41SJames Wright
162f965f5c6SJames Wright PetscCall(MatCeedSetContextReal(user->mat_ijacobian, "ijacobian time shift", shift));
16391c97f41SJames Wright
16491c97f41SJames Wright if (J_is_matceed || J_is_mffd) {
165d6bf345cSJed Brown PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
166d6bf345cSJed Brown PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
16791c97f41SJames Wright } else PetscCall(MatCeedAssembleCOO(user->mat_ijacobian, J));
16891c97f41SJames Wright
16991c97f41SJames Wright if (J_pre_is_matceed && J != J_pre) {
17091c97f41SJames Wright PetscCall(MatAssemblyBegin(J_pre, MAT_FINAL_ASSEMBLY));
17191c97f41SJames Wright PetscCall(MatAssemblyEnd(J_pre, MAT_FINAL_ASSEMBLY));
17291c97f41SJames Wright } else if (!J_pre_is_matceed && !J_pre_is_mffd && J != J_pre) {
17391c97f41SJames Wright PetscCall(MatCeedAssembleCOO(user->mat_ijacobian, J_pre));
174d6bf345cSJed Brown }
175ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
176e334ad8fSJed Brown }
177e334ad8fSJed Brown
WriteOutput(User user,Vec Q,PetscInt step_no,PetscScalar time)1782b730f8bSJeremy L Thompson PetscErrorCode WriteOutput(User user, Vec Q, PetscInt step_no, PetscScalar time) {
17977841947SLeila Ghaffari Vec Q_loc;
18077841947SLeila Ghaffari char file_path[PETSC_MAX_PATH_LEN];
18177841947SLeila Ghaffari PetscViewer viewer;
18277841947SLeila Ghaffari
183f17d818dSJames Wright PetscFunctionBeginUser;
18437cbb16aSJed Brown if (user->app_ctx->checkpoint_vtk) {
18577841947SLeila Ghaffari // Set up output
186f14660a4SJames Wright PetscCall(DMGetLocalVector(user->dm, &Q_loc));
187f14660a4SJames Wright PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec"));
188f14660a4SJames Wright PetscCall(VecZeroEntries(Q_loc));
189f14660a4SJames Wright PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
19077841947SLeila Ghaffari
19177841947SLeila Ghaffari // Output
19237cbb16aSJed Brown PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no));
193f14660a4SJames Wright
1942b730f8bSJeremy L Thompson PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path, FILE_MODE_WRITE, &viewer));
195f14660a4SJames Wright PetscCall(VecView(Q_loc, viewer));
196f14660a4SJames Wright PetscCall(PetscViewerDestroy(&viewer));
19777841947SLeila Ghaffari if (user->dm_viz) {
19877841947SLeila Ghaffari Vec Q_refined, Q_refined_loc;
19977841947SLeila Ghaffari char file_path_refined[PETSC_MAX_PATH_LEN];
20077841947SLeila Ghaffari PetscViewer viewer_refined;
20177841947SLeila Ghaffari
202f14660a4SJames Wright PetscCall(DMGetGlobalVector(user->dm_viz, &Q_refined));
203f14660a4SJames Wright PetscCall(DMGetLocalVector(user->dm_viz, &Q_refined_loc));
204f14660a4SJames Wright PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined"));
205f14660a4SJames Wright
206f14660a4SJames Wright PetscCall(MatInterpolate(user->interp_viz, Q, Q_refined));
207f14660a4SJames Wright PetscCall(VecZeroEntries(Q_refined_loc));
2082b730f8bSJeremy L Thompson PetscCall(DMGlobalToLocal(user->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc));
209f14660a4SJames Wright
2101a8516d0SJames Wright PetscCall(PetscSNPrintf(file_path_refined, sizeof file_path_refined, "%s/nsrefined-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir,
2111a8516d0SJames Wright step_no));
212f14660a4SJames Wright
2132b730f8bSJeremy L Thompson PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined), file_path_refined, FILE_MODE_WRITE, &viewer_refined));
214f14660a4SJames Wright PetscCall(VecView(Q_refined_loc, viewer_refined));
215f14660a4SJames Wright PetscCall(DMRestoreLocalVector(user->dm_viz, &Q_refined_loc));
216f14660a4SJames Wright PetscCall(DMRestoreGlobalVector(user->dm_viz, &Q_refined));
217f14660a4SJames Wright PetscCall(PetscViewerDestroy(&viewer_refined));
21877841947SLeila Ghaffari }
219f14660a4SJames Wright PetscCall(DMRestoreLocalVector(user->dm, &Q_loc));
22037cbb16aSJed Brown }
22177841947SLeila Ghaffari
22277841947SLeila Ghaffari // Save data in a binary file for continuation of simulations
22369293791SJames Wright if (user->app_ctx->add_stepnum2bin) {
22437cbb16aSJed Brown PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", user->app_ctx->output_dir, step_no));
22569293791SJames Wright } else {
2262b730f8bSJeremy L Thompson PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", user->app_ctx->output_dir));
22769293791SJames Wright }
2282b730f8bSJeremy L Thompson PetscCall(PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, &viewer));
229f14660a4SJames Wright
2300de6a49fSJames Wright PetscInt32 token = PetscDefined(USE_64BIT_INDICES) ? FLUIDS_FILE_TOKEN_64 : FLUIDS_FILE_TOKEN_32;
2310de6a49fSJames Wright PetscCall(PetscViewerBinaryWrite(viewer, &token, 1, PETSC_INT32));
2324de8550aSJed Brown PetscCall(PetscViewerBinaryWrite(viewer, &step_no, 1, PETSC_INT));
2334de8550aSJed Brown time /= user->units->second; // Dimensionalize time back
2344de8550aSJed Brown PetscCall(PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL));
235f14660a4SJames Wright PetscCall(VecView(Q, viewer));
236f14660a4SJames Wright PetscCall(PetscViewerDestroy(&viewer));
237ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
238f14660a4SJames Wright }
239f14660a4SJames Wright
240ca69d878SAdeleke O. Bankole // CSV Monitor
TSMonitor_WallForce(TS ts,PetscInt step_no,PetscReal time,Vec Q,void * ctx)241ca69d878SAdeleke O. Bankole PetscErrorCode TSMonitor_WallForce(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
242ca69d878SAdeleke O. Bankole User user = ctx;
243ca69d878SAdeleke O. Bankole Vec G_loc;
244ca69d878SAdeleke O. Bankole PetscInt num_wall = user->app_ctx->wall_forces.num_wall, dim = 3;
245ca69d878SAdeleke O. Bankole const PetscInt *walls = user->app_ctx->wall_forces.walls;
246ca69d878SAdeleke O. Bankole PetscViewer viewer = user->app_ctx->wall_forces.viewer;
247ca69d878SAdeleke O. Bankole PetscViewerFormat format = user->app_ctx->wall_forces.viewer_format;
248ca69d878SAdeleke O. Bankole PetscScalar *reaction_force;
249ca69d878SAdeleke O. Bankole PetscBool iascii;
250ca69d878SAdeleke O. Bankole
251ca69d878SAdeleke O. Bankole PetscFunctionBeginUser;
252ee4ca9cbSJames Wright if (!viewer) PetscFunctionReturn(PETSC_SUCCESS);
253ca69d878SAdeleke O. Bankole PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
254ca69d878SAdeleke O. Bankole PetscCall(PetscMalloc1(num_wall * dim, &reaction_force));
255ca69d878SAdeleke O. Bankole PetscCall(Surface_Forces_NS(user->dm, G_loc, num_wall, walls, reaction_force));
256ca69d878SAdeleke O. Bankole PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc));
257ca69d878SAdeleke O. Bankole
258ca69d878SAdeleke O. Bankole PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
259ca69d878SAdeleke O. Bankole
260ca69d878SAdeleke O. Bankole if (iascii) {
261ca69d878SAdeleke O. Bankole if (format == PETSC_VIEWER_ASCII_CSV && !user->app_ctx->wall_forces.header_written) {
262ca69d878SAdeleke O. Bankole PetscCall(PetscViewerASCIIPrintf(viewer, "Step,Time,Wall,ForceX,ForceY,ForceZ\n"));
263ca69d878SAdeleke O. Bankole user->app_ctx->wall_forces.header_written = PETSC_TRUE;
264ca69d878SAdeleke O. Bankole }
265ca69d878SAdeleke O. Bankole for (PetscInt w = 0; w < num_wall; w++) {
266ca69d878SAdeleke O. Bankole PetscInt wall = walls[w];
267ca69d878SAdeleke O. Bankole if (format == PETSC_VIEWER_ASCII_CSV) {
268ca69d878SAdeleke O. Bankole PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT ",%g,%" PetscInt_FMT ",%g,%g,%g\n", step_no, time, wall,
269ca69d878SAdeleke O. Bankole reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
270ca69d878SAdeleke O. Bankole
271ca69d878SAdeleke O. Bankole } else {
272ca69d878SAdeleke O. Bankole PetscCall(PetscViewerASCIIPrintf(viewer, "Wall %" PetscInt_FMT " Forces: Force_x = %12g, Force_y = %12g, Force_z = %12g\n", wall,
273ca69d878SAdeleke O. Bankole reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2]));
274ca69d878SAdeleke O. Bankole }
275ca69d878SAdeleke O. Bankole }
276ca69d878SAdeleke O. Bankole }
277ca69d878SAdeleke O. Bankole PetscCall(PetscFree(reaction_force));
278ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
279ca69d878SAdeleke O. Bankole }
280ca69d878SAdeleke O. Bankole
281f14660a4SJames Wright // User provided TS Monitor
TSMonitor_NS(TS ts,PetscInt step_no,PetscReal time,Vec Q,void * ctx)2822b730f8bSJeremy L Thompson PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) {
283f14660a4SJames Wright User user = ctx;
284f14660a4SJames Wright
285f17d818dSJames Wright PetscFunctionBeginUser;
28637cbb16aSJed Brown // Print every 'checkpoint_interval' steps
287894de27cSJames Wright if (user->app_ctx->checkpoint_interval <= 0 || step_no % user->app_ctx->checkpoint_interval != 0 ||
28849aac155SJeremy L Thompson (user->app_ctx->cont_steps == step_no && step_no != 0)) {
289ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
29049aac155SJeremy L Thompson }
291f14660a4SJames Wright
292f14660a4SJames Wright PetscCall(WriteOutput(user, Q, step_no, time));
293ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
29477841947SLeila Ghaffari }
29577841947SLeila Ghaffari
29677841947SLeila Ghaffari // TS: Create, setup, and solve
TSSolve_NS(DM dm,User user,AppCtx app_ctx,Physics phys,ProblemData problem,Vec * Q,PetscScalar * f_time,TS * ts)29724078fc4SJames Wright PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, ProblemData problem, Vec *Q, PetscScalar *f_time, TS *ts) {
29877841947SLeila Ghaffari MPI_Comm comm = user->comm;
29977841947SLeila Ghaffari TSAdapt adapt;
30077841947SLeila Ghaffari PetscScalar final_time;
30177841947SLeila Ghaffari
302f17d818dSJames Wright PetscFunctionBeginUser;
3032b730f8bSJeremy L Thompson PetscCall(TSCreate(comm, ts));
3042b730f8bSJeremy L Thompson PetscCall(TSSetDM(*ts, dm));
3051a7db67cSJames Wright PetscCall(TSSetApplicationContext(*ts, user));
30677841947SLeila Ghaffari if (phys->implicit) {
3072b730f8bSJeremy L Thompson PetscCall(TSSetType(*ts, TSBDF));
30877841947SLeila Ghaffari if (user->op_ifunction) {
3092b730f8bSJeremy L Thompson PetscCall(TSSetIFunction(*ts, NULL, IFunction_NS, &user));
31077841947SLeila Ghaffari } else { // Implicit integrators can fall back to using an RHSFunction
3112b730f8bSJeremy L Thompson PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
31277841947SLeila Ghaffari }
31391c97f41SJames Wright if (user->mat_ijacobian) {
3142b730f8bSJeremy L Thompson PetscCall(DMTSSetIJacobian(dm, FormIJacobian_NS, &user));
315e334ad8fSJed Brown }
31677841947SLeila Ghaffari } else {
3179ad5e8e4SJames Wright PetscCheck(user->op_rhs_ctx, comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction");
3182b730f8bSJeremy L Thompson PetscCall(TSSetType(*ts, TSRK));
3192b730f8bSJeremy L Thompson PetscCall(TSRKSetType(*ts, TSRK5F));
3202b730f8bSJeremy L Thompson PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user));
32177841947SLeila Ghaffari }
3222b730f8bSJeremy L Thompson PetscCall(TSSetMaxTime(*ts, 500. * user->units->second));
3232b730f8bSJeremy L Thompson PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER));
32451b00d91SJames Wright if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(TSSetErrorIfStepFails(*ts, PETSC_FALSE));
3252b730f8bSJeremy L Thompson PetscCall(TSSetTimeStep(*ts, 1.e-2 * user->units->second));
3262b730f8bSJeremy L Thompson PetscCall(TSGetAdapt(*ts, &adapt));
3272b730f8bSJeremy L Thompson PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * user->units->second, 1.e2 * user->units->second));
3282b730f8bSJeremy L Thompson PetscCall(TSSetFromOptions(*ts));
32991c97f41SJames Wright if (user->mat_ijacobian) {
33091c97f41SJames Wright if (app_ctx->amat_type && !strcmp(app_ctx->amat_type, MATSHELL)) {
33191c97f41SJames Wright SNES snes;
33291c97f41SJames Wright KSP ksp;
3337f2a9303SJames Wright Mat Pmat, Amat;
33491c97f41SJames Wright
33591c97f41SJames Wright PetscCall(TSGetSNES(*ts, &snes));
33691c97f41SJames Wright PetscCall(SNESGetKSP(snes, &ksp));
3377f2a9303SJames Wright PetscCall(CreateSolveOperatorsFromMatCeed(ksp, user->mat_ijacobian, PETSC_FALSE, &Amat, &Pmat));
33891c97f41SJames Wright PetscCall(TSSetIJacobian(*ts, user->mat_ijacobian, Pmat, NULL, NULL));
3397f2a9303SJames Wright PetscCall(MatDestroy(&Amat));
3407f2a9303SJames Wright PetscCall(MatDestroy(&Pmat));
34191c97f41SJames Wright }
34291c97f41SJames Wright }
34315637395SJames Wright user->time_bc_set = -1.0; // require all BCs be updated
344d55646a4SJames Wright if (app_ctx->cont_steps) { // continue from previous timestep data
34577841947SLeila Ghaffari PetscInt count;
34677841947SLeila Ghaffari PetscViewer viewer;
3472b730f8bSJeremy L Thompson
3484de8550aSJed Brown if (app_ctx->cont_time <= 0) { // Legacy files did not include step number and time
3492b730f8bSJeremy L Thompson PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_time_file, FILE_MODE_READ, &viewer));
3504de8550aSJed Brown PetscCall(PetscViewerBinaryRead(viewer, &app_ctx->cont_time, 1, &count, PETSC_REAL));
3512b730f8bSJeremy L Thompson PetscCall(PetscViewerDestroy(&viewer));
3524de8550aSJed Brown PetscCheck(app_ctx->cont_steps != -1, comm, PETSC_ERR_ARG_INCOMP,
3534de8550aSJed Brown "-continue step number not specified, but checkpoint file does not contain a step number (likely written by older code version)");
3544de8550aSJed Brown }
3554de8550aSJed Brown PetscCall(TSSetTime(*ts, app_ctx->cont_time * user->units->second));
356d8e0aecdSJed Brown PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps));
35777841947SLeila Ghaffari }
3588fb33541SJames Wright if (app_ctx->test_type == TESTTYPE_NONE) {
3592b730f8bSJeremy L Thompson PetscCall(TSMonitorSet(*ts, TSMonitor_NS, user, NULL));
3608fb33541SJames Wright }
361ca69d878SAdeleke O. Bankole if (app_ctx->wall_forces.viewer) {
362ca69d878SAdeleke O. Bankole PetscCall(TSMonitorSet(*ts, TSMonitor_WallForce, user, NULL));
363ca69d878SAdeleke O. Bankole }
364b7d66439SJames Wright if (app_ctx->turb_spanstats_enable) {
365f5452247SJames Wright PetscCall(TSMonitorSet(*ts, TSMonitor_TurbulenceStatistics, user, NULL));
366495b9769SJames Wright CeedScalar previous_time = app_ctx->cont_time * user->units->second;
367a424bcd0SJames Wright PetscCallCeed(user->ceed,
368a424bcd0SJames Wright CeedOperatorSetContextDouble(user->spanstats.op_stats_collect_ctx->op, user->spanstats.previous_time_label, &previous_time));
369a175e481SJames Wright }
3704e9802d1SJames Wright if (app_ctx->diff_filter_monitor) PetscCall(TSMonitorSet(*ts, TSMonitor_DifferentialFilter, user, NULL));
37177841947SLeila Ghaffari
3727bc7b61fSJames Wright if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(PrintRunInfo(user, user->phys, problem, *ts));
37377841947SLeila Ghaffari // Solve
374d8e0aecdSJed Brown PetscReal start_time;
375d8e0aecdSJed Brown PetscInt start_step;
3762b730f8bSJeremy L Thompson PetscCall(TSGetTime(*ts, &start_time));
377d8e0aecdSJed Brown PetscCall(TSGetStepNumber(*ts, &start_step));
3787b39487dSJeremy L Thompson
379711a423aSJed Brown PetscCall(PetscLogDefaultBegin()); // So we can use PetscLogStageGetPerfInfo without -log_view
3807b39487dSJeremy L Thompson PetscPreLoadBegin(PETSC_FALSE, "Fluids Solve");
3817b39487dSJeremy L Thompson PetscCall(TSSetTime(*ts, start_time));
382d8e0aecdSJed Brown PetscCall(TSSetStepNumber(*ts, start_step));
3837b39487dSJeremy L Thompson if (PetscPreLoadingOn) {
3847b39487dSJeremy L Thompson // LCOV_EXCL_START
3857b39487dSJeremy L Thompson SNES snes;
3867b39487dSJeremy L Thompson Vec Q_preload;
3877b39487dSJeremy L Thompson PetscReal rtol;
3887b39487dSJeremy L Thompson PetscCall(VecDuplicate(*Q, &Q_preload));
3897b39487dSJeremy L Thompson PetscCall(VecCopy(*Q, Q_preload));
3907b39487dSJeremy L Thompson PetscCall(TSGetSNES(*ts, &snes));
3917b39487dSJeremy L Thompson PetscCall(SNESGetTolerances(snes, NULL, &rtol, NULL, NULL, NULL));
3922b730f8bSJeremy L Thompson PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, .99, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
393fc4c3d69SJames Wright PetscCall(TSSetSolution(*ts, Q_preload));
3947b39487dSJeremy L Thompson PetscCall(TSStep(*ts));
3952b730f8bSJeremy L Thompson PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, rtol, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
3967b39487dSJeremy L Thompson PetscCall(VecDestroy(&Q_preload));
3977b39487dSJeremy L Thompson // LCOV_EXCL_STOP
3987b39487dSJeremy L Thompson } else {
3992b730f8bSJeremy L Thompson PetscCall(PetscBarrier((PetscObject)*ts));
4002b730f8bSJeremy L Thompson PetscCall(TSSolve(*ts, *Q));
4017b39487dSJeremy L Thompson }
4027b39487dSJeremy L Thompson PetscPreLoadEnd();
4037b39487dSJeremy L Thompson
4047b39487dSJeremy L Thompson PetscCall(TSGetSolveTime(*ts, &final_time));
40577841947SLeila Ghaffari *f_time = final_time;
4067b39487dSJeremy L Thompson
4078fb33541SJames Wright if (app_ctx->test_type == TESTTYPE_NONE) {
408f14660a4SJames Wright PetscInt step_no;
409f14660a4SJames Wright PetscCall(TSGetStepNumber(*ts, &step_no));
410a175e481SJames Wright if (user->app_ctx->checkpoint_interval > 0 || user->app_ctx->checkpoint_interval == -1) {
411f14660a4SJames Wright PetscCall(WriteOutput(user, *Q, step_no, final_time));
412f14660a4SJames Wright }
413f14660a4SJames Wright
41475d1798cSJames Wright PetscLogStage stage_id;
415711a423aSJed Brown PetscEventPerfInfo stage_perf;
4167b39487dSJeremy L Thompson
4177b39487dSJeremy L Thompson PetscCall(PetscLogStageGetId("Fluids Solve", &stage_id));
418711a423aSJed Brown PetscCall(PetscLogStageGetPerfInfo(stage_id, &stage_perf));
419711a423aSJed Brown PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time taken for solution (sec): %g\n", stage_perf.time));
42077841947SLeila Ghaffari }
421ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
42277841947SLeila Ghaffari }
423