1*dc936754SJeremy 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 18f8e2d240SJames Wright // @brief Create KSP to solve the inverse mass operator for explicit time stepping schemes 19d6ca8aaaSJames Wright PetscErrorCode CreateKSPMassOperator(User user, CeedData ceed_data) { 20f8e2d240SJames Wright Ceed ceed = user->ceed; 21f8e2d240SJames Wright DM dm = user->dm; 22a515125bSLeila Ghaffari CeedQFunction qf_mass; 23a515125bSLeila Ghaffari CeedOperator op_mass; 24a515125bSLeila Ghaffari CeedInt num_comp_q, q_data_size; 25a515125bSLeila Ghaffari 2606f41313SJames Wright PetscFunctionBeginUser; 27b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q)); 28b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &q_data_size)); 29a515125bSLeila Ghaffari 309f59f36eSJames Wright PetscCall(CreateMassQFunction(ceed, num_comp_q, q_data_size, &qf_mass)); 31b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass)); 32b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 3358e1cbfdSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data)); 34b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 35a515125bSLeila Ghaffari 36f8e2d240SJames Wright { // -- Setup KSP for mass operator 373170c09fSJames Wright Mat mat_mass; 38614823f4SJames Wright Vec Zeros_loc; 39f8e2d240SJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)dm); 40a515125bSLeila Ghaffari 41614823f4SJames Wright PetscCall(DMCreateLocalVector(dm, &Zeros_loc)); 42614823f4SJames Wright PetscCall(VecZeroEntries(Zeros_loc)); 433170c09fSJames Wright PetscCall(MatCeedCreate(dm, dm, op_mass, NULL, &mat_mass)); 443170c09fSJames Wright PetscCall(MatCeedSetLocalVectors(mat_mass, Zeros_loc, NULL)); 45a515125bSLeila Ghaffari 46f8e2d240SJames Wright PetscCall(KSPCreate(comm, &user->mass_ksp)); 47f8e2d240SJames Wright PetscCall(KSPSetOptionsPrefix(user->mass_ksp, "mass_")); 48f8e2d240SJames Wright { // lumped by default 49f8e2d240SJames Wright PC pc; 50f8e2d240SJames Wright PetscCall(KSPGetPC(user->mass_ksp, &pc)); 51f8e2d240SJames Wright PetscCall(PCSetType(pc, PCJACOBI)); 52f8e2d240SJames Wright PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM)); 53f8e2d240SJames Wright PetscCall(KSPSetType(user->mass_ksp, KSPPREONLY)); 54f8e2d240SJames Wright } 553170c09fSJames Wright PetscCall(KSPSetFromOptions_WithMatCeed(user->mass_ksp, mat_mass)); 56614823f4SJames Wright PetscCall(VecDestroy(&Zeros_loc)); 57614823f4SJames Wright PetscCall(MatDestroy(&mat_mass)); 58f8e2d240SJames Wright } 59a515125bSLeila Ghaffari 60b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 61b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass)); 62d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 63a515125bSLeila Ghaffari } 64a515125bSLeila Ghaffari 65c996854bSJames Wright // Insert Boundary values if it's a new time 66c996854bSJames Wright PetscErrorCode UpdateBoundaryValues(User user, Vec Q_loc, PetscReal t) { 67c996854bSJames Wright PetscFunctionBeginUser; 68c996854bSJames Wright if (user->time_bc_set != t) { 69c996854bSJames Wright PetscCall(DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL)); 70c996854bSJames Wright user->time_bc_set = t; 71c996854bSJames Wright } 72d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 73c996854bSJames Wright } 74c996854bSJames Wright 75a515125bSLeila Ghaffari // RHS (Explicit time-stepper) function setup 76a515125bSLeila Ghaffari // This is the RHS of the ODE, given as u_t = G(t,u) 77a515125bSLeila Ghaffari // This function takes in a state vector Q and writes into G 78a515125bSLeila Ghaffari PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) { 79a515125bSLeila Ghaffari User user = *(User *)user_data; 80701e5830SJames Wright Ceed ceed = user->ceed; 81fd969b44SJames Wright PetscScalar dt; 82da5fe0e4SJames Wright Vec Q_loc = user->Q_loc; 83a515125bSLeila Ghaffari 8406f41313SJames Wright PetscFunctionBeginUser; 85e2f84137SJeremy L Thompson // Update time dependent data 86c996854bSJames Wright PetscCall(UpdateBoundaryValues(user, Q_loc, t)); 87701e5830SJames Wright if (user->phys->solution_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_rhs_ctx->op, user->phys->solution_time_label, &t)); 882b916ea7SJeremy L Thompson PetscCall(TSGetTimeStep(ts, &dt)); 89701e5830SJames Wright if (user->phys->timestep_size_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_rhs_ctx->op, user->phys->timestep_size_label, &dt)); 90a515125bSLeila Ghaffari 91da5fe0e4SJames Wright PetscCall(ApplyCeedOperatorGlobalToGlobal(Q, G, user->op_rhs_ctx)); 92a515125bSLeila Ghaffari 93f3fcf8f4SJames Wright // Inverse of the lumped mass matrix 94f8e2d240SJames Wright PetscCall(KSPSolve(user->mass_ksp, G, G)); 95d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 96a515125bSLeila Ghaffari } 97a515125bSLeila Ghaffari 98c5e9980aSAdeleke O. Bankole // Surface forces function setup 99c5e9980aSAdeleke O. Bankole static PetscErrorCode Surface_Forces_NS(DM dm, Vec G_loc, PetscInt num_walls, const PetscInt walls[], PetscScalar *reaction_force) { 100c5e9980aSAdeleke O. Bankole DMLabel face_label; 101c5e9980aSAdeleke O. Bankole const PetscScalar *g; 1022004e3acSAdeleke O. Bankole PetscInt dof, dim = 3; 103c5e9980aSAdeleke O. Bankole MPI_Comm comm; 1042004e3acSAdeleke O. Bankole PetscSection s; 105c5e9980aSAdeleke O. Bankole 106c5e9980aSAdeleke O. Bankole PetscFunctionBeginUser; 107c5e9980aSAdeleke O. Bankole PetscCall(PetscArrayzero(reaction_force, num_walls * dim)); 108c5e9980aSAdeleke O. Bankole PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 109c5e9980aSAdeleke O. Bankole PetscCall(DMGetLabel(dm, "Face Sets", &face_label)); 110c5e9980aSAdeleke O. Bankole PetscCall(VecGetArrayRead(G_loc, &g)); 111c5e9980aSAdeleke O. Bankole for (PetscInt w = 0; w < num_walls; w++) { 112c5e9980aSAdeleke O. Bankole const PetscInt wall = walls[w]; 113c5e9980aSAdeleke O. Bankole IS wall_is; 1142004e3acSAdeleke O. Bankole PetscCall(DMGetLocalSection(dm, &s)); 115c5e9980aSAdeleke O. Bankole PetscCall(DMLabelGetStratumIS(face_label, wall, &wall_is)); 116c5e9980aSAdeleke O. Bankole if (wall_is) { // There exist such points on this process 117c5e9980aSAdeleke O. Bankole PetscInt num_points; 1182004e3acSAdeleke O. Bankole PetscInt num_comp = 0; 119c5e9980aSAdeleke O. Bankole const PetscInt *points; 1202004e3acSAdeleke O. Bankole PetscCall(PetscSectionGetFieldComponents(s, 0, &num_comp)); 121c5e9980aSAdeleke O. Bankole PetscCall(ISGetSize(wall_is, &num_points)); 122c5e9980aSAdeleke O. Bankole PetscCall(ISGetIndices(wall_is, &points)); 123c5e9980aSAdeleke O. Bankole for (PetscInt i = 0; i < num_points; i++) { 124c5e9980aSAdeleke O. Bankole const PetscInt p = points[i]; 125c5e9980aSAdeleke O. Bankole const StateConservative *r; 126c5e9980aSAdeleke O. Bankole PetscCall(DMPlexPointLocalRead(dm, p, g, &r)); 1272004e3acSAdeleke O. Bankole PetscCall(PetscSectionGetDof(s, p, &dof)); 1282004e3acSAdeleke O. Bankole for (PetscInt node = 0; node < dof / num_comp; node++) { 129c5e9980aSAdeleke O. Bankole for (PetscInt j = 0; j < 3; j++) { 1302004e3acSAdeleke O. Bankole reaction_force[w * dim + j] -= r[node].momentum[j]; 1312004e3acSAdeleke O. Bankole } 132c5e9980aSAdeleke O. Bankole } 133c5e9980aSAdeleke O. Bankole } 134c5e9980aSAdeleke O. Bankole PetscCall(ISRestoreIndices(wall_is, &points)); 135c5e9980aSAdeleke O. Bankole } 136c5e9980aSAdeleke O. Bankole PetscCall(ISDestroy(&wall_is)); 137c5e9980aSAdeleke O. Bankole } 138c5e9980aSAdeleke O. Bankole PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, reaction_force, dim * num_walls, MPIU_SCALAR, MPI_SUM, comm)); 139c5e9980aSAdeleke O. Bankole // Restore Vectors 140c5e9980aSAdeleke O. Bankole PetscCall(VecRestoreArrayRead(G_loc, &g)); 141d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 142c5e9980aSAdeleke O. Bankole } 143c5e9980aSAdeleke O. Bankole 144a515125bSLeila Ghaffari // Implicit time-stepper function setup 1452b916ea7SJeremy L Thompson PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) { 146a515125bSLeila Ghaffari User user = *(User *)user_data; 147701e5830SJames Wright Ceed ceed = user->ceed; 148fd969b44SJames Wright PetscScalar dt; 149e2f84137SJeremy L Thompson Vec Q_loc = user->Q_loc, Q_dot_loc = user->Q_dot_loc, G_loc; 150a515125bSLeila Ghaffari PetscMemType q_mem_type, q_dot_mem_type, g_mem_type; 151a515125bSLeila Ghaffari 15206f41313SJames Wright PetscFunctionBeginUser; 153e2f84137SJeremy L Thompson // Get local vectors 154c5e9980aSAdeleke O. Bankole PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); 155e2f84137SJeremy L Thompson 156e2f84137SJeremy L Thompson // Update time dependent data 157c996854bSJames Wright PetscCall(UpdateBoundaryValues(user, Q_loc, t)); 158701e5830SJames Wright if (user->phys->solution_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_ifunction, user->phys->solution_time_label, &t)); 1592b916ea7SJeremy L Thompson PetscCall(TSGetTimeStep(ts, &dt)); 160701e5830SJames Wright if (user->phys->timestep_size_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_ifunction, user->phys->timestep_size_label, &dt)); 161a515125bSLeila Ghaffari 162a515125bSLeila Ghaffari // Global-to-local 16306108310SJames Wright PetscCall(DMGlobalToLocalBegin(user->dm, Q, INSERT_VALUES, Q_loc)); 16406108310SJames Wright PetscCall(DMGlobalToLocalBegin(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc)); 16506108310SJames Wright PetscCall(DMGlobalToLocalEnd(user->dm, Q, INSERT_VALUES, Q_loc)); 16606108310SJames Wright PetscCall(DMGlobalToLocalEnd(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc)); 167a515125bSLeila Ghaffari 168a515125bSLeila Ghaffari // Place PETSc vectors in CEED vectors 169a7dac1d5SJames Wright PetscCall(VecReadPetscToCeed(Q_loc, &q_mem_type, user->q_ceed)); 170a7dac1d5SJames Wright PetscCall(VecReadPetscToCeed(Q_dot_loc, &q_dot_mem_type, user->q_dot_ceed)); 171a7dac1d5SJames Wright PetscCall(VecPetscToCeed(G_loc, &g_mem_type, user->g_ceed)); 172a515125bSLeila Ghaffari 173a515125bSLeila Ghaffari // Apply CEED operator 1747eedc94cSJames Wright PetscCall(PetscLogEventBegin(FLUIDS_CeedOperatorApply, Q, G, 0, 0)); 1757eedc94cSJames Wright PetscCall(PetscLogGpuTimeBegin()); 176b4c37c5cSJames Wright PetscCallCeed(user->ceed, CeedOperatorApply(user->op_ifunction, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE)); 1777eedc94cSJames Wright PetscCall(PetscLogGpuTimeEnd()); 1787eedc94cSJames Wright PetscCall(PetscLogEventEnd(FLUIDS_CeedOperatorApply, Q, G, 0, 0)); 179a515125bSLeila Ghaffari 180a515125bSLeila Ghaffari // Restore vectors 181a7dac1d5SJames Wright PetscCall(VecReadCeedToPetsc(user->q_ceed, q_mem_type, Q_loc)); 182a7dac1d5SJames Wright PetscCall(VecReadCeedToPetsc(user->q_dot_ceed, q_dot_mem_type, Q_dot_loc)); 183a7dac1d5SJames Wright PetscCall(VecCeedToPetsc(user->g_ceed, g_mem_type, G_loc)); 184a515125bSLeila Ghaffari 18501ab89c1SJames Wright if (user->app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) { 186ad494f68SJames Wright PetscCall(SgsDDApplyIFunction(user, Q_loc, G_loc)); 18701ab89c1SJames Wright } 1889c678832SJames Wright 189a515125bSLeila Ghaffari // Local-to-Global 1902b916ea7SJeremy L Thompson PetscCall(VecZeroEntries(G)); 1912b916ea7SJeremy L Thompson PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G)); 192a515125bSLeila Ghaffari 193a515125bSLeila Ghaffari // Restore vectors 194c5e9980aSAdeleke O. Bankole PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); 195d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 196a515125bSLeila Ghaffari } 197a515125bSLeila Ghaffari 1982b916ea7SJeremy L Thompson PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, PetscReal shift, Mat J, Mat J_pre, void *user_data) { 199f0b65372SJed Brown User user = *(User *)user_data; 200b4c37c5cSJames Wright Ceed ceed = user->ceed; 201ebfabadfSJames Wright PetscBool J_is_matceed, J_is_mffd, J_pre_is_matceed, J_pre_is_mffd; 20206f41313SJames Wright 203f0b65372SJed Brown PetscFunctionBeginUser; 20404855949SJed Brown PetscCall(PetscObjectTypeCompare((PetscObject)J, MATMFFD, &J_is_mffd)); 205ebfabadfSJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)J, MATCEED, &J_is_matceed)); 206ebfabadfSJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATMFFD, &J_pre_is_mffd)); 207ebfabadfSJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATCEED, &J_pre_is_matceed)); 208ebfabadfSJames Wright if (user->phys->ijacobian_time_shift_label) { 209ebfabadfSJames Wright CeedOperator op_ijacobian; 210ebfabadfSJames Wright 211ebfabadfSJames Wright PetscCall(MatCeedGetCeedOperators(user->mat_ijacobian, &op_ijacobian, NULL)); 212ebfabadfSJames Wright PetscCallCeed(ceed, CeedOperatorSetContextDouble(op_ijacobian, user->phys->ijacobian_time_shift_label, &shift)); 213f0b65372SJed Brown } 214ebfabadfSJames Wright 215ebfabadfSJames Wright if (J_is_matceed || J_is_mffd) { 21604855949SJed Brown PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 21704855949SJed Brown PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 218ebfabadfSJames Wright } else PetscCall(MatCeedAssembleCOO(user->mat_ijacobian, J)); 219ebfabadfSJames Wright 220ebfabadfSJames Wright if (J_pre_is_matceed && J != J_pre) { 221ebfabadfSJames Wright PetscCall(MatAssemblyBegin(J_pre, MAT_FINAL_ASSEMBLY)); 222ebfabadfSJames Wright PetscCall(MatAssemblyEnd(J_pre, MAT_FINAL_ASSEMBLY)); 223ebfabadfSJames Wright } else if (!J_pre_is_matceed && !J_pre_is_mffd && J != J_pre) { 224ebfabadfSJames Wright PetscCall(MatCeedAssembleCOO(user->mat_ijacobian, J_pre)); 22504855949SJed Brown } 226d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 227f0b65372SJed Brown } 228f0b65372SJed Brown 2292b916ea7SJeremy L Thompson PetscErrorCode WriteOutput(User user, Vec Q, PetscInt step_no, PetscScalar time) { 230a515125bSLeila Ghaffari Vec Q_loc; 231a515125bSLeila Ghaffari char file_path[PETSC_MAX_PATH_LEN]; 232a515125bSLeila Ghaffari PetscViewer viewer; 233a515125bSLeila Ghaffari 23406f41313SJames Wright PetscFunctionBeginUser; 235852e5969SJed Brown if (user->app_ctx->checkpoint_vtk) { 236a515125bSLeila Ghaffari // Set up output 2377538d537SJames Wright PetscCall(DMGetLocalVector(user->dm, &Q_loc)); 2387538d537SJames Wright PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec")); 2397538d537SJames Wright PetscCall(VecZeroEntries(Q_loc)); 2407538d537SJames Wright PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc)); 241a515125bSLeila Ghaffari 242a515125bSLeila Ghaffari // Output 243852e5969SJed Brown PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no)); 2447538d537SJames Wright 2452b916ea7SJeremy L Thompson PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path, FILE_MODE_WRITE, &viewer)); 2467538d537SJames Wright PetscCall(VecView(Q_loc, viewer)); 2477538d537SJames Wright PetscCall(PetscViewerDestroy(&viewer)); 248a515125bSLeila Ghaffari if (user->dm_viz) { 249a515125bSLeila Ghaffari Vec Q_refined, Q_refined_loc; 250a515125bSLeila Ghaffari char file_path_refined[PETSC_MAX_PATH_LEN]; 251a515125bSLeila Ghaffari PetscViewer viewer_refined; 252a515125bSLeila Ghaffari 2537538d537SJames Wright PetscCall(DMGetGlobalVector(user->dm_viz, &Q_refined)); 2547538d537SJames Wright PetscCall(DMGetLocalVector(user->dm_viz, &Q_refined_loc)); 2557538d537SJames Wright PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined")); 2567538d537SJames Wright 2577538d537SJames Wright PetscCall(MatInterpolate(user->interp_viz, Q, Q_refined)); 2587538d537SJames Wright PetscCall(VecZeroEntries(Q_refined_loc)); 2592b916ea7SJeremy L Thompson PetscCall(DMGlobalToLocal(user->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc)); 2607538d537SJames Wright 261852e5969SJed Brown PetscCall( 262852e5969SJed Brown PetscSNPrintf(file_path_refined, sizeof file_path_refined, "%s/nsrefined-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no)); 2637538d537SJames Wright 2642b916ea7SJeremy L Thompson PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined), file_path_refined, FILE_MODE_WRITE, &viewer_refined)); 2657538d537SJames Wright PetscCall(VecView(Q_refined_loc, viewer_refined)); 2667538d537SJames Wright PetscCall(DMRestoreLocalVector(user->dm_viz, &Q_refined_loc)); 2677538d537SJames Wright PetscCall(DMRestoreGlobalVector(user->dm_viz, &Q_refined)); 2687538d537SJames Wright PetscCall(PetscViewerDestroy(&viewer_refined)); 269a515125bSLeila Ghaffari } 2707538d537SJames Wright PetscCall(DMRestoreLocalVector(user->dm, &Q_loc)); 271852e5969SJed Brown } 272a515125bSLeila Ghaffari 273a515125bSLeila Ghaffari // Save data in a binary file for continuation of simulations 27491a36801SJames Wright if (user->app_ctx->add_stepnum2bin) { 275852e5969SJed Brown PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", user->app_ctx->output_dir, step_no)); 27691a36801SJames Wright } else { 2772b916ea7SJeremy L Thompson PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", user->app_ctx->output_dir)); 27891a36801SJames Wright } 2792b916ea7SJeremy L Thompson PetscCall(PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, &viewer)); 2807538d537SJames Wright 281e1233009SJames Wright PetscInt32 token = PetscDefined(USE_64BIT_INDICES) ? FLUIDS_FILE_TOKEN_64 : FLUIDS_FILE_TOKEN_32; 282e1233009SJames Wright PetscCall(PetscViewerBinaryWrite(viewer, &token, 1, PETSC_INT32)); 2839293eaa1SJed Brown PetscCall(PetscViewerBinaryWrite(viewer, &step_no, 1, PETSC_INT)); 2849293eaa1SJed Brown time /= user->units->second; // Dimensionalize time back 2859293eaa1SJed Brown PetscCall(PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL)); 2867538d537SJames Wright PetscCall(VecView(Q, viewer)); 2877538d537SJames Wright PetscCall(PetscViewerDestroy(&viewer)); 288d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 2897538d537SJames Wright } 2907538d537SJames Wright 291c5e9980aSAdeleke O. Bankole // CSV Monitor 292c5e9980aSAdeleke O. Bankole PetscErrorCode TSMonitor_WallForce(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) { 293c5e9980aSAdeleke O. Bankole User user = ctx; 294c5e9980aSAdeleke O. Bankole Vec G_loc; 295c5e9980aSAdeleke O. Bankole PetscInt num_wall = user->app_ctx->wall_forces.num_wall, dim = 3; 296c5e9980aSAdeleke O. Bankole const PetscInt *walls = user->app_ctx->wall_forces.walls; 297c5e9980aSAdeleke O. Bankole PetscViewer viewer = user->app_ctx->wall_forces.viewer; 298c5e9980aSAdeleke O. Bankole PetscViewerFormat format = user->app_ctx->wall_forces.viewer_format; 299c5e9980aSAdeleke O. Bankole PetscScalar *reaction_force; 300c5e9980aSAdeleke O. Bankole PetscBool iascii; 301c5e9980aSAdeleke O. Bankole 302c5e9980aSAdeleke O. Bankole PetscFunctionBeginUser; 303d949ddfcSJames Wright if (!viewer) PetscFunctionReturn(PETSC_SUCCESS); 304c5e9980aSAdeleke O. Bankole PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); 305c5e9980aSAdeleke O. Bankole PetscCall(PetscMalloc1(num_wall * dim, &reaction_force)); 306c5e9980aSAdeleke O. Bankole PetscCall(Surface_Forces_NS(user->dm, G_loc, num_wall, walls, reaction_force)); 307c5e9980aSAdeleke O. Bankole PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); 308c5e9980aSAdeleke O. Bankole 309c5e9980aSAdeleke O. Bankole PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 310c5e9980aSAdeleke O. Bankole 311c5e9980aSAdeleke O. Bankole if (iascii) { 312c5e9980aSAdeleke O. Bankole if (format == PETSC_VIEWER_ASCII_CSV && !user->app_ctx->wall_forces.header_written) { 313c5e9980aSAdeleke O. Bankole PetscCall(PetscViewerASCIIPrintf(viewer, "Step,Time,Wall,ForceX,ForceY,ForceZ\n")); 314c5e9980aSAdeleke O. Bankole user->app_ctx->wall_forces.header_written = PETSC_TRUE; 315c5e9980aSAdeleke O. Bankole } 316c5e9980aSAdeleke O. Bankole for (PetscInt w = 0; w < num_wall; w++) { 317c5e9980aSAdeleke O. Bankole PetscInt wall = walls[w]; 318c5e9980aSAdeleke O. Bankole if (format == PETSC_VIEWER_ASCII_CSV) { 319c5e9980aSAdeleke O. Bankole PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT ",%g,%" PetscInt_FMT ",%g,%g,%g\n", step_no, time, wall, 320c5e9980aSAdeleke O. Bankole reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2])); 321c5e9980aSAdeleke O. Bankole 322c5e9980aSAdeleke O. Bankole } else { 323c5e9980aSAdeleke O. Bankole PetscCall(PetscViewerASCIIPrintf(viewer, "Wall %" PetscInt_FMT " Forces: Force_x = %12g, Force_y = %12g, Force_z = %12g\n", wall, 324c5e9980aSAdeleke O. Bankole reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2])); 325c5e9980aSAdeleke O. Bankole } 326c5e9980aSAdeleke O. Bankole } 327c5e9980aSAdeleke O. Bankole } 328c5e9980aSAdeleke O. Bankole PetscCall(PetscFree(reaction_force)); 329d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 330c5e9980aSAdeleke O. Bankole } 331c5e9980aSAdeleke O. Bankole 3327538d537SJames Wright // User provided TS Monitor 3332b916ea7SJeremy L Thompson PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) { 3347538d537SJames Wright User user = ctx; 3357538d537SJames Wright 33606f41313SJames Wright PetscFunctionBeginUser; 337852e5969SJed Brown // Print every 'checkpoint_interval' steps 338c539088bSJames Wright if (user->app_ctx->checkpoint_interval <= 0 || step_no % user->app_ctx->checkpoint_interval != 0 || 339e419654dSJeremy L Thompson (user->app_ctx->cont_steps == step_no && step_no != 0)) { 340d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 341e419654dSJeremy L Thompson } 3427538d537SJames Wright 3437538d537SJames Wright PetscCall(WriteOutput(user, Q, step_no, time)); 344d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 345a515125bSLeila Ghaffari } 346a515125bSLeila Ghaffari 347a515125bSLeila Ghaffari // TS: Create, setup, and solve 3482b916ea7SJeremy L Thompson PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, Vec *Q, PetscScalar *f_time, TS *ts) { 349a515125bSLeila Ghaffari MPI_Comm comm = user->comm; 350a515125bSLeila Ghaffari TSAdapt adapt; 351a515125bSLeila Ghaffari PetscScalar final_time; 352a515125bSLeila Ghaffari 35306f41313SJames Wright PetscFunctionBeginUser; 3542b916ea7SJeremy L Thompson PetscCall(TSCreate(comm, ts)); 3552b916ea7SJeremy L Thompson PetscCall(TSSetDM(*ts, dm)); 356632a41e1SJames Wright PetscCall(TSSetApplicationContext(*ts, user)); 357a515125bSLeila Ghaffari if (phys->implicit) { 3582b916ea7SJeremy L Thompson PetscCall(TSSetType(*ts, TSBDF)); 359a515125bSLeila Ghaffari if (user->op_ifunction) { 3602b916ea7SJeremy L Thompson PetscCall(TSSetIFunction(*ts, NULL, IFunction_NS, &user)); 361a515125bSLeila Ghaffari } else { // Implicit integrators can fall back to using an RHSFunction 3622b916ea7SJeremy L Thompson PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user)); 363a515125bSLeila Ghaffari } 364ebfabadfSJames Wright if (user->mat_ijacobian) { 3652b916ea7SJeremy L Thompson PetscCall(DMTSSetIJacobian(dm, FormIJacobian_NS, &user)); 366f0b65372SJed Brown } 367a515125bSLeila Ghaffari } else { 368da5fe0e4SJames Wright PetscCheck(user->op_rhs_ctx, comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction"); 3692b916ea7SJeremy L Thompson PetscCall(TSSetType(*ts, TSRK)); 3702b916ea7SJeremy L Thompson PetscCall(TSRKSetType(*ts, TSRK5F)); 3712b916ea7SJeremy L Thompson PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user)); 372a515125bSLeila Ghaffari } 3732b916ea7SJeremy L Thompson PetscCall(TSSetMaxTime(*ts, 500. * user->units->second)); 3742b916ea7SJeremy L Thompson PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER)); 37522387d3aSJames Wright if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(TSSetErrorIfStepFails(*ts, PETSC_FALSE)); 3762b916ea7SJeremy L Thompson PetscCall(TSSetTimeStep(*ts, 1.e-2 * user->units->second)); 3772b916ea7SJeremy L Thompson PetscCall(TSGetAdapt(*ts, &adapt)); 3782b916ea7SJeremy L Thompson PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * user->units->second, 1.e2 * user->units->second)); 3792b916ea7SJeremy L Thompson PetscCall(TSSetFromOptions(*ts)); 380ebfabadfSJames Wright if (user->mat_ijacobian) { 381ebfabadfSJames Wright if (app_ctx->amat_type && !strcmp(app_ctx->amat_type, MATSHELL)) { 382ebfabadfSJames Wright SNES snes; 383ebfabadfSJames Wright KSP ksp; 3843170c09fSJames Wright Mat Pmat, Amat; 385ebfabadfSJames Wright 386ebfabadfSJames Wright PetscCall(TSGetSNES(*ts, &snes)); 387ebfabadfSJames Wright PetscCall(SNESGetKSP(snes, &ksp)); 3883170c09fSJames Wright PetscCall(CreateSolveOperatorsFromMatCeed(ksp, user->mat_ijacobian, PETSC_FALSE, &Amat, &Pmat)); 389ebfabadfSJames Wright PetscCall(TSSetIJacobian(*ts, user->mat_ijacobian, Pmat, NULL, NULL)); 3903170c09fSJames Wright PetscCall(MatDestroy(&Amat)); 3913170c09fSJames Wright PetscCall(MatDestroy(&Pmat)); 392ebfabadfSJames Wright } 393ebfabadfSJames Wright } 39491f639d2SJames Wright user->time_bc_set = -1.0; // require all BCs be updated 395c26b555cSJames Wright if (app_ctx->cont_steps) { // continue from previous timestep data 396a515125bSLeila Ghaffari PetscInt count; 397a515125bSLeila Ghaffari PetscViewer viewer; 3982b916ea7SJeremy L Thompson 3999293eaa1SJed Brown if (app_ctx->cont_time <= 0) { // Legacy files did not include step number and time 4002b916ea7SJeremy L Thompson PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_time_file, FILE_MODE_READ, &viewer)); 4019293eaa1SJed Brown PetscCall(PetscViewerBinaryRead(viewer, &app_ctx->cont_time, 1, &count, PETSC_REAL)); 4022b916ea7SJeremy L Thompson PetscCall(PetscViewerDestroy(&viewer)); 4039293eaa1SJed Brown PetscCheck(app_ctx->cont_steps != -1, comm, PETSC_ERR_ARG_INCOMP, 4049293eaa1SJed Brown "-continue step number not specified, but checkpoint file does not contain a step number (likely written by older code version)"); 4059293eaa1SJed Brown } 4069293eaa1SJed Brown PetscCall(TSSetTime(*ts, app_ctx->cont_time * user->units->second)); 40774a6f4ddSJed Brown PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps)); 408a515125bSLeila Ghaffari } 4090e1e9333SJames Wright if (app_ctx->test_type == TESTTYPE_NONE) { 4102b916ea7SJeremy L Thompson PetscCall(TSMonitorSet(*ts, TSMonitor_NS, user, NULL)); 4110e1e9333SJames Wright } 412c5e9980aSAdeleke O. Bankole if (app_ctx->wall_forces.viewer) { 413c5e9980aSAdeleke O. Bankole PetscCall(TSMonitorSet(*ts, TSMonitor_WallForce, user, NULL)); 414c5e9980aSAdeleke O. Bankole } 415c931fa59SJames Wright if (app_ctx->turb_spanstats_enable) { 41691933550SJames Wright PetscCall(TSMonitorSet(*ts, TSMonitor_TurbulenceStatistics, user, NULL)); 417b8daee98SJames Wright CeedScalar previous_time = app_ctx->cont_time * user->units->second; 418b4c37c5cSJames Wright PetscCallCeed(user->ceed, 419b4c37c5cSJames Wright CeedOperatorSetContextDouble(user->spanstats.op_stats_collect_ctx->op, user->spanstats.previous_time_label, &previous_time)); 420b0488d1fSJames Wright } 42188b07121SJames Wright if (app_ctx->diff_filter_monitor) PetscCall(TSMonitorSet(*ts, TSMonitor_DifferentialFilter, user, NULL)); 422a515125bSLeila Ghaffari 4231c17f66aSJames Wright if (app_ctx->sgs_train_enable) { 4241c17f66aSJames Wright PetscCall(TSMonitorSet(*ts, TSMonitor_SGS_DD_Training, user, NULL)); 4251c17f66aSJames Wright PetscCall(TSSetPostStep(*ts, TSPostStep_SGS_DD_Training)); 4261c17f66aSJames Wright } 427a515125bSLeila Ghaffari // Solve 42874a6f4ddSJed Brown PetscReal start_time; 42974a6f4ddSJed Brown PetscInt start_step; 4302b916ea7SJeremy L Thompson PetscCall(TSGetTime(*ts, &start_time)); 43174a6f4ddSJed Brown PetscCall(TSGetStepNumber(*ts, &start_step)); 43291982731SJeremy L Thompson 433df4304b5SJed Brown PetscCall(PetscLogDefaultBegin()); // So we can use PetscLogStageGetPerfInfo without -log_view 43491982731SJeremy L Thompson PetscPreLoadBegin(PETSC_FALSE, "Fluids Solve"); 43591982731SJeremy L Thompson PetscCall(TSSetTime(*ts, start_time)); 43674a6f4ddSJed Brown PetscCall(TSSetStepNumber(*ts, start_step)); 43791982731SJeremy L Thompson if (PetscPreLoadingOn) { 43891982731SJeremy L Thompson // LCOV_EXCL_START 43991982731SJeremy L Thompson SNES snes; 44091982731SJeremy L Thompson Vec Q_preload; 44191982731SJeremy L Thompson PetscReal rtol; 44291982731SJeremy L Thompson PetscCall(VecDuplicate(*Q, &Q_preload)); 44391982731SJeremy L Thompson PetscCall(VecCopy(*Q, Q_preload)); 44491982731SJeremy L Thompson PetscCall(TSGetSNES(*ts, &snes)); 44591982731SJeremy L Thompson PetscCall(SNESGetTolerances(snes, NULL, &rtol, NULL, NULL, NULL)); 4462b916ea7SJeremy L Thompson PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, .99, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 44722c1b34eSJames Wright PetscCall(TSSetSolution(*ts, Q_preload)); 44891982731SJeremy L Thompson PetscCall(TSStep(*ts)); 4492b916ea7SJeremy L Thompson PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, rtol, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 45091982731SJeremy L Thompson PetscCall(VecDestroy(&Q_preload)); 45191982731SJeremy L Thompson // LCOV_EXCL_STOP 45291982731SJeremy L Thompson } else { 4532b916ea7SJeremy L Thompson PetscCall(PetscBarrier((PetscObject)*ts)); 4542b916ea7SJeremy L Thompson PetscCall(TSSolve(*ts, *Q)); 45591982731SJeremy L Thompson } 45691982731SJeremy L Thompson PetscPreLoadEnd(); 45791982731SJeremy L Thompson 45891982731SJeremy L Thompson PetscCall(TSGetSolveTime(*ts, &final_time)); 459a515125bSLeila Ghaffari *f_time = final_time; 46091982731SJeremy L Thompson 4610e1e9333SJames Wright if (app_ctx->test_type == TESTTYPE_NONE) { 4627538d537SJames Wright PetscInt step_no; 4637538d537SJames Wright PetscCall(TSGetStepNumber(*ts, &step_no)); 464b0488d1fSJames Wright if (user->app_ctx->checkpoint_interval > 0 || user->app_ctx->checkpoint_interval == -1) { 4657538d537SJames Wright PetscCall(WriteOutput(user, *Q, step_no, final_time)); 4667538d537SJames Wright } 4677538d537SJames Wright 4687eedc94cSJames Wright PetscLogStage stage_id; 469df4304b5SJed Brown PetscEventPerfInfo stage_perf; 47091982731SJeremy L Thompson 47191982731SJeremy L Thompson PetscCall(PetscLogStageGetId("Fluids Solve", &stage_id)); 472df4304b5SJed Brown PetscCall(PetscLogStageGetPerfInfo(stage_id, &stage_perf)); 473df4304b5SJed Brown PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time taken for solution (sec): %g\n", stage_perf.time)); 474a515125bSLeila Ghaffari } 475d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 476a515125bSLeila Ghaffari } 477