1dc936754SJeremy L Thompson // Copyright (c) 2017-2024, 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 18c996854bSJames Wright // Insert Boundary values if it's a new time 19c996854bSJames Wright PetscErrorCode UpdateBoundaryValues(User user, Vec Q_loc, PetscReal t) { 20c996854bSJames Wright PetscFunctionBeginUser; 21c996854bSJames Wright if (user->time_bc_set != t) { 22c996854bSJames Wright PetscCall(DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL)); 23c996854bSJames Wright user->time_bc_set = t; 24c996854bSJames Wright } 25d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 26c996854bSJames Wright } 27c996854bSJames Wright 28a515125bSLeila Ghaffari // RHS (Explicit time-stepper) function setup 29a515125bSLeila Ghaffari // This is the RHS of the ODE, given as u_t = G(t,u) 30a515125bSLeila Ghaffari // This function takes in a state vector Q and writes into G 31a515125bSLeila Ghaffari PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) { 32a515125bSLeila Ghaffari User user = *(User *)user_data; 33701e5830SJames Wright Ceed ceed = user->ceed; 34fd969b44SJames Wright PetscScalar dt; 35da5fe0e4SJames Wright Vec Q_loc = user->Q_loc; 36*a78efa86SJames Wright PetscMemType q_mem_type; 37a515125bSLeila Ghaffari 3806f41313SJames Wright PetscFunctionBeginUser; 39e2f84137SJeremy L Thompson // Update time dependent data 40c996854bSJames Wright PetscCall(UpdateBoundaryValues(user, Q_loc, t)); 41701e5830SJames Wright if (user->phys->solution_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_rhs_ctx->op, user->phys->solution_time_label, &t)); 422b916ea7SJeremy L Thompson PetscCall(TSGetTimeStep(ts, &dt)); 43701e5830SJames Wright if (user->phys->timestep_size_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_rhs_ctx->op, user->phys->timestep_size_label, &dt)); 44a515125bSLeila Ghaffari 45da5fe0e4SJames Wright PetscCall(ApplyCeedOperatorGlobalToGlobal(Q, G, user->op_rhs_ctx)); 46a515125bSLeila Ghaffari 47*a78efa86SJames Wright PetscCall(VecReadPetscToCeed(Q_loc, &q_mem_type, user->q_ceed)); 48*a78efa86SJames Wright 49*a78efa86SJames Wright // Inverse of the mass matrix 50f8e2d240SJames Wright PetscCall(KSPSolve(user->mass_ksp, G, G)); 51*a78efa86SJames Wright 52*a78efa86SJames Wright PetscCall(VecReadCeedToPetsc(user->q_ceed, q_mem_type, Q_loc)); 53d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 54a515125bSLeila Ghaffari } 55a515125bSLeila Ghaffari 56c5e9980aSAdeleke O. Bankole // Surface forces function setup 57c5e9980aSAdeleke O. Bankole static PetscErrorCode Surface_Forces_NS(DM dm, Vec G_loc, PetscInt num_walls, const PetscInt walls[], PetscScalar *reaction_force) { 58c5e9980aSAdeleke O. Bankole DMLabel face_label; 59c5e9980aSAdeleke O. Bankole const PetscScalar *g; 602004e3acSAdeleke O. Bankole PetscInt dof, dim = 3; 61c5e9980aSAdeleke O. Bankole MPI_Comm comm; 622004e3acSAdeleke O. Bankole PetscSection s; 63c5e9980aSAdeleke O. Bankole 64c5e9980aSAdeleke O. Bankole PetscFunctionBeginUser; 65c5e9980aSAdeleke O. Bankole PetscCall(PetscArrayzero(reaction_force, num_walls * dim)); 66c5e9980aSAdeleke O. Bankole PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 67c5e9980aSAdeleke O. Bankole PetscCall(DMGetLabel(dm, "Face Sets", &face_label)); 68c5e9980aSAdeleke O. Bankole PetscCall(VecGetArrayRead(G_loc, &g)); 69c5e9980aSAdeleke O. Bankole for (PetscInt w = 0; w < num_walls; w++) { 70c5e9980aSAdeleke O. Bankole const PetscInt wall = walls[w]; 71c5e9980aSAdeleke O. Bankole IS wall_is; 722004e3acSAdeleke O. Bankole PetscCall(DMGetLocalSection(dm, &s)); 73c5e9980aSAdeleke O. Bankole PetscCall(DMLabelGetStratumIS(face_label, wall, &wall_is)); 74c5e9980aSAdeleke O. Bankole if (wall_is) { // There exist such points on this process 75c5e9980aSAdeleke O. Bankole PetscInt num_points; 762004e3acSAdeleke O. Bankole PetscInt num_comp = 0; 77c5e9980aSAdeleke O. Bankole const PetscInt *points; 782004e3acSAdeleke O. Bankole PetscCall(PetscSectionGetFieldComponents(s, 0, &num_comp)); 79c5e9980aSAdeleke O. Bankole PetscCall(ISGetSize(wall_is, &num_points)); 80c5e9980aSAdeleke O. Bankole PetscCall(ISGetIndices(wall_is, &points)); 81c5e9980aSAdeleke O. Bankole for (PetscInt i = 0; i < num_points; i++) { 82c5e9980aSAdeleke O. Bankole const PetscInt p = points[i]; 83c5e9980aSAdeleke O. Bankole const StateConservative *r; 84c5e9980aSAdeleke O. Bankole PetscCall(DMPlexPointLocalRead(dm, p, g, &r)); 852004e3acSAdeleke O. Bankole PetscCall(PetscSectionGetDof(s, p, &dof)); 862004e3acSAdeleke O. Bankole for (PetscInt node = 0; node < dof / num_comp; node++) { 87c5e9980aSAdeleke O. Bankole for (PetscInt j = 0; j < 3; j++) { 882004e3acSAdeleke O. Bankole reaction_force[w * dim + j] -= r[node].momentum[j]; 892004e3acSAdeleke O. Bankole } 90c5e9980aSAdeleke O. Bankole } 91c5e9980aSAdeleke O. Bankole } 92c5e9980aSAdeleke O. Bankole PetscCall(ISRestoreIndices(wall_is, &points)); 93c5e9980aSAdeleke O. Bankole } 94c5e9980aSAdeleke O. Bankole PetscCall(ISDestroy(&wall_is)); 95c5e9980aSAdeleke O. Bankole } 96c5e9980aSAdeleke O. Bankole PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, reaction_force, dim * num_walls, MPIU_SCALAR, MPI_SUM, comm)); 97c5e9980aSAdeleke O. Bankole // Restore Vectors 98c5e9980aSAdeleke O. Bankole PetscCall(VecRestoreArrayRead(G_loc, &g)); 99d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 100c5e9980aSAdeleke O. Bankole } 101c5e9980aSAdeleke O. Bankole 102a515125bSLeila Ghaffari // Implicit time-stepper function setup 1032b916ea7SJeremy L Thompson PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) { 104a515125bSLeila Ghaffari User user = *(User *)user_data; 105701e5830SJames Wright Ceed ceed = user->ceed; 106fd969b44SJames Wright PetscScalar dt; 107e2f84137SJeremy L Thompson Vec Q_loc = user->Q_loc, Q_dot_loc = user->Q_dot_loc, G_loc; 108a515125bSLeila Ghaffari PetscMemType q_mem_type, q_dot_mem_type, g_mem_type; 109a515125bSLeila Ghaffari 11006f41313SJames Wright PetscFunctionBeginUser; 111e2f84137SJeremy L Thompson // Get local vectors 112c5e9980aSAdeleke O. Bankole PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); 113e2f84137SJeremy L Thompson 114e2f84137SJeremy L Thompson // Update time dependent data 115c996854bSJames Wright PetscCall(UpdateBoundaryValues(user, Q_loc, t)); 116701e5830SJames Wright if (user->phys->solution_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_ifunction, user->phys->solution_time_label, &t)); 1172b916ea7SJeremy L Thompson PetscCall(TSGetTimeStep(ts, &dt)); 118701e5830SJames Wright if (user->phys->timestep_size_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_ifunction, user->phys->timestep_size_label, &dt)); 119a515125bSLeila Ghaffari 120a515125bSLeila Ghaffari // Global-to-local 12106108310SJames Wright PetscCall(DMGlobalToLocalBegin(user->dm, Q, INSERT_VALUES, Q_loc)); 12206108310SJames Wright PetscCall(DMGlobalToLocalBegin(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc)); 12306108310SJames Wright PetscCall(DMGlobalToLocalEnd(user->dm, Q, INSERT_VALUES, Q_loc)); 12406108310SJames Wright PetscCall(DMGlobalToLocalEnd(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc)); 125a515125bSLeila Ghaffari 126a515125bSLeila Ghaffari // Place PETSc vectors in CEED vectors 127a7dac1d5SJames Wright PetscCall(VecReadPetscToCeed(Q_loc, &q_mem_type, user->q_ceed)); 128a7dac1d5SJames Wright PetscCall(VecReadPetscToCeed(Q_dot_loc, &q_dot_mem_type, user->q_dot_ceed)); 129a7dac1d5SJames Wright PetscCall(VecPetscToCeed(G_loc, &g_mem_type, user->g_ceed)); 130a515125bSLeila Ghaffari 131a515125bSLeila Ghaffari // Apply CEED operator 1327eedc94cSJames Wright PetscCall(PetscLogEventBegin(FLUIDS_CeedOperatorApply, Q, G, 0, 0)); 1337eedc94cSJames Wright PetscCall(PetscLogGpuTimeBegin()); 134b4c37c5cSJames Wright PetscCallCeed(user->ceed, CeedOperatorApply(user->op_ifunction, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE)); 1357eedc94cSJames Wright PetscCall(PetscLogGpuTimeEnd()); 1367eedc94cSJames Wright PetscCall(PetscLogEventEnd(FLUIDS_CeedOperatorApply, Q, G, 0, 0)); 137a515125bSLeila Ghaffari 138a515125bSLeila Ghaffari // Restore vectors 139a7dac1d5SJames Wright PetscCall(VecReadCeedToPetsc(user->q_ceed, q_mem_type, Q_loc)); 140a7dac1d5SJames Wright PetscCall(VecReadCeedToPetsc(user->q_dot_ceed, q_dot_mem_type, Q_dot_loc)); 141a7dac1d5SJames Wright PetscCall(VecCeedToPetsc(user->g_ceed, g_mem_type, G_loc)); 142a515125bSLeila Ghaffari 14301ab89c1SJames Wright if (user->app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) { 144ad494f68SJames Wright PetscCall(SgsDDApplyIFunction(user, Q_loc, G_loc)); 14501ab89c1SJames Wright } 1469c678832SJames Wright 147a515125bSLeila Ghaffari // Local-to-Global 1482b916ea7SJeremy L Thompson PetscCall(VecZeroEntries(G)); 1492b916ea7SJeremy L Thompson PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G)); 150a515125bSLeila Ghaffari 151a515125bSLeila Ghaffari // Restore vectors 152c5e9980aSAdeleke O. Bankole PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); 153d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 154a515125bSLeila Ghaffari } 155a515125bSLeila Ghaffari 1562b916ea7SJeremy L Thompson PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, PetscReal shift, Mat J, Mat J_pre, void *user_data) { 157f0b65372SJed Brown User user = *(User *)user_data; 158b4c37c5cSJames Wright Ceed ceed = user->ceed; 159ebfabadfSJames Wright PetscBool J_is_matceed, J_is_mffd, J_pre_is_matceed, J_pre_is_mffd; 16006f41313SJames Wright 161f0b65372SJed Brown PetscFunctionBeginUser; 16204855949SJed Brown PetscCall(PetscObjectTypeCompare((PetscObject)J, MATMFFD, &J_is_mffd)); 163ebfabadfSJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)J, MATCEED, &J_is_matceed)); 164ebfabadfSJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATMFFD, &J_pre_is_mffd)); 165ebfabadfSJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATCEED, &J_pre_is_matceed)); 166ebfabadfSJames Wright if (user->phys->ijacobian_time_shift_label) { 167ebfabadfSJames Wright CeedOperator op_ijacobian; 168ebfabadfSJames Wright 169ebfabadfSJames Wright PetscCall(MatCeedGetCeedOperators(user->mat_ijacobian, &op_ijacobian, NULL)); 170ebfabadfSJames Wright PetscCallCeed(ceed, CeedOperatorSetContextDouble(op_ijacobian, user->phys->ijacobian_time_shift_label, &shift)); 171f0b65372SJed Brown } 172ebfabadfSJames Wright 173ebfabadfSJames Wright if (J_is_matceed || J_is_mffd) { 17404855949SJed Brown PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 17504855949SJed Brown PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 176ebfabadfSJames Wright } else PetscCall(MatCeedAssembleCOO(user->mat_ijacobian, J)); 177ebfabadfSJames Wright 178ebfabadfSJames Wright if (J_pre_is_matceed && J != J_pre) { 179ebfabadfSJames Wright PetscCall(MatAssemblyBegin(J_pre, MAT_FINAL_ASSEMBLY)); 180ebfabadfSJames Wright PetscCall(MatAssemblyEnd(J_pre, MAT_FINAL_ASSEMBLY)); 181ebfabadfSJames Wright } else if (!J_pre_is_matceed && !J_pre_is_mffd && J != J_pre) { 182ebfabadfSJames Wright PetscCall(MatCeedAssembleCOO(user->mat_ijacobian, J_pre)); 18304855949SJed Brown } 184d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 185f0b65372SJed Brown } 186f0b65372SJed Brown 1872b916ea7SJeremy L Thompson PetscErrorCode WriteOutput(User user, Vec Q, PetscInt step_no, PetscScalar time) { 188a515125bSLeila Ghaffari Vec Q_loc; 189a515125bSLeila Ghaffari char file_path[PETSC_MAX_PATH_LEN]; 190a515125bSLeila Ghaffari PetscViewer viewer; 191a515125bSLeila Ghaffari 19206f41313SJames Wright PetscFunctionBeginUser; 193852e5969SJed Brown if (user->app_ctx->checkpoint_vtk) { 194a515125bSLeila Ghaffari // Set up output 1957538d537SJames Wright PetscCall(DMGetLocalVector(user->dm, &Q_loc)); 1967538d537SJames Wright PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec")); 1977538d537SJames Wright PetscCall(VecZeroEntries(Q_loc)); 1987538d537SJames Wright PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc)); 199a515125bSLeila Ghaffari 200a515125bSLeila Ghaffari // Output 201852e5969SJed Brown PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no)); 2027538d537SJames Wright 2032b916ea7SJeremy L Thompson PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path, FILE_MODE_WRITE, &viewer)); 2047538d537SJames Wright PetscCall(VecView(Q_loc, viewer)); 2057538d537SJames Wright PetscCall(PetscViewerDestroy(&viewer)); 206a515125bSLeila Ghaffari if (user->dm_viz) { 207a515125bSLeila Ghaffari Vec Q_refined, Q_refined_loc; 208a515125bSLeila Ghaffari char file_path_refined[PETSC_MAX_PATH_LEN]; 209a515125bSLeila Ghaffari PetscViewer viewer_refined; 210a515125bSLeila Ghaffari 2117538d537SJames Wright PetscCall(DMGetGlobalVector(user->dm_viz, &Q_refined)); 2127538d537SJames Wright PetscCall(DMGetLocalVector(user->dm_viz, &Q_refined_loc)); 2137538d537SJames Wright PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined")); 2147538d537SJames Wright 2157538d537SJames Wright PetscCall(MatInterpolate(user->interp_viz, Q, Q_refined)); 2167538d537SJames Wright PetscCall(VecZeroEntries(Q_refined_loc)); 2172b916ea7SJeremy L Thompson PetscCall(DMGlobalToLocal(user->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc)); 2187538d537SJames Wright 219852e5969SJed Brown PetscCall( 220852e5969SJed Brown PetscSNPrintf(file_path_refined, sizeof file_path_refined, "%s/nsrefined-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no)); 2217538d537SJames Wright 2222b916ea7SJeremy L Thompson PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined), file_path_refined, FILE_MODE_WRITE, &viewer_refined)); 2237538d537SJames Wright PetscCall(VecView(Q_refined_loc, viewer_refined)); 2247538d537SJames Wright PetscCall(DMRestoreLocalVector(user->dm_viz, &Q_refined_loc)); 2257538d537SJames Wright PetscCall(DMRestoreGlobalVector(user->dm_viz, &Q_refined)); 2267538d537SJames Wright PetscCall(PetscViewerDestroy(&viewer_refined)); 227a515125bSLeila Ghaffari } 2287538d537SJames Wright PetscCall(DMRestoreLocalVector(user->dm, &Q_loc)); 229852e5969SJed Brown } 230a515125bSLeila Ghaffari 231a515125bSLeila Ghaffari // Save data in a binary file for continuation of simulations 23291a36801SJames Wright if (user->app_ctx->add_stepnum2bin) { 233852e5969SJed Brown PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", user->app_ctx->output_dir, step_no)); 23491a36801SJames Wright } else { 2352b916ea7SJeremy L Thompson PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", user->app_ctx->output_dir)); 23691a36801SJames Wright } 2372b916ea7SJeremy L Thompson PetscCall(PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, &viewer)); 2387538d537SJames Wright 239e1233009SJames Wright PetscInt32 token = PetscDefined(USE_64BIT_INDICES) ? FLUIDS_FILE_TOKEN_64 : FLUIDS_FILE_TOKEN_32; 240e1233009SJames Wright PetscCall(PetscViewerBinaryWrite(viewer, &token, 1, PETSC_INT32)); 2419293eaa1SJed Brown PetscCall(PetscViewerBinaryWrite(viewer, &step_no, 1, PETSC_INT)); 2429293eaa1SJed Brown time /= user->units->second; // Dimensionalize time back 2439293eaa1SJed Brown PetscCall(PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL)); 2447538d537SJames Wright PetscCall(VecView(Q, viewer)); 2457538d537SJames Wright PetscCall(PetscViewerDestroy(&viewer)); 246d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 2477538d537SJames Wright } 2487538d537SJames Wright 249c5e9980aSAdeleke O. Bankole // CSV Monitor 250c5e9980aSAdeleke O. Bankole PetscErrorCode TSMonitor_WallForce(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) { 251c5e9980aSAdeleke O. Bankole User user = ctx; 252c5e9980aSAdeleke O. Bankole Vec G_loc; 253c5e9980aSAdeleke O. Bankole PetscInt num_wall = user->app_ctx->wall_forces.num_wall, dim = 3; 254c5e9980aSAdeleke O. Bankole const PetscInt *walls = user->app_ctx->wall_forces.walls; 255c5e9980aSAdeleke O. Bankole PetscViewer viewer = user->app_ctx->wall_forces.viewer; 256c5e9980aSAdeleke O. Bankole PetscViewerFormat format = user->app_ctx->wall_forces.viewer_format; 257c5e9980aSAdeleke O. Bankole PetscScalar *reaction_force; 258c5e9980aSAdeleke O. Bankole PetscBool iascii; 259c5e9980aSAdeleke O. Bankole 260c5e9980aSAdeleke O. Bankole PetscFunctionBeginUser; 261d949ddfcSJames Wright if (!viewer) PetscFunctionReturn(PETSC_SUCCESS); 262c5e9980aSAdeleke O. Bankole PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); 263c5e9980aSAdeleke O. Bankole PetscCall(PetscMalloc1(num_wall * dim, &reaction_force)); 264c5e9980aSAdeleke O. Bankole PetscCall(Surface_Forces_NS(user->dm, G_loc, num_wall, walls, reaction_force)); 265c5e9980aSAdeleke O. Bankole PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); 266c5e9980aSAdeleke O. Bankole 267c5e9980aSAdeleke O. Bankole PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 268c5e9980aSAdeleke O. Bankole 269c5e9980aSAdeleke O. Bankole if (iascii) { 270c5e9980aSAdeleke O. Bankole if (format == PETSC_VIEWER_ASCII_CSV && !user->app_ctx->wall_forces.header_written) { 271c5e9980aSAdeleke O. Bankole PetscCall(PetscViewerASCIIPrintf(viewer, "Step,Time,Wall,ForceX,ForceY,ForceZ\n")); 272c5e9980aSAdeleke O. Bankole user->app_ctx->wall_forces.header_written = PETSC_TRUE; 273c5e9980aSAdeleke O. Bankole } 274c5e9980aSAdeleke O. Bankole for (PetscInt w = 0; w < num_wall; w++) { 275c5e9980aSAdeleke O. Bankole PetscInt wall = walls[w]; 276c5e9980aSAdeleke O. Bankole if (format == PETSC_VIEWER_ASCII_CSV) { 277c5e9980aSAdeleke O. Bankole PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT ",%g,%" PetscInt_FMT ",%g,%g,%g\n", step_no, time, wall, 278c5e9980aSAdeleke O. Bankole reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2])); 279c5e9980aSAdeleke O. Bankole 280c5e9980aSAdeleke O. Bankole } else { 281c5e9980aSAdeleke O. Bankole PetscCall(PetscViewerASCIIPrintf(viewer, "Wall %" PetscInt_FMT " Forces: Force_x = %12g, Force_y = %12g, Force_z = %12g\n", wall, 282c5e9980aSAdeleke O. Bankole reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2])); 283c5e9980aSAdeleke O. Bankole } 284c5e9980aSAdeleke O. Bankole } 285c5e9980aSAdeleke O. Bankole } 286c5e9980aSAdeleke O. Bankole PetscCall(PetscFree(reaction_force)); 287d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 288c5e9980aSAdeleke O. Bankole } 289c5e9980aSAdeleke O. Bankole 2907538d537SJames Wright // User provided TS Monitor 2912b916ea7SJeremy L Thompson PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) { 2927538d537SJames Wright User user = ctx; 2937538d537SJames Wright 29406f41313SJames Wright PetscFunctionBeginUser; 295852e5969SJed Brown // Print every 'checkpoint_interval' steps 296c539088bSJames Wright if (user->app_ctx->checkpoint_interval <= 0 || step_no % user->app_ctx->checkpoint_interval != 0 || 297e419654dSJeremy L Thompson (user->app_ctx->cont_steps == step_no && step_no != 0)) { 298d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 299e419654dSJeremy L Thompson } 3007538d537SJames Wright 3017538d537SJames Wright PetscCall(WriteOutput(user, Q, step_no, time)); 302d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 303a515125bSLeila Ghaffari } 304a515125bSLeila Ghaffari 305a515125bSLeila Ghaffari // TS: Create, setup, and solve 3062b916ea7SJeremy L Thompson PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, Vec *Q, PetscScalar *f_time, TS *ts) { 307a515125bSLeila Ghaffari MPI_Comm comm = user->comm; 308a515125bSLeila Ghaffari TSAdapt adapt; 309a515125bSLeila Ghaffari PetscScalar final_time; 310a515125bSLeila Ghaffari 31106f41313SJames Wright PetscFunctionBeginUser; 3122b916ea7SJeremy L Thompson PetscCall(TSCreate(comm, ts)); 3132b916ea7SJeremy L Thompson PetscCall(TSSetDM(*ts, dm)); 314632a41e1SJames Wright PetscCall(TSSetApplicationContext(*ts, user)); 315a515125bSLeila Ghaffari if (phys->implicit) { 3162b916ea7SJeremy L Thompson PetscCall(TSSetType(*ts, TSBDF)); 317a515125bSLeila Ghaffari if (user->op_ifunction) { 3182b916ea7SJeremy L Thompson PetscCall(TSSetIFunction(*ts, NULL, IFunction_NS, &user)); 319a515125bSLeila Ghaffari } else { // Implicit integrators can fall back to using an RHSFunction 3202b916ea7SJeremy L Thompson PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user)); 321a515125bSLeila Ghaffari } 322ebfabadfSJames Wright if (user->mat_ijacobian) { 3232b916ea7SJeremy L Thompson PetscCall(DMTSSetIJacobian(dm, FormIJacobian_NS, &user)); 324f0b65372SJed Brown } 325a515125bSLeila Ghaffari } else { 326da5fe0e4SJames Wright PetscCheck(user->op_rhs_ctx, comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction"); 3272b916ea7SJeremy L Thompson PetscCall(TSSetType(*ts, TSRK)); 3282b916ea7SJeremy L Thompson PetscCall(TSRKSetType(*ts, TSRK5F)); 3292b916ea7SJeremy L Thompson PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user)); 330a515125bSLeila Ghaffari } 3312b916ea7SJeremy L Thompson PetscCall(TSSetMaxTime(*ts, 500. * user->units->second)); 3322b916ea7SJeremy L Thompson PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER)); 33322387d3aSJames Wright if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(TSSetErrorIfStepFails(*ts, PETSC_FALSE)); 3342b916ea7SJeremy L Thompson PetscCall(TSSetTimeStep(*ts, 1.e-2 * user->units->second)); 3352b916ea7SJeremy L Thompson PetscCall(TSGetAdapt(*ts, &adapt)); 3362b916ea7SJeremy L Thompson PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * user->units->second, 1.e2 * user->units->second)); 3372b916ea7SJeremy L Thompson PetscCall(TSSetFromOptions(*ts)); 338ebfabadfSJames Wright if (user->mat_ijacobian) { 339ebfabadfSJames Wright if (app_ctx->amat_type && !strcmp(app_ctx->amat_type, MATSHELL)) { 340ebfabadfSJames Wright SNES snes; 341ebfabadfSJames Wright KSP ksp; 3423170c09fSJames Wright Mat Pmat, Amat; 343ebfabadfSJames Wright 344ebfabadfSJames Wright PetscCall(TSGetSNES(*ts, &snes)); 345ebfabadfSJames Wright PetscCall(SNESGetKSP(snes, &ksp)); 3463170c09fSJames Wright PetscCall(CreateSolveOperatorsFromMatCeed(ksp, user->mat_ijacobian, PETSC_FALSE, &Amat, &Pmat)); 347ebfabadfSJames Wright PetscCall(TSSetIJacobian(*ts, user->mat_ijacobian, Pmat, NULL, NULL)); 3483170c09fSJames Wright PetscCall(MatDestroy(&Amat)); 3493170c09fSJames Wright PetscCall(MatDestroy(&Pmat)); 350ebfabadfSJames Wright } 351ebfabadfSJames Wright } 35291f639d2SJames Wright user->time_bc_set = -1.0; // require all BCs be updated 353c26b555cSJames Wright if (app_ctx->cont_steps) { // continue from previous timestep data 354a515125bSLeila Ghaffari PetscInt count; 355a515125bSLeila Ghaffari PetscViewer viewer; 3562b916ea7SJeremy L Thompson 3579293eaa1SJed Brown if (app_ctx->cont_time <= 0) { // Legacy files did not include step number and time 3582b916ea7SJeremy L Thompson PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_time_file, FILE_MODE_READ, &viewer)); 3599293eaa1SJed Brown PetscCall(PetscViewerBinaryRead(viewer, &app_ctx->cont_time, 1, &count, PETSC_REAL)); 3602b916ea7SJeremy L Thompson PetscCall(PetscViewerDestroy(&viewer)); 3619293eaa1SJed Brown PetscCheck(app_ctx->cont_steps != -1, comm, PETSC_ERR_ARG_INCOMP, 3629293eaa1SJed Brown "-continue step number not specified, but checkpoint file does not contain a step number (likely written by older code version)"); 3639293eaa1SJed Brown } 3649293eaa1SJed Brown PetscCall(TSSetTime(*ts, app_ctx->cont_time * user->units->second)); 36574a6f4ddSJed Brown PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps)); 366a515125bSLeila Ghaffari } 3670e1e9333SJames Wright if (app_ctx->test_type == TESTTYPE_NONE) { 3682b916ea7SJeremy L Thompson PetscCall(TSMonitorSet(*ts, TSMonitor_NS, user, NULL)); 3690e1e9333SJames Wright } 370c5e9980aSAdeleke O. Bankole if (app_ctx->wall_forces.viewer) { 371c5e9980aSAdeleke O. Bankole PetscCall(TSMonitorSet(*ts, TSMonitor_WallForce, user, NULL)); 372c5e9980aSAdeleke O. Bankole } 373c931fa59SJames Wright if (app_ctx->turb_spanstats_enable) { 37491933550SJames Wright PetscCall(TSMonitorSet(*ts, TSMonitor_TurbulenceStatistics, user, NULL)); 375b8daee98SJames Wright CeedScalar previous_time = app_ctx->cont_time * user->units->second; 376b4c37c5cSJames Wright PetscCallCeed(user->ceed, 377b4c37c5cSJames Wright CeedOperatorSetContextDouble(user->spanstats.op_stats_collect_ctx->op, user->spanstats.previous_time_label, &previous_time)); 378b0488d1fSJames Wright } 37988b07121SJames Wright if (app_ctx->diff_filter_monitor) PetscCall(TSMonitorSet(*ts, TSMonitor_DifferentialFilter, user, NULL)); 380a515125bSLeila Ghaffari 3811c17f66aSJames Wright if (app_ctx->sgs_train_enable) { 3821c17f66aSJames Wright PetscCall(TSMonitorSet(*ts, TSMonitor_SGS_DD_Training, user, NULL)); 3831c17f66aSJames Wright PetscCall(TSSetPostStep(*ts, TSPostStep_SGS_DD_Training)); 3841c17f66aSJames Wright } 385a515125bSLeila Ghaffari // Solve 38674a6f4ddSJed Brown PetscReal start_time; 38774a6f4ddSJed Brown PetscInt start_step; 3882b916ea7SJeremy L Thompson PetscCall(TSGetTime(*ts, &start_time)); 38974a6f4ddSJed Brown PetscCall(TSGetStepNumber(*ts, &start_step)); 39091982731SJeremy L Thompson 391df4304b5SJed Brown PetscCall(PetscLogDefaultBegin()); // So we can use PetscLogStageGetPerfInfo without -log_view 39291982731SJeremy L Thompson PetscPreLoadBegin(PETSC_FALSE, "Fluids Solve"); 39391982731SJeremy L Thompson PetscCall(TSSetTime(*ts, start_time)); 39474a6f4ddSJed Brown PetscCall(TSSetStepNumber(*ts, start_step)); 39591982731SJeremy L Thompson if (PetscPreLoadingOn) { 39691982731SJeremy L Thompson // LCOV_EXCL_START 39791982731SJeremy L Thompson SNES snes; 39891982731SJeremy L Thompson Vec Q_preload; 39991982731SJeremy L Thompson PetscReal rtol; 40091982731SJeremy L Thompson PetscCall(VecDuplicate(*Q, &Q_preload)); 40191982731SJeremy L Thompson PetscCall(VecCopy(*Q, Q_preload)); 40291982731SJeremy L Thompson PetscCall(TSGetSNES(*ts, &snes)); 40391982731SJeremy L Thompson PetscCall(SNESGetTolerances(snes, NULL, &rtol, NULL, NULL, NULL)); 4042b916ea7SJeremy L Thompson PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, .99, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 40522c1b34eSJames Wright PetscCall(TSSetSolution(*ts, Q_preload)); 40691982731SJeremy L Thompson PetscCall(TSStep(*ts)); 4072b916ea7SJeremy L Thompson PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, rtol, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 40891982731SJeremy L Thompson PetscCall(VecDestroy(&Q_preload)); 40991982731SJeremy L Thompson // LCOV_EXCL_STOP 41091982731SJeremy L Thompson } else { 4112b916ea7SJeremy L Thompson PetscCall(PetscBarrier((PetscObject)*ts)); 4122b916ea7SJeremy L Thompson PetscCall(TSSolve(*ts, *Q)); 41391982731SJeremy L Thompson } 41491982731SJeremy L Thompson PetscPreLoadEnd(); 41591982731SJeremy L Thompson 41691982731SJeremy L Thompson PetscCall(TSGetSolveTime(*ts, &final_time)); 417a515125bSLeila Ghaffari *f_time = final_time; 41891982731SJeremy L Thompson 4190e1e9333SJames Wright if (app_ctx->test_type == TESTTYPE_NONE) { 4207538d537SJames Wright PetscInt step_no; 4217538d537SJames Wright PetscCall(TSGetStepNumber(*ts, &step_no)); 422b0488d1fSJames Wright if (user->app_ctx->checkpoint_interval > 0 || user->app_ctx->checkpoint_interval == -1) { 4237538d537SJames Wright PetscCall(WriteOutput(user, *Q, step_no, final_time)); 4247538d537SJames Wright } 4257538d537SJames Wright 4267eedc94cSJames Wright PetscLogStage stage_id; 427df4304b5SJed Brown PetscEventPerfInfo stage_perf; 42891982731SJeremy L Thompson 42991982731SJeremy L Thompson PetscCall(PetscLogStageGetId("Fluids Solve", &stage_id)); 430df4304b5SJed Brown PetscCall(PetscLogStageGetPerfInfo(stage_id, &stage_perf)); 431df4304b5SJed Brown PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time taken for solution (sec): %g\n", stage_perf.time)); 432a515125bSLeila Ghaffari } 433d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 434a515125bSLeila Ghaffari } 435