1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 3a515125bSLeila Ghaffari 4a515125bSLeila Ghaffari /// @file 5a515125bSLeila Ghaffari /// Time-stepping functions for Navier-Stokes example using PETSc 6a515125bSLeila Ghaffari 7e419654dSJeremy L Thompson #include <ceed.h> 8e419654dSJeremy L Thompson #include <petscdmplex.h> 9e419654dSJeremy L Thompson #include <petscts.h> 10e419654dSJeremy L Thompson 11149fb536SJames Wright #include <navierstokes.h> 12c5e9980aSAdeleke O. Bankole #include "../qfunctions/newtonian_state.h" 13a515125bSLeila Ghaffari 144cb09fddSJames Wright // @brief Insert Boundary values if it's a new time 15*0c373b74SJames Wright PetscErrorCode UpdateBoundaryValues(Honee honee, Vec Q_loc, PetscReal t) { 16c996854bSJames Wright PetscFunctionBeginUser; 17*0c373b74SJames Wright if (honee->time_bc_set != t) { 18*0c373b74SJames Wright PetscCall(DMPlexInsertBoundaryValues(honee->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL)); 19*0c373b74SJames Wright honee->time_bc_set = t; 20c996854bSJames Wright } 21d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 22c996854bSJames Wright } 23c996854bSJames Wright 24a515125bSLeila Ghaffari // RHS (Explicit time-stepper) function setup 25a515125bSLeila Ghaffari // This is the RHS of the ODE, given as u_t = G(t,u) 26a515125bSLeila Ghaffari // This function takes in a state vector Q and writes into G 27a515125bSLeila Ghaffari PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) { 28*0c373b74SJames Wright Honee honee = *(Honee *)user_data; 29*0c373b74SJames Wright Ceed ceed = honee->ceed; 30fd969b44SJames Wright PetscScalar dt; 31*0c373b74SJames Wright Vec Q_loc = honee->Q_loc; 32a78efa86SJames Wright PetscMemType q_mem_type; 33a515125bSLeila Ghaffari 3406f41313SJames Wright PetscFunctionBeginUser; 35e2f84137SJeremy L Thompson // Update time dependent data 36*0c373b74SJames Wright PetscCall(UpdateBoundaryValues(honee, Q_loc, t)); 37*0c373b74SJames Wright if (honee->phys->solution_time_label) 38*0c373b74SJames Wright PetscCallCeed(ceed, CeedOperatorSetContextDouble(honee->op_rhs_ctx->op, honee->phys->solution_time_label, &t)); 392b916ea7SJeremy L Thompson PetscCall(TSGetTimeStep(ts, &dt)); 40*0c373b74SJames Wright if (honee->phys->timestep_size_label) 41*0c373b74SJames Wright PetscCallCeed(ceed, CeedOperatorSetContextDouble(honee->op_rhs_ctx->op, honee->phys->timestep_size_label, &dt)); 42a515125bSLeila Ghaffari 43*0c373b74SJames Wright PetscCall(DMGlobalToLocal(honee->dm, Q, INSERT_VALUES, Q_loc)); 44*0c373b74SJames Wright if (honee->app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) PetscCall(DivDiffFluxProjectionApply(honee->diff_flux_proj, Q_loc)); 45*0c373b74SJames Wright PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, G, honee->op_rhs_ctx)); 46a515125bSLeila Ghaffari 47*0c373b74SJames Wright PetscCall(VecReadPetscToCeed(Q_loc, &q_mem_type, honee->q_ceed)); 48a78efa86SJames Wright 49a78efa86SJames Wright // Inverse of the mass matrix 50*0c373b74SJames Wright PetscCall(KSPSolve(honee->mass_ksp, G, G)); 51a78efa86SJames Wright 52*0c373b74SJames Wright PetscCall(VecReadCeedToPetsc(honee->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; 5962ab3c0bSJames Wright const PetscScalar *g_array; 6062ab3c0bSJames Wright PetscInt dim = 3; 6162ab3c0bSJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)dm); 6262ab3c0bSJames Wright PetscSection section; 63c5e9980aSAdeleke O. Bankole 64c5e9980aSAdeleke O. Bankole PetscFunctionBeginUser; 65c5e9980aSAdeleke O. Bankole PetscCall(DMGetLabel(dm, "Face Sets", &face_label)); 6662ab3c0bSJames Wright PetscCall(VecGetArrayRead(G_loc, &g_array)); 67c5e9980aSAdeleke O. Bankole for (PetscInt w = 0; w < num_walls; w++) { 6862ab3c0bSJames Wright const PetscInt wall = walls[w], *points; 69c5e9980aSAdeleke O. Bankole IS wall_is; 7062ab3c0bSJames Wright PetscInt num_points, num_comp = 0; 7162ab3c0bSJames Wright 72c5e9980aSAdeleke O. Bankole PetscCall(DMLabelGetStratumIS(face_label, wall, &wall_is)); 7362ab3c0bSJames Wright if (!wall_is) continue; // No wall points on this process, skip 7462ab3c0bSJames Wright 7562ab3c0bSJames Wright PetscCall(DMGetLocalSection(dm, §ion)); 7662ab3c0bSJames Wright PetscCall(PetscSectionGetFieldComponents(section, 0, &num_comp)); 77c5e9980aSAdeleke O. Bankole PetscCall(ISGetSize(wall_is, &num_points)); 78c5e9980aSAdeleke O. Bankole PetscCall(ISGetIndices(wall_is, &points)); 79c5e9980aSAdeleke O. Bankole for (PetscInt i = 0; i < num_points; i++) { 80c5e9980aSAdeleke O. Bankole const PetscInt p = points[i]; 81c5e9980aSAdeleke O. Bankole const StateConservative *r; 8262ab3c0bSJames Wright PetscInt dof; 8362ab3c0bSJames Wright 8462ab3c0bSJames Wright PetscCall(DMPlexPointLocalRead(dm, p, g_array, &r)); 8562ab3c0bSJames Wright PetscCall(PetscSectionGetDof(section, p, &dof)); 862004e3acSAdeleke O. Bankole for (PetscInt node = 0; node < dof / num_comp; node++) { 8762ab3c0bSJames Wright for (PetscInt j = 0; j < dim; 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 PetscCall(ISDestroy(&wall_is)); 94c5e9980aSAdeleke O. Bankole } 95c5e9980aSAdeleke O. Bankole PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, reaction_force, dim * num_walls, MPIU_SCALAR, MPI_SUM, comm)); 9662ab3c0bSJames Wright PetscCall(VecRestoreArrayRead(G_loc, &g_array)); 97d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 98c5e9980aSAdeleke O. Bankole } 99c5e9980aSAdeleke O. Bankole 100a515125bSLeila Ghaffari // Implicit time-stepper function setup 1012b916ea7SJeremy L Thompson PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) { 102*0c373b74SJames Wright Honee honee = *(Honee *)user_data; 103*0c373b74SJames Wright Ceed ceed = honee->ceed; 104fd969b44SJames Wright PetscScalar dt; 105*0c373b74SJames Wright Vec Q_loc = honee->Q_loc, Q_dot_loc = honee->Q_dot_loc, G_loc; 106a515125bSLeila Ghaffari PetscMemType q_mem_type, q_dot_mem_type, g_mem_type; 107a515125bSLeila Ghaffari 10806f41313SJames Wright PetscFunctionBeginUser; 109e2f84137SJeremy L Thompson // Get local vectors 110*0c373b74SJames Wright PetscCall(DMGetNamedLocalVector(honee->dm, "ResidualLocal", &G_loc)); 111e2f84137SJeremy L Thompson 112e2f84137SJeremy L Thompson // Update time dependent data 113*0c373b74SJames Wright PetscCall(UpdateBoundaryValues(honee, Q_loc, t)); 114*0c373b74SJames Wright if (honee->phys->solution_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(honee->op_ifunction, honee->phys->solution_time_label, &t)); 1152b916ea7SJeremy L Thompson PetscCall(TSGetTimeStep(ts, &dt)); 116*0c373b74SJames Wright if (honee->phys->timestep_size_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(honee->op_ifunction, honee->phys->timestep_size_label, &dt)); 117a515125bSLeila Ghaffari 118a515125bSLeila Ghaffari // Global-to-local 119*0c373b74SJames Wright PetscCall(DMGlobalToLocalBegin(honee->dm, Q, INSERT_VALUES, Q_loc)); 120*0c373b74SJames Wright PetscCall(DMGlobalToLocalBegin(honee->dm, Q_dot, INSERT_VALUES, Q_dot_loc)); 121*0c373b74SJames Wright PetscCall(DMGlobalToLocalEnd(honee->dm, Q, INSERT_VALUES, Q_loc)); 122*0c373b74SJames Wright if (honee->app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) PetscCall(DivDiffFluxProjectionApply(honee->diff_flux_proj, Q_loc)); 123*0c373b74SJames Wright PetscCall(DMGlobalToLocalEnd(honee->dm, Q_dot, INSERT_VALUES, Q_dot_loc)); 124a515125bSLeila Ghaffari 125a515125bSLeila Ghaffari // Place PETSc vectors in CEED vectors 126*0c373b74SJames Wright PetscCall(VecReadPetscToCeed(Q_loc, &q_mem_type, honee->q_ceed)); 127*0c373b74SJames Wright PetscCall(VecReadPetscToCeed(Q_dot_loc, &q_dot_mem_type, honee->q_dot_ceed)); 128*0c373b74SJames Wright PetscCall(VecPetscToCeed(G_loc, &g_mem_type, honee->g_ceed)); 129a515125bSLeila Ghaffari 130a515125bSLeila Ghaffari // Apply CEED operator 1317eedc94cSJames Wright PetscCall(PetscLogEventBegin(FLUIDS_CeedOperatorApply, Q, G, 0, 0)); 1327eedc94cSJames Wright PetscCall(PetscLogGpuTimeBegin()); 133*0c373b74SJames Wright PetscCallCeed(honee->ceed, CeedOperatorApply(honee->op_ifunction, honee->q_ceed, honee->g_ceed, CEED_REQUEST_IMMEDIATE)); 1347eedc94cSJames Wright PetscCall(PetscLogGpuTimeEnd()); 1357eedc94cSJames Wright PetscCall(PetscLogEventEnd(FLUIDS_CeedOperatorApply, Q, G, 0, 0)); 136a515125bSLeila Ghaffari 137a515125bSLeila Ghaffari // Restore vectors 138*0c373b74SJames Wright PetscCall(VecReadCeedToPetsc(honee->q_ceed, q_mem_type, Q_loc)); 139*0c373b74SJames Wright PetscCall(VecReadCeedToPetsc(honee->q_dot_ceed, q_dot_mem_type, Q_dot_loc)); 140*0c373b74SJames Wright PetscCall(VecCeedToPetsc(honee->g_ceed, g_mem_type, G_loc)); 141a515125bSLeila Ghaffari 142*0c373b74SJames Wright if (honee->app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) { 143*0c373b74SJames Wright PetscCall(SgsDDApplyIFunction(honee, Q_loc, G_loc)); 14401ab89c1SJames Wright } 1459c678832SJames Wright 146a515125bSLeila Ghaffari // Local-to-Global 1472b916ea7SJeremy L Thompson PetscCall(VecZeroEntries(G)); 148*0c373b74SJames Wright PetscCall(DMLocalToGlobal(honee->dm, G_loc, ADD_VALUES, G)); 149a515125bSLeila Ghaffari 150a515125bSLeila Ghaffari // Restore vectors 151*0c373b74SJames Wright PetscCall(DMRestoreNamedLocalVector(honee->dm, "ResidualLocal", &G_loc)); 152d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 153a515125bSLeila Ghaffari } 154a515125bSLeila Ghaffari 1552b916ea7SJeremy L Thompson PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, PetscReal shift, Mat J, Mat J_pre, void *user_data) { 156*0c373b74SJames Wright Honee honee = *(Honee *)user_data; 157ebfabadfSJames Wright PetscBool J_is_matceed, J_is_mffd, J_pre_is_matceed, J_pre_is_mffd; 15806f41313SJames Wright 159f0b65372SJed Brown PetscFunctionBeginUser; 16004855949SJed Brown PetscCall(PetscObjectTypeCompare((PetscObject)J, MATMFFD, &J_is_mffd)); 161ebfabadfSJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)J, MATCEED, &J_is_matceed)); 162ebfabadfSJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATMFFD, &J_pre_is_mffd)); 163ebfabadfSJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATCEED, &J_pre_is_matceed)); 164ebfabadfSJames Wright 165*0c373b74SJames Wright PetscCall(MatCeedSetContextReal(honee->mat_ijacobian, "ijacobian time shift", shift)); 166ebfabadfSJames Wright 167ebfabadfSJames Wright if (J_is_matceed || J_is_mffd) { 16804855949SJed Brown PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 16904855949SJed Brown PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 170*0c373b74SJames Wright } else PetscCall(MatCeedAssembleCOO(honee->mat_ijacobian, J)); 171ebfabadfSJames Wright 172ebfabadfSJames Wright if (J_pre_is_matceed && J != J_pre) { 173ebfabadfSJames Wright PetscCall(MatAssemblyBegin(J_pre, MAT_FINAL_ASSEMBLY)); 174ebfabadfSJames Wright PetscCall(MatAssemblyEnd(J_pre, MAT_FINAL_ASSEMBLY)); 175ebfabadfSJames Wright } else if (!J_pre_is_matceed && !J_pre_is_mffd && J != J_pre) { 176*0c373b74SJames Wright PetscCall(MatCeedAssembleCOO(honee->mat_ijacobian, J_pre)); 17704855949SJed Brown } 178d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 179f0b65372SJed Brown } 180f0b65372SJed Brown 181*0c373b74SJames Wright PetscErrorCode WriteOutput(Honee honee, Vec Q, PetscInt step_no, PetscScalar time) { 182a515125bSLeila Ghaffari Vec Q_loc; 183a515125bSLeila Ghaffari char file_path[PETSC_MAX_PATH_LEN]; 184a515125bSLeila Ghaffari PetscViewer viewer; 185a515125bSLeila Ghaffari 18606f41313SJames Wright PetscFunctionBeginUser; 187*0c373b74SJames Wright if (honee->app_ctx->checkpoint_vtk) { 188a515125bSLeila Ghaffari // Set up output 189*0c373b74SJames Wright PetscCall(DMGetLocalVector(honee->dm, &Q_loc)); 1907538d537SJames Wright PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec")); 1917538d537SJames Wright PetscCall(VecZeroEntries(Q_loc)); 192*0c373b74SJames Wright PetscCall(DMGlobalToLocal(honee->dm, Q, INSERT_VALUES, Q_loc)); 193a515125bSLeila Ghaffari 194a515125bSLeila Ghaffari // Output 195*0c373b74SJames Wright PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-%03" PetscInt_FMT ".vtu", honee->app_ctx->output_dir, step_no)); 1967538d537SJames Wright 1972b916ea7SJeremy L Thompson PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path, FILE_MODE_WRITE, &viewer)); 1987538d537SJames Wright PetscCall(VecView(Q_loc, viewer)); 1997538d537SJames Wright PetscCall(PetscViewerDestroy(&viewer)); 200*0c373b74SJames Wright if (honee->dm_viz) { 201a515125bSLeila Ghaffari Vec Q_refined, Q_refined_loc; 202a515125bSLeila Ghaffari char file_path_refined[PETSC_MAX_PATH_LEN]; 203a515125bSLeila Ghaffari PetscViewer viewer_refined; 204a515125bSLeila Ghaffari 205*0c373b74SJames Wright PetscCall(DMGetGlobalVector(honee->dm_viz, &Q_refined)); 206*0c373b74SJames Wright PetscCall(DMGetLocalVector(honee->dm_viz, &Q_refined_loc)); 2077538d537SJames Wright PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined")); 2087538d537SJames Wright 209*0c373b74SJames Wright PetscCall(MatInterpolate(honee->interp_viz, Q, Q_refined)); 2107538d537SJames Wright PetscCall(VecZeroEntries(Q_refined_loc)); 211*0c373b74SJames Wright PetscCall(DMGlobalToLocal(honee->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc)); 2127538d537SJames Wright 213852e5969SJed Brown PetscCall( 214*0c373b74SJames Wright PetscSNPrintf(file_path_refined, sizeof file_path_refined, "%s/nsrefined-%03" PetscInt_FMT ".vtu", honee->app_ctx->output_dir, step_no)); 2157538d537SJames Wright 2162b916ea7SJeremy L Thompson PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined), file_path_refined, FILE_MODE_WRITE, &viewer_refined)); 2177538d537SJames Wright PetscCall(VecView(Q_refined_loc, viewer_refined)); 218*0c373b74SJames Wright PetscCall(DMRestoreLocalVector(honee->dm_viz, &Q_refined_loc)); 219*0c373b74SJames Wright PetscCall(DMRestoreGlobalVector(honee->dm_viz, &Q_refined)); 2207538d537SJames Wright PetscCall(PetscViewerDestroy(&viewer_refined)); 221a515125bSLeila Ghaffari } 222*0c373b74SJames Wright PetscCall(DMRestoreLocalVector(honee->dm, &Q_loc)); 223852e5969SJed Brown } 224a515125bSLeila Ghaffari 225a515125bSLeila Ghaffari // Save data in a binary file for continuation of simulations 226*0c373b74SJames Wright if (honee->app_ctx->add_stepnum2bin) { 227*0c373b74SJames Wright PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", honee->app_ctx->output_dir, step_no)); 22891a36801SJames Wright } else { 229*0c373b74SJames Wright PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", honee->app_ctx->output_dir)); 23091a36801SJames Wright } 231*0c373b74SJames Wright PetscCall(PetscViewerBinaryOpen(honee->comm, file_path, FILE_MODE_WRITE, &viewer)); 2327538d537SJames Wright 233*0c373b74SJames Wright time /= honee->units->second; // Dimensionalize time back 234b237916aSJames Wright PetscCall(HoneeWriteBinaryVec(viewer, Q, time, step_no)); 2357538d537SJames Wright PetscCall(PetscViewerDestroy(&viewer)); 236d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 2377538d537SJames Wright } 2387538d537SJames Wright 239c5e9980aSAdeleke O. Bankole // CSV Monitor 240c5e9980aSAdeleke O. Bankole PetscErrorCode TSMonitor_WallForce(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) { 241*0c373b74SJames Wright Honee honee = ctx; 242c5e9980aSAdeleke O. Bankole Vec G_loc; 243*0c373b74SJames Wright PetscInt num_wall = honee->app_ctx->wall_forces.num_wall, dim = 3; 244*0c373b74SJames Wright const PetscInt *walls = honee->app_ctx->wall_forces.walls; 245*0c373b74SJames Wright PetscViewer viewer = honee->app_ctx->wall_forces.viewer; 246*0c373b74SJames Wright PetscViewerFormat format = honee->app_ctx->wall_forces.viewer_format; 247c5e9980aSAdeleke O. Bankole PetscScalar *reaction_force; 24862ab3c0bSJames Wright PetscBool is_ascii; 249c5e9980aSAdeleke O. Bankole 250c5e9980aSAdeleke O. Bankole PetscFunctionBeginUser; 251d949ddfcSJames Wright if (!viewer) PetscFunctionReturn(PETSC_SUCCESS); 252*0c373b74SJames Wright PetscCall(DMGetNamedLocalVector(honee->dm, "ResidualLocal", &G_loc)); 25362ab3c0bSJames Wright PetscCall(PetscCalloc1(num_wall * dim, &reaction_force)); 254*0c373b74SJames Wright PetscCall(Surface_Forces_NS(honee->dm, G_loc, num_wall, walls, reaction_force)); 255*0c373b74SJames Wright PetscCall(DMRestoreNamedLocalVector(honee->dm, "ResidualLocal", &G_loc)); 256c5e9980aSAdeleke O. Bankole 25762ab3c0bSJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &is_ascii)); 258c5e9980aSAdeleke O. Bankole 25962ab3c0bSJames Wright if (is_ascii) { 260*0c373b74SJames Wright if (format == PETSC_VIEWER_ASCII_CSV && !honee->app_ctx->wall_forces.header_written) { 261c5e9980aSAdeleke O. Bankole PetscCall(PetscViewerASCIIPrintf(viewer, "Step,Time,Wall,ForceX,ForceY,ForceZ\n")); 262*0c373b74SJames Wright honee->app_ctx->wall_forces.header_written = PETSC_TRUE; 263c5e9980aSAdeleke O. Bankole } 264c5e9980aSAdeleke O. Bankole for (PetscInt w = 0; w < num_wall; w++) { 265c5e9980aSAdeleke O. Bankole PetscInt wall = walls[w]; 266c5e9980aSAdeleke O. Bankole if (format == PETSC_VIEWER_ASCII_CSV) { 267c5e9980aSAdeleke O. Bankole PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT ",%g,%" PetscInt_FMT ",%g,%g,%g\n", step_no, time, wall, 268c5e9980aSAdeleke O. Bankole reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2])); 269c5e9980aSAdeleke O. Bankole 270c5e9980aSAdeleke O. Bankole } else { 271c5e9980aSAdeleke O. Bankole PetscCall(PetscViewerASCIIPrintf(viewer, "Wall %" PetscInt_FMT " Forces: Force_x = %12g, Force_y = %12g, Force_z = %12g\n", wall, 272c5e9980aSAdeleke O. Bankole reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2])); 273c5e9980aSAdeleke O. Bankole } 274c5e9980aSAdeleke O. Bankole } 275c5e9980aSAdeleke O. Bankole } 276c5e9980aSAdeleke O. Bankole PetscCall(PetscFree(reaction_force)); 277d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 278c5e9980aSAdeleke O. Bankole } 279c5e9980aSAdeleke O. Bankole 2807538d537SJames Wright // User provided TS Monitor 2812b916ea7SJeremy L Thompson PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) { 282*0c373b74SJames Wright Honee honee = ctx; 2837538d537SJames Wright 28406f41313SJames Wright PetscFunctionBeginUser; 285852e5969SJed Brown // Print every 'checkpoint_interval' steps 286*0c373b74SJames Wright if (honee->app_ctx->checkpoint_interval <= 0 || step_no % honee->app_ctx->checkpoint_interval != 0 || 287*0c373b74SJames Wright (honee->app_ctx->cont_steps == step_no && step_no != 0)) { 288d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 289e419654dSJeremy L Thompson } 2907538d537SJames Wright 291*0c373b74SJames Wright PetscCall(WriteOutput(honee, Q, step_no, time)); 292d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 293a515125bSLeila Ghaffari } 294a515125bSLeila Ghaffari 295a515125bSLeila Ghaffari // TS: Create, setup, and solve 296*0c373b74SJames Wright PetscErrorCode TSSolve_NS(DM dm, Honee honee, AppCtx app_ctx, Physics phys, ProblemData problem, Vec *Q, PetscScalar *f_time, TS *ts) { 297*0c373b74SJames Wright MPI_Comm comm = honee->comm; 298a515125bSLeila Ghaffari TSAdapt adapt; 299a515125bSLeila Ghaffari PetscScalar final_time; 300a515125bSLeila Ghaffari 30106f41313SJames Wright PetscFunctionBeginUser; 3022b916ea7SJeremy L Thompson PetscCall(TSCreate(comm, ts)); 3032b916ea7SJeremy L Thompson PetscCall(TSSetDM(*ts, dm)); 304*0c373b74SJames Wright PetscCall(TSSetApplicationContext(*ts, honee)); 305a515125bSLeila Ghaffari if (phys->implicit) { 3062b916ea7SJeremy L Thompson PetscCall(TSSetType(*ts, TSBDF)); 307*0c373b74SJames Wright if (honee->op_ifunction) { 308*0c373b74SJames Wright PetscCall(TSSetIFunction(*ts, NULL, IFunction_NS, &honee)); 309a515125bSLeila Ghaffari } else { // Implicit integrators can fall back to using an RHSFunction 310*0c373b74SJames Wright PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &honee)); 311a515125bSLeila Ghaffari } 312*0c373b74SJames Wright if (honee->mat_ijacobian) { 313*0c373b74SJames Wright PetscCall(DMTSSetIJacobian(dm, FormIJacobian_NS, &honee)); 314f0b65372SJed Brown } 315a515125bSLeila Ghaffari } else { 316*0c373b74SJames Wright PetscCheck(honee->op_rhs_ctx, comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction"); 3172b916ea7SJeremy L Thompson PetscCall(TSSetType(*ts, TSRK)); 3182b916ea7SJeremy L Thompson PetscCall(TSRKSetType(*ts, TSRK5F)); 319*0c373b74SJames Wright PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &honee)); 320a515125bSLeila Ghaffari } 321*0c373b74SJames Wright PetscCall(TSSetMaxTime(*ts, 500. * honee->units->second)); 3222b916ea7SJeremy L Thompson PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER)); 32322387d3aSJames Wright if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(TSSetErrorIfStepFails(*ts, PETSC_FALSE)); 324*0c373b74SJames Wright PetscCall(TSSetTimeStep(*ts, 1.e-2 * honee->units->second)); 3252b916ea7SJeremy L Thompson PetscCall(TSGetAdapt(*ts, &adapt)); 326*0c373b74SJames Wright PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * honee->units->second, 1.e2 * honee->units->second)); 3272b916ea7SJeremy L Thompson PetscCall(TSSetFromOptions(*ts)); 328*0c373b74SJames Wright if (honee->mat_ijacobian) { 329ebfabadfSJames Wright if (app_ctx->amat_type && !strcmp(app_ctx->amat_type, MATSHELL)) { 330ebfabadfSJames Wright SNES snes; 331ebfabadfSJames Wright KSP ksp; 3323170c09fSJames Wright Mat Pmat, Amat; 333ebfabadfSJames Wright 334ebfabadfSJames Wright PetscCall(TSGetSNES(*ts, &snes)); 335ebfabadfSJames Wright PetscCall(SNESGetKSP(snes, &ksp)); 336*0c373b74SJames Wright PetscCall(CreateSolveOperatorsFromMatCeed(ksp, honee->mat_ijacobian, PETSC_FALSE, &Amat, &Pmat)); 337*0c373b74SJames Wright PetscCall(TSSetIJacobian(*ts, honee->mat_ijacobian, Pmat, NULL, NULL)); 3383170c09fSJames Wright PetscCall(MatDestroy(&Amat)); 3393170c09fSJames Wright PetscCall(MatDestroy(&Pmat)); 340ebfabadfSJames Wright } 341ebfabadfSJames Wright } 342*0c373b74SJames Wright honee->time_bc_set = -1.0; // require all BCs be updated 343c26b555cSJames Wright if (app_ctx->cont_steps) { // continue from previous timestep data 344a515125bSLeila Ghaffari PetscInt count; 345a515125bSLeila Ghaffari PetscViewer viewer; 3462b916ea7SJeremy L Thompson 3479293eaa1SJed Brown if (app_ctx->cont_time <= 0) { // Legacy files did not include step number and time 3482b916ea7SJeremy L Thompson PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_time_file, FILE_MODE_READ, &viewer)); 3499293eaa1SJed Brown PetscCall(PetscViewerBinaryRead(viewer, &app_ctx->cont_time, 1, &count, PETSC_REAL)); 3502b916ea7SJeremy L Thompson PetscCall(PetscViewerDestroy(&viewer)); 3519293eaa1SJed Brown PetscCheck(app_ctx->cont_steps != -1, comm, PETSC_ERR_ARG_INCOMP, 3529293eaa1SJed Brown "-continue step number not specified, but checkpoint file does not contain a step number (likely written by older code version)"); 3539293eaa1SJed Brown } 354*0c373b74SJames Wright PetscCall(TSSetTime(*ts, app_ctx->cont_time * honee->units->second)); 35574a6f4ddSJed Brown PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps)); 356a515125bSLeila Ghaffari } 3570e1e9333SJames Wright if (app_ctx->test_type == TESTTYPE_NONE) { 358*0c373b74SJames Wright PetscCall(TSMonitorSet(*ts, TSMonitor_NS, honee, NULL)); 3590e1e9333SJames Wright } 360c5e9980aSAdeleke O. Bankole if (app_ctx->wall_forces.viewer) { 361*0c373b74SJames Wright PetscCall(TSMonitorSet(*ts, TSMonitor_WallForce, honee, NULL)); 362c5e9980aSAdeleke O. Bankole } 363c931fa59SJames Wright if (app_ctx->turb_spanstats_enable) { 364*0c373b74SJames Wright PetscCall(TSMonitorSet(*ts, TSMonitor_TurbulenceStatistics, honee, NULL)); 365*0c373b74SJames Wright CeedScalar previous_time = app_ctx->cont_time * honee->units->second; 366*0c373b74SJames Wright PetscCallCeed(honee->ceed, 367*0c373b74SJames Wright CeedOperatorSetContextDouble(honee->spanstats.op_stats_collect_ctx->op, honee->spanstats.previous_time_label, &previous_time)); 368b0488d1fSJames Wright } 369*0c373b74SJames Wright if (app_ctx->diff_filter_monitor) PetscCall(TSMonitorSet(*ts, TSMonitor_DifferentialFilter, honee, NULL)); 370a515125bSLeila Ghaffari 3711c17f66aSJames Wright if (app_ctx->sgs_train_enable) { 372*0c373b74SJames Wright PetscCall(TSMonitorSet(*ts, TSMonitor_SGS_DD_Training, honee, NULL)); 3731c17f66aSJames Wright PetscCall(TSSetPostStep(*ts, TSPostStep_SGS_DD_Training)); 3741c17f66aSJames Wright } 375e539bbbaSJames Wright 376*0c373b74SJames Wright if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(PrintRunInfo(honee, honee->phys, problem, *ts)); 377a515125bSLeila Ghaffari // Solve 37874a6f4ddSJed Brown PetscReal start_time; 37974a6f4ddSJed Brown PetscInt start_step; 3802b916ea7SJeremy L Thompson PetscCall(TSGetTime(*ts, &start_time)); 38174a6f4ddSJed Brown PetscCall(TSGetStepNumber(*ts, &start_step)); 38291982731SJeremy L Thompson 383df4304b5SJed Brown PetscCall(PetscLogDefaultBegin()); // So we can use PetscLogStageGetPerfInfo without -log_view 38491982731SJeremy L Thompson PetscPreLoadBegin(PETSC_FALSE, "Fluids Solve"); 38591982731SJeremy L Thompson PetscCall(TSSetTime(*ts, start_time)); 38674a6f4ddSJed Brown PetscCall(TSSetStepNumber(*ts, start_step)); 38791982731SJeremy L Thompson if (PetscPreLoadingOn) { 38891982731SJeremy L Thompson // LCOV_EXCL_START 38991982731SJeremy L Thompson SNES snes; 39091982731SJeremy L Thompson Vec Q_preload; 39191982731SJeremy L Thompson PetscReal rtol; 39291982731SJeremy L Thompson PetscCall(VecDuplicate(*Q, &Q_preload)); 39391982731SJeremy L Thompson PetscCall(VecCopy(*Q, Q_preload)); 39491982731SJeremy L Thompson PetscCall(TSGetSNES(*ts, &snes)); 39591982731SJeremy L Thompson PetscCall(SNESGetTolerances(snes, NULL, &rtol, NULL, NULL, NULL)); 3962b916ea7SJeremy L Thompson PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, .99, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 39722c1b34eSJames Wright PetscCall(TSSetSolution(*ts, Q_preload)); 39891982731SJeremy L Thompson PetscCall(TSStep(*ts)); 3992b916ea7SJeremy L Thompson PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, rtol, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 40091982731SJeremy L Thompson PetscCall(VecDestroy(&Q_preload)); 40191982731SJeremy L Thompson // LCOV_EXCL_STOP 40291982731SJeremy L Thompson } else { 4032b916ea7SJeremy L Thompson PetscCall(PetscBarrier((PetscObject)*ts)); 4042b916ea7SJeremy L Thompson PetscCall(TSSolve(*ts, *Q)); 40591982731SJeremy L Thompson } 40691982731SJeremy L Thompson PetscPreLoadEnd(); 40791982731SJeremy L Thompson 40891982731SJeremy L Thompson PetscCall(TSGetSolveTime(*ts, &final_time)); 409a515125bSLeila Ghaffari *f_time = final_time; 41091982731SJeremy L Thompson 4110e1e9333SJames Wright if (app_ctx->test_type == TESTTYPE_NONE) { 4127538d537SJames Wright PetscInt step_no; 4137538d537SJames Wright PetscCall(TSGetStepNumber(*ts, &step_no)); 414*0c373b74SJames Wright if (honee->app_ctx->checkpoint_interval > 0 || honee->app_ctx->checkpoint_interval == -1) { 415*0c373b74SJames Wright PetscCall(WriteOutput(honee, *Q, step_no, final_time)); 4167538d537SJames Wright } 4177538d537SJames Wright 4187eedc94cSJames Wright PetscLogStage stage_id; 419df4304b5SJed Brown PetscEventPerfInfo stage_perf; 42091982731SJeremy L Thompson 42191982731SJeremy L Thompson PetscCall(PetscLogStageGetId("Fluids Solve", &stage_id)); 422df4304b5SJed Brown PetscCall(PetscLogStageGetPerfInfo(stage_id, &stage_perf)); 423df4304b5SJed Brown PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time taken for solution (sec): %g\n", stage_perf.time)); 424a515125bSLeila Ghaffari } 425d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 426a515125bSLeila Ghaffari } 427