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 5ea615d4cSJames Wright /// Time-stepping functions for HONEE 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 150c373b74SJames Wright PetscErrorCode UpdateBoundaryValues(Honee honee, Vec Q_loc, PetscReal t) { 16c996854bSJames Wright PetscFunctionBeginUser; 170c373b74SJames Wright if (honee->time_bc_set != t) { 180c373b74SJames Wright PetscCall(DMPlexInsertBoundaryValues(honee->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL)); 190c373b74SJames 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) { 280c373b74SJames Wright Honee honee = *(Honee *)user_data; 290c373b74SJames Wright Ceed ceed = honee->ceed; 30fd969b44SJames Wright PetscScalar dt; 31dc3c760aSJames Wright Vec Q_loc = honee->Q_loc, R; 32a78efa86SJames Wright PetscMemType q_mem_type; 33a515125bSLeila Ghaffari 3406f41313SJames Wright PetscFunctionBeginUser; 35e2f84137SJeremy L Thompson // Update time dependent data 360c373b74SJames Wright PetscCall(UpdateBoundaryValues(honee, Q_loc, t)); 370c373b74SJames Wright if (honee->phys->solution_time_label) 380c373b74SJames Wright PetscCallCeed(ceed, CeedOperatorSetContextDouble(honee->op_rhs_ctx->op, honee->phys->solution_time_label, &t)); 392b916ea7SJeremy L Thompson PetscCall(TSGetTimeStep(ts, &dt)); 400c373b74SJames Wright if (honee->phys->timestep_size_label) 410c373b74SJames Wright PetscCallCeed(ceed, CeedOperatorSetContextDouble(honee->op_rhs_ctx->op, honee->phys->timestep_size_label, &dt)); 42a515125bSLeila Ghaffari 43dc3c760aSJames Wright PetscCall(DMGetNamedGlobalVector(honee->dm, "RHS Residual", &R)); 440c373b74SJames Wright PetscCall(DMGlobalToLocal(honee->dm, Q, INSERT_VALUES, Q_loc)); 450c373b74SJames Wright if (honee->app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) PetscCall(DivDiffFluxProjectionApply(honee->diff_flux_proj, Q_loc)); 46dc3c760aSJames Wright PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, R, honee->op_rhs_ctx)); 47a78efa86SJames Wright 48a78efa86SJames Wright // Inverse of the mass matrix 49dc3c760aSJames Wright PetscCall(VecReadPetscToCeed(Q_loc, &q_mem_type, honee->q_ceed)); 50dc3c760aSJames Wright PetscCall(KSPSolve(honee->mass_ksp, R, G)); 510c373b74SJames Wright PetscCall(VecReadCeedToPetsc(honee->q_ceed, q_mem_type, Q_loc)); 52dc3c760aSJames Wright 53dc3c760aSJames Wright PetscCall(DMRestoreNamedGlobalVector(honee->dm, "RHS Residual", &R)); 54d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 55a515125bSLeila Ghaffari } 56a515125bSLeila Ghaffari 57c5e9980aSAdeleke O. Bankole // Surface forces function setup 58c5e9980aSAdeleke O. Bankole static PetscErrorCode Surface_Forces_NS(DM dm, Vec G_loc, PetscInt num_walls, const PetscInt walls[], PetscScalar *reaction_force) { 59c5e9980aSAdeleke O. Bankole DMLabel face_label; 6062ab3c0bSJames Wright const PetscScalar *g_array; 6162ab3c0bSJames Wright PetscInt dim = 3; 6262ab3c0bSJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)dm); 6362ab3c0bSJames Wright PetscSection section; 64c5e9980aSAdeleke O. Bankole 65c5e9980aSAdeleke O. Bankole PetscFunctionBeginUser; 66c5e9980aSAdeleke O. Bankole PetscCall(DMGetLabel(dm, "Face Sets", &face_label)); 6762ab3c0bSJames Wright PetscCall(VecGetArrayRead(G_loc, &g_array)); 68c5e9980aSAdeleke O. Bankole for (PetscInt w = 0; w < num_walls; w++) { 6962ab3c0bSJames Wright const PetscInt wall = walls[w], *points; 70c5e9980aSAdeleke O. Bankole IS wall_is; 7162ab3c0bSJames Wright PetscInt num_points, num_comp = 0; 7262ab3c0bSJames Wright 73c5e9980aSAdeleke O. Bankole PetscCall(DMLabelGetStratumIS(face_label, wall, &wall_is)); 7462ab3c0bSJames Wright if (!wall_is) continue; // No wall points on this process, skip 7562ab3c0bSJames Wright 7662ab3c0bSJames Wright PetscCall(DMGetLocalSection(dm, §ion)); 7762ab3c0bSJames Wright PetscCall(PetscSectionGetFieldComponents(section, 0, &num_comp)); 78c5e9980aSAdeleke O. Bankole PetscCall(ISGetSize(wall_is, &num_points)); 79c5e9980aSAdeleke O. Bankole PetscCall(ISGetIndices(wall_is, &points)); 80c5e9980aSAdeleke O. Bankole for (PetscInt i = 0; i < num_points; i++) { 81c5e9980aSAdeleke O. Bankole const PetscInt p = points[i]; 82c5e9980aSAdeleke O. Bankole const StateConservative *r; 8362ab3c0bSJames Wright PetscInt dof; 8462ab3c0bSJames Wright 8562ab3c0bSJames Wright PetscCall(DMPlexPointLocalRead(dm, p, g_array, &r)); 8662ab3c0bSJames Wright PetscCall(PetscSectionGetDof(section, p, &dof)); 872004e3acSAdeleke O. Bankole for (PetscInt node = 0; node < dof / num_comp; node++) { 8862ab3c0bSJames Wright for (PetscInt j = 0; j < dim; j++) { 892004e3acSAdeleke O. Bankole reaction_force[w * dim + j] -= r[node].momentum[j]; 902004e3acSAdeleke O. Bankole } 91c5e9980aSAdeleke O. Bankole } 92c5e9980aSAdeleke O. Bankole } 93c5e9980aSAdeleke O. Bankole PetscCall(ISRestoreIndices(wall_is, &points)); 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)); 9762ab3c0bSJames Wright PetscCall(VecRestoreArrayRead(G_loc, &g_array)); 98d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 99c5e9980aSAdeleke O. Bankole } 100c5e9980aSAdeleke O. Bankole 101a515125bSLeila Ghaffari // Implicit time-stepper function setup 1022b916ea7SJeremy L Thompson PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) { 1030c373b74SJames Wright Honee honee = *(Honee *)user_data; 1040c373b74SJames Wright Ceed ceed = honee->ceed; 105fd969b44SJames Wright PetscScalar dt; 1060c373b74SJames Wright Vec Q_loc = honee->Q_loc, Q_dot_loc = honee->Q_dot_loc, G_loc; 107a515125bSLeila Ghaffari PetscMemType q_mem_type, q_dot_mem_type, g_mem_type; 108a515125bSLeila Ghaffari 10906f41313SJames Wright PetscFunctionBeginUser; 110*dabd2275SJames Wright PetscCall(DMGlobalToLocalBegin(honee->dm, Q_dot, INSERT_VALUES, Q_dot_loc)); 1110c373b74SJames Wright PetscCall(DMGetNamedLocalVector(honee->dm, "ResidualLocal", &G_loc)); 112e2f84137SJeremy L Thompson 113e2f84137SJeremy L Thompson // Update time dependent data 1140c373b74SJames Wright PetscCall(UpdateBoundaryValues(honee, Q_loc, t)); 1150c373b74SJames Wright if (honee->phys->solution_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(honee->op_ifunction, honee->phys->solution_time_label, &t)); 1162b916ea7SJeremy L Thompson PetscCall(TSGetTimeStep(ts, &dt)); 1170c373b74SJames Wright if (honee->phys->timestep_size_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(honee->op_ifunction, honee->phys->timestep_size_label, &dt)); 118a515125bSLeila Ghaffari 119a515125bSLeila Ghaffari // Global-to-local 1200c373b74SJames Wright PetscCall(DMGlobalToLocalBegin(honee->dm, Q, INSERT_VALUES, Q_loc)); 1210c373b74SJames Wright PetscCall(DMGlobalToLocalEnd(honee->dm, Q, INSERT_VALUES, Q_loc)); 1220c373b74SJames Wright if (honee->app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) PetscCall(DivDiffFluxProjectionApply(honee->diff_flux_proj, Q_loc)); 1230c373b74SJames Wright PetscCall(DMGlobalToLocalEnd(honee->dm, Q_dot, INSERT_VALUES, Q_dot_loc)); 124a515125bSLeila Ghaffari 125a515125bSLeila Ghaffari // Place PETSc vectors in CEED vectors 1260c373b74SJames Wright PetscCall(VecReadPetscToCeed(Q_loc, &q_mem_type, honee->q_ceed)); 1270c373b74SJames Wright PetscCall(VecReadPetscToCeed(Q_dot_loc, &q_dot_mem_type, honee->q_dot_ceed)); 1280c373b74SJames Wright PetscCall(VecPetscToCeed(G_loc, &g_mem_type, honee->g_ceed)); 129a515125bSLeila Ghaffari 130a515125bSLeila Ghaffari // Apply CEED operator 131ea615d4cSJames Wright PetscCall(PetscLogEventBegin(HONEE_CeedOperatorApply, Q, G, 0, 0)); 1327eedc94cSJames Wright PetscCall(PetscLogGpuTimeBegin()); 1330c373b74SJames Wright PetscCallCeed(honee->ceed, CeedOperatorApply(honee->op_ifunction, honee->q_ceed, honee->g_ceed, CEED_REQUEST_IMMEDIATE)); 1347eedc94cSJames Wright PetscCall(PetscLogGpuTimeEnd()); 135ea615d4cSJames Wright PetscCall(PetscLogEventEnd(HONEE_CeedOperatorApply, Q, G, 0, 0)); 136a515125bSLeila Ghaffari 137a515125bSLeila Ghaffari // Restore vectors 1380c373b74SJames Wright PetscCall(VecReadCeedToPetsc(honee->q_ceed, q_mem_type, Q_loc)); 1390c373b74SJames Wright PetscCall(VecReadCeedToPetsc(honee->q_dot_ceed, q_dot_mem_type, Q_dot_loc)); 1400c373b74SJames Wright PetscCall(VecCeedToPetsc(honee->g_ceed, g_mem_type, G_loc)); 141a515125bSLeila Ghaffari 1420c373b74SJames Wright if (honee->app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) { 1430c373b74SJames Wright PetscCall(SgsDDApplyIFunction(honee, Q_loc, G_loc)); 14401ab89c1SJames Wright } 1459c678832SJames Wright 146a515125bSLeila Ghaffari // Local-to-Global 1472b916ea7SJeremy L Thompson PetscCall(VecZeroEntries(G)); 1480c373b74SJames Wright PetscCall(DMLocalToGlobal(honee->dm, G_loc, ADD_VALUES, G)); 149a515125bSLeila Ghaffari 150a515125bSLeila Ghaffari // Restore vectors 1510c373b74SJames 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) { 1560c373b74SJames 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 1650c373b74SJames 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)); 1700c373b74SJames 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) { 1760c373b74SJames Wright PetscCall(MatCeedAssembleCOO(honee->mat_ijacobian, J_pre)); 17704855949SJed Brown } 178d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 179f0b65372SJed Brown } 180f0b65372SJed Brown 1810c373b74SJames 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; 1870c373b74SJames Wright if (honee->app_ctx->checkpoint_vtk) { 188a515125bSLeila Ghaffari // Set up output 1890c373b74SJames Wright PetscCall(DMGetLocalVector(honee->dm, &Q_loc)); 1907538d537SJames Wright PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec")); 1917538d537SJames Wright PetscCall(VecZeroEntries(Q_loc)); 1920c373b74SJames Wright PetscCall(DMGlobalToLocal(honee->dm, Q, INSERT_VALUES, Q_loc)); 193a515125bSLeila Ghaffari 194a515125bSLeila Ghaffari // Output 1950c373b74SJames 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)); 2000c373b74SJames 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 2050c373b74SJames Wright PetscCall(DMGetGlobalVector(honee->dm_viz, &Q_refined)); 2060c373b74SJames Wright PetscCall(DMGetLocalVector(honee->dm_viz, &Q_refined_loc)); 2077538d537SJames Wright PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined")); 2087538d537SJames Wright 2090c373b74SJames Wright PetscCall(MatInterpolate(honee->interp_viz, Q, Q_refined)); 2107538d537SJames Wright PetscCall(VecZeroEntries(Q_refined_loc)); 2110c373b74SJames Wright PetscCall(DMGlobalToLocal(honee->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc)); 2127538d537SJames Wright 213852e5969SJed Brown PetscCall( 2140c373b74SJames 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)); 2180c373b74SJames Wright PetscCall(DMRestoreLocalVector(honee->dm_viz, &Q_refined_loc)); 2190c373b74SJames Wright PetscCall(DMRestoreGlobalVector(honee->dm_viz, &Q_refined)); 2207538d537SJames Wright PetscCall(PetscViewerDestroy(&viewer_refined)); 221a515125bSLeila Ghaffari } 2220c373b74SJames Wright PetscCall(DMRestoreLocalVector(honee->dm, &Q_loc)); 223852e5969SJed Brown } 224a515125bSLeila Ghaffari 225a515125bSLeila Ghaffari // Save data in a binary file for continuation of simulations 2260c373b74SJames Wright if (honee->app_ctx->add_stepnum2bin) { 2270c373b74SJames Wright PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", honee->app_ctx->output_dir, step_no)); 22891a36801SJames Wright } else { 2290c373b74SJames Wright PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", honee->app_ctx->output_dir)); 23091a36801SJames Wright } 2310c373b74SJames Wright PetscCall(PetscViewerBinaryOpen(honee->comm, file_path, FILE_MODE_WRITE, &viewer)); 2327538d537SJames Wright 2330c373b74SJames 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) { 2410c373b74SJames Wright Honee honee = ctx; 242c5e9980aSAdeleke O. Bankole Vec G_loc; 2430c373b74SJames Wright PetscInt num_wall = honee->app_ctx->wall_forces.num_wall, dim = 3; 2440c373b74SJames Wright const PetscInt *walls = honee->app_ctx->wall_forces.walls; 2450c373b74SJames Wright PetscViewer viewer = honee->app_ctx->wall_forces.viewer; 2460c373b74SJames 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); 2520c373b74SJames Wright PetscCall(DMGetNamedLocalVector(honee->dm, "ResidualLocal", &G_loc)); 25362ab3c0bSJames Wright PetscCall(PetscCalloc1(num_wall * dim, &reaction_force)); 2540c373b74SJames Wright PetscCall(Surface_Forces_NS(honee->dm, G_loc, num_wall, walls, reaction_force)); 2550c373b74SJames 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) { 2600c373b74SJames 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")); 2620c373b74SJames 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) { 2820c373b74SJames Wright Honee honee = ctx; 2837538d537SJames Wright 28406f41313SJames Wright PetscFunctionBeginUser; 285852e5969SJed Brown // Print every 'checkpoint_interval' steps 2860c373b74SJames Wright if (honee->app_ctx->checkpoint_interval <= 0 || step_no % honee->app_ctx->checkpoint_interval != 0 || 2870c373b74SJames Wright (honee->app_ctx->cont_steps == step_no && step_no != 0)) { 288d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 289e419654dSJeremy L Thompson } 2907538d537SJames Wright 2910c373b74SJames Wright PetscCall(WriteOutput(honee, Q, step_no, time)); 292d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 293a515125bSLeila Ghaffari } 294a515125bSLeila Ghaffari 2958b774af8SJames Wright PetscErrorCode TSPostStep_CheckStep(TS ts) { 2968b774af8SJames Wright Honee honee; 2978b774af8SJames Wright PetscReal norm; 2988b774af8SJames Wright PetscInt step; 2998b774af8SJames Wright Vec Q; 3008b774af8SJames Wright 3018b774af8SJames Wright PetscFunctionBeginUser; 3028b774af8SJames Wright PetscCall(TSGetApplicationContext(ts, &honee)); 3038b774af8SJames Wright PetscCall(TSGetStepNumber(ts, &step)); 3048b774af8SJames Wright PetscCall(TSGetSolution(ts, &Q)); 3058b774af8SJames Wright if (step % honee->app_ctx->check_step_interval) PetscFunctionReturn(PETSC_SUCCESS); 3068b774af8SJames Wright PetscCall(VecNorm(Q, NORM_1, &norm)); 307827052a7SJames Wright if (PetscIsInfOrNanReal(norm)) { 3088b774af8SJames Wright PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "Solution diverged: Nans found in solution\n")); 3098b774af8SJames Wright PetscCall(TSSetConvergedReason(ts, TS_DIVERGED_NONLINEAR_SOLVE)); 3108b774af8SJames Wright } 3118b774af8SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3128b774af8SJames Wright } 3138b774af8SJames Wright 314354560d1SJames Wright PetscErrorCode TSPostStep_MaxWallTime(TS ts) { 315354560d1SJames Wright Honee honee; 316354560d1SJames Wright PetscInt step; 317354560d1SJames Wright PetscMPIInt rank; 318354560d1SJames Wright MPI_Comm comm; 319354560d1SJames Wright PetscBool is_wall_time_exceeded = PETSC_FALSE; 320354560d1SJames Wright 321354560d1SJames Wright PetscFunctionBeginUser; 322354560d1SJames Wright PetscCall(TSGetApplicationContext(ts, &honee)); 323354560d1SJames Wright PetscCall(TSGetStepNumber(ts, &step)); 324354560d1SJames Wright if (step % honee->max_wall_time_interval) PetscFunctionReturn(PETSC_SUCCESS); 325354560d1SJames Wright PetscCall(PetscObjectGetComm((PetscObject)ts, &comm)); 326354560d1SJames Wright PetscCallMPI(MPI_Comm_rank(comm, &rank)); 327354560d1SJames Wright if (rank == 0) is_wall_time_exceeded = time(NULL) > honee->max_wall_time ? PETSC_TRUE : PETSC_FALSE; 328354560d1SJames Wright // Must broadcast to avoid race condition 329354560d1SJames Wright PetscCallMPI(MPI_Bcast(&is_wall_time_exceeded, 1, MPIU_BOOL, 0, comm)); 330354560d1SJames Wright if (PetscUnlikely(is_wall_time_exceeded)) { 331354560d1SJames Wright PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "Stopping TSSolve: Set max wall time exceeded\n")); 332354560d1SJames Wright PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_USER)); 333354560d1SJames Wright } 334354560d1SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 335354560d1SJames Wright } 336354560d1SJames Wright 3374c6ae86eSJames Wright /** 3384c6ae86eSJames Wright @brief TSPostStep for HONEE 3394c6ae86eSJames Wright 3404c6ae86eSJames Wright `TSSetPostStep()` only accepts a single function argument, so this function groups together all post-step 3414c6ae86eSJames Wright functionality needed for HONEE features 3424c6ae86eSJames Wright 3434c6ae86eSJames Wright @param[in] ts TS object 3444c6ae86eSJames Wright **/ 3454c6ae86eSJames Wright PetscErrorCode TSPostStep_Honee(TS ts) { 3464c6ae86eSJames Wright Honee honee; 3474c6ae86eSJames Wright 3484c6ae86eSJames Wright PetscFunctionBeginUser; 3494c6ae86eSJames Wright PetscCall(TSGetApplicationContext(ts, &honee)); 3504c6ae86eSJames Wright if (honee->max_wall_time > 0) PetscCall(TSPostStep_MaxWallTime(ts)); 3514c6ae86eSJames Wright if (honee->app_ctx->sgs_train_enable) PetscCall(TSPostStep_SGS_DD_Training(ts)); 3524c6ae86eSJames Wright if (honee->app_ctx->check_step_interval > 0) PetscCall(TSPostStep_CheckStep(ts)); 3534c6ae86eSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3544c6ae86eSJames Wright } 3554c6ae86eSJames Wright 356a515125bSLeila Ghaffari // TS: Create, setup, and solve 3570c373b74SJames Wright PetscErrorCode TSSolve_NS(DM dm, Honee honee, AppCtx app_ctx, Physics phys, ProblemData problem, Vec *Q, PetscScalar *f_time, TS *ts) { 3580c373b74SJames Wright MPI_Comm comm = honee->comm; 359a515125bSLeila Ghaffari TSAdapt adapt; 360a515125bSLeila Ghaffari PetscScalar final_time; 361a515125bSLeila Ghaffari 36206f41313SJames Wright PetscFunctionBeginUser; 3632b916ea7SJeremy L Thompson PetscCall(TSCreate(comm, ts)); 3642b916ea7SJeremy L Thompson PetscCall(TSSetDM(*ts, dm)); 3650c373b74SJames Wright PetscCall(TSSetApplicationContext(*ts, honee)); 366a515125bSLeila Ghaffari if (phys->implicit) { 3672b916ea7SJeremy L Thompson PetscCall(TSSetType(*ts, TSBDF)); 3680c373b74SJames Wright if (honee->op_ifunction) { 3690c373b74SJames Wright PetscCall(TSSetIFunction(*ts, NULL, IFunction_NS, &honee)); 370a515125bSLeila Ghaffari } else { // Implicit integrators can fall back to using an RHSFunction 3710c373b74SJames Wright PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &honee)); 372a515125bSLeila Ghaffari } 3730c373b74SJames Wright if (honee->mat_ijacobian) { 3740c373b74SJames Wright PetscCall(DMTSSetIJacobian(dm, FormIJacobian_NS, &honee)); 375f0b65372SJed Brown } 376a515125bSLeila Ghaffari } else { 3770c373b74SJames Wright PetscCheck(honee->op_rhs_ctx, comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction"); 3782b916ea7SJeremy L Thompson PetscCall(TSSetType(*ts, TSRK)); 3792b916ea7SJeremy L Thompson PetscCall(TSRKSetType(*ts, TSRK5F)); 3800c373b74SJames Wright PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &honee)); 381a515125bSLeila Ghaffari } 3820c373b74SJames Wright PetscCall(TSSetMaxTime(*ts, 500. * honee->units->second)); 3832b916ea7SJeremy L Thompson PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER)); 38422387d3aSJames Wright if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(TSSetErrorIfStepFails(*ts, PETSC_FALSE)); 3850c373b74SJames Wright PetscCall(TSSetTimeStep(*ts, 1.e-2 * honee->units->second)); 3862b916ea7SJeremy L Thompson PetscCall(TSGetAdapt(*ts, &adapt)); 3870c373b74SJames Wright PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * honee->units->second, 1.e2 * honee->units->second)); 3882b916ea7SJeremy L Thompson PetscCall(TSSetFromOptions(*ts)); 3890c373b74SJames Wright if (honee->mat_ijacobian) { 390ebfabadfSJames Wright if (app_ctx->amat_type && !strcmp(app_ctx->amat_type, MATSHELL)) { 391ebfabadfSJames Wright SNES snes; 392ebfabadfSJames Wright KSP ksp; 3933170c09fSJames Wright Mat Pmat, Amat; 394ebfabadfSJames Wright 395ebfabadfSJames Wright PetscCall(TSGetSNES(*ts, &snes)); 396ebfabadfSJames Wright PetscCall(SNESGetKSP(snes, &ksp)); 3970c373b74SJames Wright PetscCall(CreateSolveOperatorsFromMatCeed(ksp, honee->mat_ijacobian, PETSC_FALSE, &Amat, &Pmat)); 3980c373b74SJames Wright PetscCall(TSSetIJacobian(*ts, honee->mat_ijacobian, Pmat, NULL, NULL)); 3993170c09fSJames Wright PetscCall(MatDestroy(&Amat)); 4003170c09fSJames Wright PetscCall(MatDestroy(&Pmat)); 401ebfabadfSJames Wright } 402ebfabadfSJames Wright } 4030c373b74SJames Wright honee->time_bc_set = -1.0; // require all BCs be updated 404c26b555cSJames Wright if (app_ctx->cont_steps) { // continue from previous timestep data 4050c373b74SJames Wright PetscCall(TSSetTime(*ts, app_ctx->cont_time * honee->units->second)); 40674a6f4ddSJed Brown PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps)); 407a515125bSLeila Ghaffari } 4084c6ae86eSJames Wright if (honee->set_poststep) PetscCall(TSSetPostStep(*ts, TSPostStep_Honee)); 4094c6ae86eSJames Wright if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(TSMonitorSet(*ts, TSMonitor_NS, honee, NULL)); 4104c6ae86eSJames Wright if (app_ctx->wall_forces.viewer) PetscCall(TSMonitorSet(*ts, TSMonitor_WallForce, honee, NULL)); 411c931fa59SJames Wright if (app_ctx->turb_spanstats_enable) { 4120c373b74SJames Wright PetscCall(TSMonitorSet(*ts, TSMonitor_TurbulenceStatistics, honee, NULL)); 4130c373b74SJames Wright CeedScalar previous_time = app_ctx->cont_time * honee->units->second; 4140c373b74SJames Wright PetscCallCeed(honee->ceed, 4150c373b74SJames Wright CeedOperatorSetContextDouble(honee->spanstats.op_stats_collect_ctx->op, honee->spanstats.previous_time_label, &previous_time)); 416b0488d1fSJames Wright } 41725125139SJames Wright PetscCall(TSMonitorSetFromOptions(*ts, "-ts_monitor_total_kinetic_energy", "Monitor total kinetic energy balance terms in the domain", NULL, 41825125139SJames Wright TSMonitor_TotalKineticEnergy, SetupMontiorTotalKineticEnergy)); 41987fd7f33SJames Wright PetscCall(TSMonitorSetFromOptions(*ts, "-ts_monitor_cfl", "Monitor element CFL statistics", NULL, TSMonitor_Cfl, SetupMontiorCfl)); 4200c373b74SJames Wright if (app_ctx->diff_filter_monitor) PetscCall(TSMonitorSet(*ts, TSMonitor_DifferentialFilter, honee, NULL)); 4214c6ae86eSJames Wright if (app_ctx->sgs_train_enable) PetscCall(TSMonitorSet(*ts, TSMonitor_SGS_DD_Training, honee, NULL)); 422e539bbbaSJames Wright 4230c373b74SJames Wright if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(PrintRunInfo(honee, honee->phys, problem, *ts)); 424a515125bSLeila Ghaffari // Solve 42574a6f4ddSJed Brown PetscReal start_time; 42674a6f4ddSJed Brown PetscInt start_step; 4272b916ea7SJeremy L Thompson PetscCall(TSGetTime(*ts, &start_time)); 42874a6f4ddSJed Brown PetscCall(TSGetStepNumber(*ts, &start_step)); 42991982731SJeremy L Thompson 430df4304b5SJed Brown PetscCall(PetscLogDefaultBegin()); // So we can use PetscLogStageGetPerfInfo without -log_view 431ea615d4cSJames Wright PetscPreLoadBegin(PETSC_FALSE, "HONEE Solve"); 43291982731SJeremy L Thompson PetscCall(TSSetTime(*ts, start_time)); 43374a6f4ddSJed Brown PetscCall(TSSetStepNumber(*ts, start_step)); 43491982731SJeremy L Thompson if (PetscPreLoadingOn) { 43591982731SJeremy L Thompson // LCOV_EXCL_START 43691982731SJeremy L Thompson SNES snes; 4374f25acd9SJames Wright KSP ksp; 43891982731SJeremy L Thompson Vec Q_preload; 4394f25acd9SJames Wright PetscReal rtol_snes, rtol_ksp; 44091982731SJeremy L Thompson PetscCall(VecDuplicate(*Q, &Q_preload)); 44191982731SJeremy L Thompson PetscCall(VecCopy(*Q, Q_preload)); 44291982731SJeremy L Thompson PetscCall(TSGetSNES(*ts, &snes)); 4434f25acd9SJames Wright PetscCall(SNESGetTolerances(snes, NULL, &rtol_snes, NULL, NULL, NULL)); 4444f25acd9SJames Wright PetscCall(SNESSetTolerances(snes, PETSC_CURRENT, .99, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT)); 4454f25acd9SJames Wright PetscCall(SNESGetKSP(snes, &ksp)); 4464f25acd9SJames Wright PetscCall(KSPGetTolerances(ksp, &rtol_ksp, NULL, NULL, NULL)); 4474f25acd9SJames Wright PetscCall(KSPSetTolerances(ksp, .99, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT)); 44822c1b34eSJames Wright PetscCall(TSSetSolution(*ts, Q_preload)); 44991982731SJeremy L Thompson PetscCall(TSStep(*ts)); 4504f25acd9SJames Wright PetscCall(SNESSetTolerances(snes, PETSC_CURRENT, rtol_snes, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT)); 4514f25acd9SJames Wright PetscCall(KSPSetTolerances(ksp, rtol_ksp, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT)); 45291982731SJeremy L Thompson PetscCall(VecDestroy(&Q_preload)); 45391982731SJeremy L Thompson // LCOV_EXCL_STOP 45491982731SJeremy L Thompson } else { 4552b916ea7SJeremy L Thompson PetscCall(PetscBarrier((PetscObject)*ts)); 4562b916ea7SJeremy L Thompson PetscCall(TSSolve(*ts, *Q)); 45791982731SJeremy L Thompson } 45891982731SJeremy L Thompson PetscPreLoadEnd(); 45991982731SJeremy L Thompson 46091982731SJeremy L Thompson PetscCall(TSGetSolveTime(*ts, &final_time)); 461a515125bSLeila Ghaffari *f_time = final_time; 46291982731SJeremy L Thompson 4630e1e9333SJames Wright if (app_ctx->test_type == TESTTYPE_NONE) { 4647538d537SJames Wright PetscInt step_no; 4657538d537SJames Wright PetscCall(TSGetStepNumber(*ts, &step_no)); 4660c373b74SJames Wright if (honee->app_ctx->checkpoint_interval > 0 || honee->app_ctx->checkpoint_interval == -1) { 4670c373b74SJames Wright PetscCall(WriteOutput(honee, *Q, step_no, final_time)); 4687538d537SJames Wright } 4697538d537SJames Wright 4707eedc94cSJames Wright PetscLogStage stage_id; 471df4304b5SJed Brown PetscEventPerfInfo stage_perf; 47291982731SJeremy L Thompson 473ea615d4cSJames Wright PetscCall(PetscLogStageGetId("HONEE Solve", &stage_id)); 474df4304b5SJed Brown PetscCall(PetscLogStageGetPerfInfo(stage_id, &stage_perf)); 475df4304b5SJed Brown PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time taken for solution (sec): %g\n", stage_perf.time)); 476a515125bSLeila Ghaffari } 477d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 478a515125bSLeila Ghaffari } 479