13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 377841947SLeila Ghaffari // 43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 577841947SLeila Ghaffari // 63d8e8822SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 777841947SLeila Ghaffari 877841947SLeila Ghaffari /// @file 977841947SLeila Ghaffari /// Time-stepping functions for Navier-Stokes example using PETSc 1077841947SLeila Ghaffari 1149aac155SJeremy L Thompson #include <ceed.h> 1249aac155SJeremy L Thompson #include <petscdmplex.h> 1349aac155SJeremy L Thompson #include <petscts.h> 1449aac155SJeremy L Thompson 1577841947SLeila Ghaffari #include "../navierstokes.h" 16ca69d878SAdeleke O. Bankole #include "../qfunctions/newtonian_state.h" 1777841947SLeila Ghaffari 1804292f7dSJames Wright // @brief Create KSP to solve the inverse mass operator for explicit time stepping schemes 19dfcf44b0SJames Wright PetscErrorCode CreateKSPMassOperator(User user, CeedData ceed_data) { 2004292f7dSJames Wright Ceed ceed = user->ceed; 2104292f7dSJames Wright DM dm = user->dm; 2277841947SLeila Ghaffari CeedQFunction qf_mass; 2377841947SLeila Ghaffari CeedOperator op_mass; 2404292f7dSJames Wright OperatorApplyContext mass_matop_ctx; 2577841947SLeila Ghaffari CeedInt num_comp_q, q_data_size; 2677841947SLeila Ghaffari 27f17d818dSJames Wright PetscFunctionBeginUser; 28a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q)); 29a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &q_data_size)); 3077841947SLeila Ghaffari 31ef080ff9SJames Wright PetscCall(CreateMassQFunction(ceed, num_comp_q, q_data_size, &qf_mass)); 32a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass)); 33a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 34356036faSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data)); 35a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 3677841947SLeila Ghaffari 3704292f7dSJames Wright { // -- Setup KSP for mass operator 38ed9ed3deSJames Wright Mat mat_mass = NULL; 390f2fa9b4SJames Wright Vec Zeros_loc; 4004292f7dSJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)dm); 4177841947SLeila Ghaffari 420f2fa9b4SJames Wright PetscCall(DMCreateLocalVector(dm, &Zeros_loc)); 430f2fa9b4SJames Wright PetscCall(VecZeroEntries(Zeros_loc)); 440f2fa9b4SJames Wright PetscCall(OperatorApplyContextCreate(dm, dm, ceed, op_mass, NULL, NULL, Zeros_loc, NULL, &mass_matop_ctx)); 4504292f7dSJames Wright PetscCall(CreateMatShell_Ceed(mass_matop_ctx, &mat_mass)); 4677841947SLeila Ghaffari 4704292f7dSJames Wright PetscCall(KSPCreate(comm, &user->mass_ksp)); 4804292f7dSJames Wright PetscCall(KSPSetOptionsPrefix(user->mass_ksp, "mass_")); 4904292f7dSJames Wright { // lumped by default 5004292f7dSJames Wright PC pc; 5104292f7dSJames Wright PetscCall(KSPGetPC(user->mass_ksp, &pc)); 5204292f7dSJames Wright PetscCall(PCSetType(pc, PCJACOBI)); 5304292f7dSJames Wright PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM)); 5404292f7dSJames Wright PetscCall(KSPSetType(user->mass_ksp, KSPPREONLY)); 5504292f7dSJames Wright } 5604292f7dSJames Wright PetscCall(KSPSetOperators(user->mass_ksp, mat_mass, mat_mass)); 5704292f7dSJames Wright PetscCall(KSPSetFromOptions(user->mass_ksp)); 580f2fa9b4SJames Wright PetscCall(VecDestroy(&Zeros_loc)); 590f2fa9b4SJames Wright PetscCall(MatDestroy(&mat_mass)); 6004292f7dSJames Wright } 6177841947SLeila Ghaffari 6277841947SLeila Ghaffari // Cleanup 63a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 64a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass)); 65ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 6677841947SLeila Ghaffari } 6777841947SLeila Ghaffari 68a0b9a424SJames Wright // Insert Boundary values if it's a new time 69a0b9a424SJames Wright PetscErrorCode UpdateBoundaryValues(User user, Vec Q_loc, PetscReal t) { 70a0b9a424SJames Wright PetscFunctionBeginUser; 71a0b9a424SJames Wright if (user->time_bc_set != t) { 72a0b9a424SJames Wright PetscCall(DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL)); 73a0b9a424SJames Wright user->time_bc_set = t; 74a0b9a424SJames Wright } 75ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 76a0b9a424SJames Wright } 77a0b9a424SJames Wright 7877841947SLeila Ghaffari // RHS (Explicit time-stepper) function setup 7977841947SLeila Ghaffari // This is the RHS of the ODE, given as u_t = G(t,u) 8077841947SLeila Ghaffari // This function takes in a state vector Q and writes into G 8177841947SLeila Ghaffari PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) { 8277841947SLeila Ghaffari User user = *(User *)user_data; 837d4c6defSJames Wright Ceed ceed = user->ceed; 84c798d55aSJames Wright PetscScalar dt; 859ad5e8e4SJames Wright Vec Q_loc = user->Q_loc; 8677841947SLeila Ghaffari 87f17d818dSJames Wright PetscFunctionBeginUser; 885e82a6e1SJeremy L Thompson // Update time dependent data 89a0b9a424SJames Wright PetscCall(UpdateBoundaryValues(user, Q_loc, t)); 907d4c6defSJames Wright if (user->phys->solution_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_rhs_ctx->op, user->phys->solution_time_label, &t)); 912b730f8bSJeremy L Thompson PetscCall(TSGetTimeStep(ts, &dt)); 927d4c6defSJames Wright if (user->phys->timestep_size_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_rhs_ctx->op, user->phys->timestep_size_label, &dt)); 9377841947SLeila Ghaffari 949ad5e8e4SJames Wright PetscCall(ApplyCeedOperatorGlobalToGlobal(Q, G, user->op_rhs_ctx)); 9577841947SLeila Ghaffari 96d6e67e47SJames Wright // Inverse of the lumped mass matrix 9704292f7dSJames Wright PetscCall(KSPSolve(user->mass_ksp, G, G)); 98ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 9977841947SLeila Ghaffari } 10077841947SLeila Ghaffari 101ca69d878SAdeleke O. Bankole // Surface forces function setup 102ca69d878SAdeleke O. Bankole static PetscErrorCode Surface_Forces_NS(DM dm, Vec G_loc, PetscInt num_walls, const PetscInt walls[], PetscScalar *reaction_force) { 103ca69d878SAdeleke O. Bankole DMLabel face_label; 104ca69d878SAdeleke O. Bankole const PetscScalar *g; 105d6734f85SAdeleke O. Bankole PetscInt dof, dim = 3; 106ca69d878SAdeleke O. Bankole MPI_Comm comm; 107d6734f85SAdeleke O. Bankole PetscSection s; 108ca69d878SAdeleke O. Bankole 109ca69d878SAdeleke O. Bankole PetscFunctionBeginUser; 110ca69d878SAdeleke O. Bankole PetscCall(PetscArrayzero(reaction_force, num_walls * dim)); 111ca69d878SAdeleke O. Bankole PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 112ca69d878SAdeleke O. Bankole PetscCall(DMGetLabel(dm, "Face Sets", &face_label)); 113ca69d878SAdeleke O. Bankole PetscCall(VecGetArrayRead(G_loc, &g)); 114ca69d878SAdeleke O. Bankole for (PetscInt w = 0; w < num_walls; w++) { 115ca69d878SAdeleke O. Bankole const PetscInt wall = walls[w]; 116ca69d878SAdeleke O. Bankole IS wall_is; 117d6734f85SAdeleke O. Bankole PetscCall(DMGetLocalSection(dm, &s)); 118ca69d878SAdeleke O. Bankole PetscCall(DMLabelGetStratumIS(face_label, wall, &wall_is)); 119ca69d878SAdeleke O. Bankole if (wall_is) { // There exist such points on this process 120ca69d878SAdeleke O. Bankole PetscInt num_points; 121d6734f85SAdeleke O. Bankole PetscInt num_comp = 0; 122ca69d878SAdeleke O. Bankole const PetscInt *points; 123d6734f85SAdeleke O. Bankole PetscCall(PetscSectionGetFieldComponents(s, 0, &num_comp)); 124ca69d878SAdeleke O. Bankole PetscCall(ISGetSize(wall_is, &num_points)); 125ca69d878SAdeleke O. Bankole PetscCall(ISGetIndices(wall_is, &points)); 126ca69d878SAdeleke O. Bankole for (PetscInt i = 0; i < num_points; i++) { 127ca69d878SAdeleke O. Bankole const PetscInt p = points[i]; 128ca69d878SAdeleke O. Bankole const StateConservative *r; 129ca69d878SAdeleke O. Bankole PetscCall(DMPlexPointLocalRead(dm, p, g, &r)); 130d6734f85SAdeleke O. Bankole PetscCall(PetscSectionGetDof(s, p, &dof)); 131d6734f85SAdeleke O. Bankole for (PetscInt node = 0; node < dof / num_comp; node++) { 132ca69d878SAdeleke O. Bankole for (PetscInt j = 0; j < 3; j++) { 133d6734f85SAdeleke O. Bankole reaction_force[w * dim + j] -= r[node].momentum[j]; 134d6734f85SAdeleke O. Bankole } 135ca69d878SAdeleke O. Bankole } 136ca69d878SAdeleke O. Bankole } 137ca69d878SAdeleke O. Bankole PetscCall(ISRestoreIndices(wall_is, &points)); 138ca69d878SAdeleke O. Bankole } 139ca69d878SAdeleke O. Bankole PetscCall(ISDestroy(&wall_is)); 140ca69d878SAdeleke O. Bankole } 141ca69d878SAdeleke O. Bankole PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, reaction_force, dim * num_walls, MPIU_SCALAR, MPI_SUM, comm)); 142ca69d878SAdeleke O. Bankole // Restore Vectors 143ca69d878SAdeleke O. Bankole PetscCall(VecRestoreArrayRead(G_loc, &g)); 144ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 145ca69d878SAdeleke O. Bankole } 146ca69d878SAdeleke O. Bankole 14777841947SLeila Ghaffari // Implicit time-stepper function setup 1482b730f8bSJeremy L Thompson PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) { 14977841947SLeila Ghaffari User user = *(User *)user_data; 1507d4c6defSJames Wright Ceed ceed = user->ceed; 151c798d55aSJames Wright PetscScalar dt; 1525e82a6e1SJeremy L Thompson Vec Q_loc = user->Q_loc, Q_dot_loc = user->Q_dot_loc, G_loc; 15377841947SLeila Ghaffari PetscMemType q_mem_type, q_dot_mem_type, g_mem_type; 15477841947SLeila Ghaffari 155f17d818dSJames Wright PetscFunctionBeginUser; 1565e82a6e1SJeremy L Thompson // Get local vectors 157ca69d878SAdeleke O. Bankole PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); 1585e82a6e1SJeremy L Thompson 1595e82a6e1SJeremy L Thompson // Update time dependent data 160a0b9a424SJames Wright PetscCall(UpdateBoundaryValues(user, Q_loc, t)); 1617d4c6defSJames Wright if (user->phys->solution_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_ifunction, user->phys->solution_time_label, &t)); 1622b730f8bSJeremy L Thompson PetscCall(TSGetTimeStep(ts, &dt)); 1637d4c6defSJames Wright if (user->phys->timestep_size_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_ifunction, user->phys->timestep_size_label, &dt)); 16477841947SLeila Ghaffari 16577841947SLeila Ghaffari // Global-to-local 166878eb0e7SJames Wright PetscCall(DMGlobalToLocalBegin(user->dm, Q, INSERT_VALUES, Q_loc)); 167878eb0e7SJames Wright PetscCall(DMGlobalToLocalBegin(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc)); 168878eb0e7SJames Wright PetscCall(DMGlobalToLocalEnd(user->dm, Q, INSERT_VALUES, Q_loc)); 169878eb0e7SJames Wright PetscCall(DMGlobalToLocalEnd(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc)); 17077841947SLeila Ghaffari 17177841947SLeila Ghaffari // Place PETSc vectors in CEED vectors 172c798d55aSJames Wright PetscCall(VecReadP2C(Q_loc, &q_mem_type, user->q_ceed)); 173c798d55aSJames Wright PetscCall(VecReadP2C(Q_dot_loc, &q_dot_mem_type, user->q_dot_ceed)); 174c798d55aSJames Wright PetscCall(VecP2C(G_loc, &g_mem_type, user->g_ceed)); 17577841947SLeila Ghaffari 17677841947SLeila Ghaffari // Apply CEED operator 17775d1798cSJames Wright PetscCall(PetscLogEventBegin(FLUIDS_CeedOperatorApply, Q, G, 0, 0)); 17875d1798cSJames Wright PetscCall(PetscLogGpuTimeBegin()); 179a424bcd0SJames Wright PetscCallCeed(user->ceed, CeedOperatorApply(user->op_ifunction, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE)); 18075d1798cSJames Wright PetscCall(PetscLogGpuTimeEnd()); 18175d1798cSJames Wright PetscCall(PetscLogEventEnd(FLUIDS_CeedOperatorApply, Q, G, 0, 0)); 18277841947SLeila Ghaffari 18377841947SLeila Ghaffari // Restore vectors 184c798d55aSJames Wright PetscCall(VecReadC2P(user->q_ceed, q_mem_type, Q_loc)); 185c798d55aSJames Wright PetscCall(VecReadC2P(user->q_dot_ceed, q_dot_mem_type, Q_dot_loc)); 186c798d55aSJames Wright PetscCall(VecC2P(user->g_ceed, g_mem_type, G_loc)); 18777841947SLeila Ghaffari 18819ffbc25SJames Wright if (user->app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) { 1898f5ab23bSJames Wright PetscCall(SgsDDApplyIFunction(user, Q_loc, G_loc)); 19019ffbc25SJames Wright } 191be34b3b7SJames Wright 19277841947SLeila Ghaffari // Local-to-Global 1932b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(G)); 1942b730f8bSJeremy L Thompson PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G)); 19577841947SLeila Ghaffari 19677841947SLeila Ghaffari // Restore vectors 197ca69d878SAdeleke O. Bankole PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); 198ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 19977841947SLeila Ghaffari } 20077841947SLeila Ghaffari 2012b730f8bSJeremy L Thompson PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, PetscReal shift, Mat J, Mat J_pre, void *user_data) { 202e334ad8fSJed Brown User user = *(User *)user_data; 203a424bcd0SJames Wright Ceed ceed = user->ceed; 204*91c97f41SJames Wright PetscBool J_is_matceed, J_is_mffd, J_pre_is_matceed, J_pre_is_mffd; 205f17d818dSJames Wright 206e334ad8fSJed Brown PetscFunctionBeginUser; 207d6bf345cSJed Brown PetscCall(PetscObjectTypeCompare((PetscObject)J, MATMFFD, &J_is_mffd)); 208*91c97f41SJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)J, MATCEED, &J_is_matceed)); 209*91c97f41SJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATMFFD, &J_pre_is_mffd)); 210*91c97f41SJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATCEED, &J_pre_is_matceed)); 211*91c97f41SJames Wright if (user->phys->ijacobian_time_shift_label) { 212*91c97f41SJames Wright CeedOperator op_ijacobian; 213*91c97f41SJames Wright 214*91c97f41SJames Wright PetscCall(MatCeedGetCeedOperators(user->mat_ijacobian, &op_ijacobian, NULL)); 215*91c97f41SJames Wright PetscCallCeed(ceed, CeedOperatorSetContextDouble(op_ijacobian, user->phys->ijacobian_time_shift_label, &shift)); 216e334ad8fSJed Brown } 217*91c97f41SJames Wright 218*91c97f41SJames Wright if (J_is_matceed || J_is_mffd) { 219d6bf345cSJed Brown PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 220d6bf345cSJed Brown PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 221*91c97f41SJames Wright } else PetscCall(MatCeedAssembleCOO(user->mat_ijacobian, J)); 222*91c97f41SJames Wright 223*91c97f41SJames Wright if (J_pre_is_matceed && J != J_pre) { 224*91c97f41SJames Wright PetscCall(MatAssemblyBegin(J_pre, MAT_FINAL_ASSEMBLY)); 225*91c97f41SJames Wright PetscCall(MatAssemblyEnd(J_pre, MAT_FINAL_ASSEMBLY)); 226*91c97f41SJames Wright } else if (!J_pre_is_matceed && !J_pre_is_mffd && J != J_pre) { 227*91c97f41SJames Wright PetscCall(MatCeedAssembleCOO(user->mat_ijacobian, J_pre)); 228d6bf345cSJed Brown } 229ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 230e334ad8fSJed Brown } 231e334ad8fSJed Brown 2322b730f8bSJeremy L Thompson PetscErrorCode WriteOutput(User user, Vec Q, PetscInt step_no, PetscScalar time) { 23377841947SLeila Ghaffari Vec Q_loc; 23477841947SLeila Ghaffari char file_path[PETSC_MAX_PATH_LEN]; 23577841947SLeila Ghaffari PetscViewer viewer; 23677841947SLeila Ghaffari 237f17d818dSJames Wright PetscFunctionBeginUser; 23837cbb16aSJed Brown if (user->app_ctx->checkpoint_vtk) { 23977841947SLeila Ghaffari // Set up output 240f14660a4SJames Wright PetscCall(DMGetLocalVector(user->dm, &Q_loc)); 241f14660a4SJames Wright PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec")); 242f14660a4SJames Wright PetscCall(VecZeroEntries(Q_loc)); 243f14660a4SJames Wright PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc)); 24477841947SLeila Ghaffari 24577841947SLeila Ghaffari // Output 24637cbb16aSJed Brown PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no)); 247f14660a4SJames Wright 2482b730f8bSJeremy L Thompson PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path, FILE_MODE_WRITE, &viewer)); 249f14660a4SJames Wright PetscCall(VecView(Q_loc, viewer)); 250f14660a4SJames Wright PetscCall(PetscViewerDestroy(&viewer)); 25177841947SLeila Ghaffari if (user->dm_viz) { 25277841947SLeila Ghaffari Vec Q_refined, Q_refined_loc; 25377841947SLeila Ghaffari char file_path_refined[PETSC_MAX_PATH_LEN]; 25477841947SLeila Ghaffari PetscViewer viewer_refined; 25577841947SLeila Ghaffari 256f14660a4SJames Wright PetscCall(DMGetGlobalVector(user->dm_viz, &Q_refined)); 257f14660a4SJames Wright PetscCall(DMGetLocalVector(user->dm_viz, &Q_refined_loc)); 258f14660a4SJames Wright PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined")); 259f14660a4SJames Wright 260f14660a4SJames Wright PetscCall(MatInterpolate(user->interp_viz, Q, Q_refined)); 261f14660a4SJames Wright PetscCall(VecZeroEntries(Q_refined_loc)); 2622b730f8bSJeremy L Thompson PetscCall(DMGlobalToLocal(user->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc)); 263f14660a4SJames Wright 26437cbb16aSJed Brown PetscCall( 26537cbb16aSJed Brown PetscSNPrintf(file_path_refined, sizeof file_path_refined, "%s/nsrefined-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no)); 266f14660a4SJames Wright 2672b730f8bSJeremy L Thompson PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined), file_path_refined, FILE_MODE_WRITE, &viewer_refined)); 268f14660a4SJames Wright PetscCall(VecView(Q_refined_loc, viewer_refined)); 269f14660a4SJames Wright PetscCall(DMRestoreLocalVector(user->dm_viz, &Q_refined_loc)); 270f14660a4SJames Wright PetscCall(DMRestoreGlobalVector(user->dm_viz, &Q_refined)); 271f14660a4SJames Wright PetscCall(PetscViewerDestroy(&viewer_refined)); 27277841947SLeila Ghaffari } 273f14660a4SJames Wright PetscCall(DMRestoreLocalVector(user->dm, &Q_loc)); 27437cbb16aSJed Brown } 27577841947SLeila Ghaffari 27677841947SLeila Ghaffari // Save data in a binary file for continuation of simulations 27769293791SJames Wright if (user->app_ctx->add_stepnum2bin) { 27837cbb16aSJed Brown PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", user->app_ctx->output_dir, step_no)); 27969293791SJames Wright } else { 2802b730f8bSJeremy L Thompson PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", user->app_ctx->output_dir)); 28169293791SJames Wright } 2822b730f8bSJeremy L Thompson PetscCall(PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, &viewer)); 283f14660a4SJames Wright 2840de6a49fSJames Wright PetscInt32 token = PetscDefined(USE_64BIT_INDICES) ? FLUIDS_FILE_TOKEN_64 : FLUIDS_FILE_TOKEN_32; 2850de6a49fSJames Wright PetscCall(PetscViewerBinaryWrite(viewer, &token, 1, PETSC_INT32)); 2864de8550aSJed Brown PetscCall(PetscViewerBinaryWrite(viewer, &step_no, 1, PETSC_INT)); 2874de8550aSJed Brown time /= user->units->second; // Dimensionalize time back 2884de8550aSJed Brown PetscCall(PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL)); 289f14660a4SJames Wright PetscCall(VecView(Q, viewer)); 290f14660a4SJames Wright PetscCall(PetscViewerDestroy(&viewer)); 291ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 292f14660a4SJames Wright } 293f14660a4SJames Wright 294ca69d878SAdeleke O. Bankole // CSV Monitor 295ca69d878SAdeleke O. Bankole PetscErrorCode TSMonitor_WallForce(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) { 296ca69d878SAdeleke O. Bankole User user = ctx; 297ca69d878SAdeleke O. Bankole Vec G_loc; 298ca69d878SAdeleke O. Bankole PetscInt num_wall = user->app_ctx->wall_forces.num_wall, dim = 3; 299ca69d878SAdeleke O. Bankole const PetscInt *walls = user->app_ctx->wall_forces.walls; 300ca69d878SAdeleke O. Bankole PetscViewer viewer = user->app_ctx->wall_forces.viewer; 301ca69d878SAdeleke O. Bankole PetscViewerFormat format = user->app_ctx->wall_forces.viewer_format; 302ca69d878SAdeleke O. Bankole PetscScalar *reaction_force; 303ca69d878SAdeleke O. Bankole PetscBool iascii; 304ca69d878SAdeleke O. Bankole 305ca69d878SAdeleke O. Bankole PetscFunctionBeginUser; 306ee4ca9cbSJames Wright if (!viewer) PetscFunctionReturn(PETSC_SUCCESS); 307ca69d878SAdeleke O. Bankole PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); 308ca69d878SAdeleke O. Bankole PetscCall(PetscMalloc1(num_wall * dim, &reaction_force)); 309ca69d878SAdeleke O. Bankole PetscCall(Surface_Forces_NS(user->dm, G_loc, num_wall, walls, reaction_force)); 310ca69d878SAdeleke O. Bankole PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); 311ca69d878SAdeleke O. Bankole 312ca69d878SAdeleke O. Bankole PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 313ca69d878SAdeleke O. Bankole 314ca69d878SAdeleke O. Bankole if (iascii) { 315ca69d878SAdeleke O. Bankole if (format == PETSC_VIEWER_ASCII_CSV && !user->app_ctx->wall_forces.header_written) { 316ca69d878SAdeleke O. Bankole PetscCall(PetscViewerASCIIPrintf(viewer, "Step,Time,Wall,ForceX,ForceY,ForceZ\n")); 317ca69d878SAdeleke O. Bankole user->app_ctx->wall_forces.header_written = PETSC_TRUE; 318ca69d878SAdeleke O. Bankole } 319ca69d878SAdeleke O. Bankole for (PetscInt w = 0; w < num_wall; w++) { 320ca69d878SAdeleke O. Bankole PetscInt wall = walls[w]; 321ca69d878SAdeleke O. Bankole if (format == PETSC_VIEWER_ASCII_CSV) { 322ca69d878SAdeleke O. Bankole PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT ",%g,%" PetscInt_FMT ",%g,%g,%g\n", step_no, time, wall, 323ca69d878SAdeleke O. Bankole reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2])); 324ca69d878SAdeleke O. Bankole 325ca69d878SAdeleke O. Bankole } else { 326ca69d878SAdeleke O. Bankole PetscCall(PetscViewerASCIIPrintf(viewer, "Wall %" PetscInt_FMT " Forces: Force_x = %12g, Force_y = %12g, Force_z = %12g\n", wall, 327ca69d878SAdeleke O. Bankole reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2])); 328ca69d878SAdeleke O. Bankole } 329ca69d878SAdeleke O. Bankole } 330ca69d878SAdeleke O. Bankole } 331ca69d878SAdeleke O. Bankole PetscCall(PetscFree(reaction_force)); 332ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 333ca69d878SAdeleke O. Bankole } 334ca69d878SAdeleke O. Bankole 335f14660a4SJames Wright // User provided TS Monitor 3362b730f8bSJeremy L Thompson PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) { 337f14660a4SJames Wright User user = ctx; 338f14660a4SJames Wright 339f17d818dSJames Wright PetscFunctionBeginUser; 34037cbb16aSJed Brown // Print every 'checkpoint_interval' steps 341894de27cSJames Wright if (user->app_ctx->checkpoint_interval <= 0 || step_no % user->app_ctx->checkpoint_interval != 0 || 34249aac155SJeremy L Thompson (user->app_ctx->cont_steps == step_no && step_no != 0)) { 343ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 34449aac155SJeremy L Thompson } 345f14660a4SJames Wright 346f14660a4SJames Wright PetscCall(WriteOutput(user, Q, step_no, time)); 347ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 34877841947SLeila Ghaffari } 34977841947SLeila Ghaffari 35077841947SLeila Ghaffari // TS: Create, setup, and solve 3512b730f8bSJeremy L Thompson PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, Vec *Q, PetscScalar *f_time, TS *ts) { 35277841947SLeila Ghaffari MPI_Comm comm = user->comm; 35377841947SLeila Ghaffari TSAdapt adapt; 35477841947SLeila Ghaffari PetscScalar final_time; 35577841947SLeila Ghaffari 356f17d818dSJames Wright PetscFunctionBeginUser; 3572b730f8bSJeremy L Thompson PetscCall(TSCreate(comm, ts)); 3582b730f8bSJeremy L Thompson PetscCall(TSSetDM(*ts, dm)); 3591a7db67cSJames Wright PetscCall(TSSetApplicationContext(*ts, user)); 36077841947SLeila Ghaffari if (phys->implicit) { 3612b730f8bSJeremy L Thompson PetscCall(TSSetType(*ts, TSBDF)); 36277841947SLeila Ghaffari if (user->op_ifunction) { 3632b730f8bSJeremy L Thompson PetscCall(TSSetIFunction(*ts, NULL, IFunction_NS, &user)); 36477841947SLeila Ghaffari } else { // Implicit integrators can fall back to using an RHSFunction 3652b730f8bSJeremy L Thompson PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user)); 36677841947SLeila Ghaffari } 367*91c97f41SJames Wright if (user->mat_ijacobian) { 3682b730f8bSJeremy L Thompson PetscCall(DMTSSetIJacobian(dm, FormIJacobian_NS, &user)); 369e334ad8fSJed Brown } 37077841947SLeila Ghaffari } else { 3719ad5e8e4SJames Wright PetscCheck(user->op_rhs_ctx, comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction"); 3722b730f8bSJeremy L Thompson PetscCall(TSSetType(*ts, TSRK)); 3732b730f8bSJeremy L Thompson PetscCall(TSRKSetType(*ts, TSRK5F)); 3742b730f8bSJeremy L Thompson PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user)); 37577841947SLeila Ghaffari } 3762b730f8bSJeremy L Thompson PetscCall(TSSetMaxTime(*ts, 500. * user->units->second)); 3772b730f8bSJeremy L Thompson PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER)); 37851b00d91SJames Wright if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(TSSetErrorIfStepFails(*ts, PETSC_FALSE)); 3792b730f8bSJeremy L Thompson PetscCall(TSSetTimeStep(*ts, 1.e-2 * user->units->second)); 3802b730f8bSJeremy L Thompson PetscCall(TSGetAdapt(*ts, &adapt)); 3812b730f8bSJeremy L Thompson PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * user->units->second, 1.e2 * user->units->second)); 3822b730f8bSJeremy L Thompson PetscCall(TSSetFromOptions(*ts)); 383*91c97f41SJames Wright if (user->mat_ijacobian) { 384*91c97f41SJames Wright if (app_ctx->amat_type && !strcmp(app_ctx->amat_type, MATSHELL)) { 385*91c97f41SJames Wright PetscBool use_matceed_pmat; 386*91c97f41SJames Wright SNES snes; 387*91c97f41SJames Wright KSP ksp; 388*91c97f41SJames Wright PC pc; 389*91c97f41SJames Wright PCType pc_type; 390*91c97f41SJames Wright Mat Pmat; 391*91c97f41SJames Wright 392*91c97f41SJames Wright PetscCall(TSGetSNES(*ts, &snes)); 393*91c97f41SJames Wright PetscCall(SNESGetKSP(snes, &ksp)); 394*91c97f41SJames Wright PetscCall(KSPGetPC(ksp, &pc)); 395*91c97f41SJames Wright PetscCall(PCGetType(pc, &pc_type)); 396*91c97f41SJames Wright PetscCall(PetscStrcmpAny(pc_type, &use_matceed_pmat, PCJACOBI, PCVPBJACOBI, PCPBJACOBI, "")); 397*91c97f41SJames Wright Pmat = use_matceed_pmat ? user->mat_ijacobian : NULL; 398*91c97f41SJames Wright PetscCall(TSSetIJacobian(*ts, user->mat_ijacobian, Pmat, NULL, NULL)); 399*91c97f41SJames Wright } 400*91c97f41SJames Wright } 40115637395SJames Wright user->time_bc_set = -1.0; // require all BCs be updated 402d55646a4SJames Wright if (app_ctx->cont_steps) { // continue from previous timestep data 40377841947SLeila Ghaffari PetscInt count; 40477841947SLeila Ghaffari PetscViewer viewer; 4052b730f8bSJeremy L Thompson 4064de8550aSJed Brown if (app_ctx->cont_time <= 0) { // Legacy files did not include step number and time 4072b730f8bSJeremy L Thompson PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_time_file, FILE_MODE_READ, &viewer)); 4084de8550aSJed Brown PetscCall(PetscViewerBinaryRead(viewer, &app_ctx->cont_time, 1, &count, PETSC_REAL)); 4092b730f8bSJeremy L Thompson PetscCall(PetscViewerDestroy(&viewer)); 4104de8550aSJed Brown PetscCheck(app_ctx->cont_steps != -1, comm, PETSC_ERR_ARG_INCOMP, 4114de8550aSJed Brown "-continue step number not specified, but checkpoint file does not contain a step number (likely written by older code version)"); 4124de8550aSJed Brown } 4134de8550aSJed Brown PetscCall(TSSetTime(*ts, app_ctx->cont_time * user->units->second)); 414d8e0aecdSJed Brown PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps)); 41577841947SLeila Ghaffari } 4168fb33541SJames Wright if (app_ctx->test_type == TESTTYPE_NONE) { 4172b730f8bSJeremy L Thompson PetscCall(TSMonitorSet(*ts, TSMonitor_NS, user, NULL)); 4188fb33541SJames Wright } 419ca69d878SAdeleke O. Bankole if (app_ctx->wall_forces.viewer) { 420ca69d878SAdeleke O. Bankole PetscCall(TSMonitorSet(*ts, TSMonitor_WallForce, user, NULL)); 421ca69d878SAdeleke O. Bankole } 422b7d66439SJames Wright if (app_ctx->turb_spanstats_enable) { 423f5452247SJames Wright PetscCall(TSMonitorSet(*ts, TSMonitor_TurbulenceStatistics, user, NULL)); 424495b9769SJames Wright CeedScalar previous_time = app_ctx->cont_time * user->units->second; 425a424bcd0SJames Wright PetscCallCeed(user->ceed, 426a424bcd0SJames Wright CeedOperatorSetContextDouble(user->spanstats.op_stats_collect_ctx->op, user->spanstats.previous_time_label, &previous_time)); 427a175e481SJames Wright } 4284e9802d1SJames Wright if (app_ctx->diff_filter_monitor) PetscCall(TSMonitorSet(*ts, TSMonitor_DifferentialFilter, user, NULL)); 42977841947SLeila Ghaffari 43028134cfaSJames Wright if (app_ctx->sgs_train_enable) { 43128134cfaSJames Wright PetscCall(TSMonitorSet(*ts, TSMonitor_SGS_DD_Training, user, NULL)); 43228134cfaSJames Wright PetscCall(TSSetPostStep(*ts, TSPostStep_SGS_DD_Training)); 43328134cfaSJames Wright } 43477841947SLeila Ghaffari // Solve 435d8e0aecdSJed Brown PetscReal start_time; 436d8e0aecdSJed Brown PetscInt start_step; 4372b730f8bSJeremy L Thompson PetscCall(TSGetTime(*ts, &start_time)); 438d8e0aecdSJed Brown PetscCall(TSGetStepNumber(*ts, &start_step)); 4397b39487dSJeremy L Thompson 440711a423aSJed Brown PetscCall(PetscLogDefaultBegin()); // So we can use PetscLogStageGetPerfInfo without -log_view 4417b39487dSJeremy L Thompson PetscPreLoadBegin(PETSC_FALSE, "Fluids Solve"); 4427b39487dSJeremy L Thompson PetscCall(TSSetTime(*ts, start_time)); 443d8e0aecdSJed Brown PetscCall(TSSetStepNumber(*ts, start_step)); 4447b39487dSJeremy L Thompson if (PetscPreLoadingOn) { 4457b39487dSJeremy L Thompson // LCOV_EXCL_START 4467b39487dSJeremy L Thompson SNES snes; 4477b39487dSJeremy L Thompson Vec Q_preload; 4487b39487dSJeremy L Thompson PetscReal rtol; 4497b39487dSJeremy L Thompson PetscCall(VecDuplicate(*Q, &Q_preload)); 4507b39487dSJeremy L Thompson PetscCall(VecCopy(*Q, Q_preload)); 4517b39487dSJeremy L Thompson PetscCall(TSGetSNES(*ts, &snes)); 4527b39487dSJeremy L Thompson PetscCall(SNESGetTolerances(snes, NULL, &rtol, NULL, NULL, NULL)); 4532b730f8bSJeremy L Thompson PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, .99, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 454fc4c3d69SJames Wright PetscCall(TSSetSolution(*ts, Q_preload)); 4557b39487dSJeremy L Thompson PetscCall(TSStep(*ts)); 4562b730f8bSJeremy L Thompson PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, rtol, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 4577b39487dSJeremy L Thompson PetscCall(VecDestroy(&Q_preload)); 4587b39487dSJeremy L Thompson // LCOV_EXCL_STOP 4597b39487dSJeremy L Thompson } else { 4602b730f8bSJeremy L Thompson PetscCall(PetscBarrier((PetscObject)*ts)); 4612b730f8bSJeremy L Thompson PetscCall(TSSolve(*ts, *Q)); 4627b39487dSJeremy L Thompson } 4637b39487dSJeremy L Thompson PetscPreLoadEnd(); 4647b39487dSJeremy L Thompson 4657b39487dSJeremy L Thompson PetscCall(TSGetSolveTime(*ts, &final_time)); 46677841947SLeila Ghaffari *f_time = final_time; 4677b39487dSJeremy L Thompson 4688fb33541SJames Wright if (app_ctx->test_type == TESTTYPE_NONE) { 469f14660a4SJames Wright PetscInt step_no; 470f14660a4SJames Wright PetscCall(TSGetStepNumber(*ts, &step_no)); 471a175e481SJames Wright if (user->app_ctx->checkpoint_interval > 0 || user->app_ctx->checkpoint_interval == -1) { 472f14660a4SJames Wright PetscCall(WriteOutput(user, *Q, step_no, final_time)); 473f14660a4SJames Wright } 474f14660a4SJames Wright 47575d1798cSJames Wright PetscLogStage stage_id; 476711a423aSJed Brown PetscEventPerfInfo stage_perf; 4777b39487dSJeremy L Thompson 4787b39487dSJeremy L Thompson PetscCall(PetscLogStageGetId("Fluids Solve", &stage_id)); 479711a423aSJed Brown PetscCall(PetscLogStageGetPerfInfo(stage_id, &stage_perf)); 480711a423aSJed Brown PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time taken for solution (sec): %g\n", stage_perf.time)); 48177841947SLeila Ghaffari } 482ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 48377841947SLeila Ghaffari } 484