1*ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2*ae2b091fSJames 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 14c996854bSJames Wright // Insert Boundary values if it's a new time 15c996854bSJames Wright PetscErrorCode UpdateBoundaryValues(User user, Vec Q_loc, PetscReal t) { 16c996854bSJames Wright PetscFunctionBeginUser; 17c996854bSJames Wright if (user->time_bc_set != t) { 18c996854bSJames Wright PetscCall(DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL)); 19c996854bSJames Wright user->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) { 28a515125bSLeila Ghaffari User user = *(User *)user_data; 29701e5830SJames Wright Ceed ceed = user->ceed; 30fd969b44SJames Wright PetscScalar dt; 31da5fe0e4SJames Wright Vec Q_loc = user->Q_loc; 32a78efa86SJames Wright PetscMemType q_mem_type; 33a515125bSLeila Ghaffari 3406f41313SJames Wright PetscFunctionBeginUser; 35e2f84137SJeremy L Thompson // Update time dependent data 36c996854bSJames Wright PetscCall(UpdateBoundaryValues(user, Q_loc, t)); 37701e5830SJames Wright if (user->phys->solution_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_rhs_ctx->op, user->phys->solution_time_label, &t)); 382b916ea7SJeremy L Thompson PetscCall(TSGetTimeStep(ts, &dt)); 39701e5830SJames Wright if (user->phys->timestep_size_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_rhs_ctx->op, user->phys->timestep_size_label, &dt)); 40a515125bSLeila Ghaffari 41da5fe0e4SJames Wright PetscCall(ApplyCeedOperatorGlobalToGlobal(Q, G, user->op_rhs_ctx)); 42a515125bSLeila Ghaffari 43a78efa86SJames Wright PetscCall(VecReadPetscToCeed(Q_loc, &q_mem_type, user->q_ceed)); 44a78efa86SJames Wright 45a78efa86SJames Wright // Inverse of the mass matrix 46f8e2d240SJames Wright PetscCall(KSPSolve(user->mass_ksp, G, G)); 47a78efa86SJames Wright 48a78efa86SJames Wright PetscCall(VecReadCeedToPetsc(user->q_ceed, q_mem_type, Q_loc)); 49d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 50a515125bSLeila Ghaffari } 51a515125bSLeila Ghaffari 52c5e9980aSAdeleke O. Bankole // Surface forces function setup 53c5e9980aSAdeleke O. Bankole static PetscErrorCode Surface_Forces_NS(DM dm, Vec G_loc, PetscInt num_walls, const PetscInt walls[], PetscScalar *reaction_force) { 54c5e9980aSAdeleke O. Bankole DMLabel face_label; 55c5e9980aSAdeleke O. Bankole const PetscScalar *g; 562004e3acSAdeleke O. Bankole PetscInt dof, dim = 3; 57c5e9980aSAdeleke O. Bankole MPI_Comm comm; 582004e3acSAdeleke O. Bankole PetscSection s; 59c5e9980aSAdeleke O. Bankole 60c5e9980aSAdeleke O. Bankole PetscFunctionBeginUser; 61c5e9980aSAdeleke O. Bankole PetscCall(PetscArrayzero(reaction_force, num_walls * dim)); 62c5e9980aSAdeleke O. Bankole PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 63c5e9980aSAdeleke O. Bankole PetscCall(DMGetLabel(dm, "Face Sets", &face_label)); 64c5e9980aSAdeleke O. Bankole PetscCall(VecGetArrayRead(G_loc, &g)); 65c5e9980aSAdeleke O. Bankole for (PetscInt w = 0; w < num_walls; w++) { 66c5e9980aSAdeleke O. Bankole const PetscInt wall = walls[w]; 67c5e9980aSAdeleke O. Bankole IS wall_is; 682004e3acSAdeleke O. Bankole PetscCall(DMGetLocalSection(dm, &s)); 69c5e9980aSAdeleke O. Bankole PetscCall(DMLabelGetStratumIS(face_label, wall, &wall_is)); 70c5e9980aSAdeleke O. Bankole if (wall_is) { // There exist such points on this process 71c5e9980aSAdeleke O. Bankole PetscInt num_points; 722004e3acSAdeleke O. Bankole PetscInt num_comp = 0; 73c5e9980aSAdeleke O. Bankole const PetscInt *points; 742004e3acSAdeleke O. Bankole PetscCall(PetscSectionGetFieldComponents(s, 0, &num_comp)); 75c5e9980aSAdeleke O. Bankole PetscCall(ISGetSize(wall_is, &num_points)); 76c5e9980aSAdeleke O. Bankole PetscCall(ISGetIndices(wall_is, &points)); 77c5e9980aSAdeleke O. Bankole for (PetscInt i = 0; i < num_points; i++) { 78c5e9980aSAdeleke O. Bankole const PetscInt p = points[i]; 79c5e9980aSAdeleke O. Bankole const StateConservative *r; 80c5e9980aSAdeleke O. Bankole PetscCall(DMPlexPointLocalRead(dm, p, g, &r)); 812004e3acSAdeleke O. Bankole PetscCall(PetscSectionGetDof(s, p, &dof)); 822004e3acSAdeleke O. Bankole for (PetscInt node = 0; node < dof / num_comp; node++) { 83c5e9980aSAdeleke O. Bankole for (PetscInt j = 0; j < 3; j++) { 842004e3acSAdeleke O. Bankole reaction_force[w * dim + j] -= r[node].momentum[j]; 852004e3acSAdeleke O. Bankole } 86c5e9980aSAdeleke O. Bankole } 87c5e9980aSAdeleke O. Bankole } 88c5e9980aSAdeleke O. Bankole PetscCall(ISRestoreIndices(wall_is, &points)); 89c5e9980aSAdeleke O. Bankole } 90c5e9980aSAdeleke O. Bankole PetscCall(ISDestroy(&wall_is)); 91c5e9980aSAdeleke O. Bankole } 92c5e9980aSAdeleke O. Bankole PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, reaction_force, dim * num_walls, MPIU_SCALAR, MPI_SUM, comm)); 93c5e9980aSAdeleke O. Bankole // Restore Vectors 94c5e9980aSAdeleke O. Bankole PetscCall(VecRestoreArrayRead(G_loc, &g)); 95d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 96c5e9980aSAdeleke O. Bankole } 97c5e9980aSAdeleke O. Bankole 98a515125bSLeila Ghaffari // Implicit time-stepper function setup 992b916ea7SJeremy L Thompson PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) { 100a515125bSLeila Ghaffari User user = *(User *)user_data; 101701e5830SJames Wright Ceed ceed = user->ceed; 102fd969b44SJames Wright PetscScalar dt; 103e2f84137SJeremy L Thompson Vec Q_loc = user->Q_loc, Q_dot_loc = user->Q_dot_loc, G_loc; 104a515125bSLeila Ghaffari PetscMemType q_mem_type, q_dot_mem_type, g_mem_type; 105a515125bSLeila Ghaffari 10606f41313SJames Wright PetscFunctionBeginUser; 107e2f84137SJeremy L Thompson // Get local vectors 108c5e9980aSAdeleke O. Bankole PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); 109e2f84137SJeremy L Thompson 110e2f84137SJeremy L Thompson // Update time dependent data 111c996854bSJames Wright PetscCall(UpdateBoundaryValues(user, Q_loc, t)); 112701e5830SJames Wright if (user->phys->solution_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_ifunction, user->phys->solution_time_label, &t)); 1132b916ea7SJeremy L Thompson PetscCall(TSGetTimeStep(ts, &dt)); 114701e5830SJames Wright if (user->phys->timestep_size_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_ifunction, user->phys->timestep_size_label, &dt)); 115a515125bSLeila Ghaffari 116a515125bSLeila Ghaffari // Global-to-local 11706108310SJames Wright PetscCall(DMGlobalToLocalBegin(user->dm, Q, INSERT_VALUES, Q_loc)); 11806108310SJames Wright PetscCall(DMGlobalToLocalBegin(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc)); 11906108310SJames Wright PetscCall(DMGlobalToLocalEnd(user->dm, Q, INSERT_VALUES, Q_loc)); 12006108310SJames Wright PetscCall(DMGlobalToLocalEnd(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc)); 121a515125bSLeila Ghaffari 122a515125bSLeila Ghaffari // Place PETSc vectors in CEED vectors 123a7dac1d5SJames Wright PetscCall(VecReadPetscToCeed(Q_loc, &q_mem_type, user->q_ceed)); 124a7dac1d5SJames Wright PetscCall(VecReadPetscToCeed(Q_dot_loc, &q_dot_mem_type, user->q_dot_ceed)); 125a7dac1d5SJames Wright PetscCall(VecPetscToCeed(G_loc, &g_mem_type, user->g_ceed)); 126a515125bSLeila Ghaffari 127a515125bSLeila Ghaffari // Apply CEED operator 1287eedc94cSJames Wright PetscCall(PetscLogEventBegin(FLUIDS_CeedOperatorApply, Q, G, 0, 0)); 1297eedc94cSJames Wright PetscCall(PetscLogGpuTimeBegin()); 130b4c37c5cSJames Wright PetscCallCeed(user->ceed, CeedOperatorApply(user->op_ifunction, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE)); 1317eedc94cSJames Wright PetscCall(PetscLogGpuTimeEnd()); 1327eedc94cSJames Wright PetscCall(PetscLogEventEnd(FLUIDS_CeedOperatorApply, Q, G, 0, 0)); 133a515125bSLeila Ghaffari 134a515125bSLeila Ghaffari // Restore vectors 135a7dac1d5SJames Wright PetscCall(VecReadCeedToPetsc(user->q_ceed, q_mem_type, Q_loc)); 136a7dac1d5SJames Wright PetscCall(VecReadCeedToPetsc(user->q_dot_ceed, q_dot_mem_type, Q_dot_loc)); 137a7dac1d5SJames Wright PetscCall(VecCeedToPetsc(user->g_ceed, g_mem_type, G_loc)); 138a515125bSLeila Ghaffari 13901ab89c1SJames Wright if (user->app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) { 140ad494f68SJames Wright PetscCall(SgsDDApplyIFunction(user, Q_loc, G_loc)); 14101ab89c1SJames Wright } 1429c678832SJames Wright 143a515125bSLeila Ghaffari // Local-to-Global 1442b916ea7SJeremy L Thompson PetscCall(VecZeroEntries(G)); 1452b916ea7SJeremy L Thompson PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G)); 146a515125bSLeila Ghaffari 147a515125bSLeila Ghaffari // Restore vectors 148c5e9980aSAdeleke O. Bankole PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); 149d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 150a515125bSLeila Ghaffari } 151a515125bSLeila Ghaffari 1522b916ea7SJeremy L Thompson PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, PetscReal shift, Mat J, Mat J_pre, void *user_data) { 153f0b65372SJed Brown User user = *(User *)user_data; 154ebfabadfSJames Wright PetscBool J_is_matceed, J_is_mffd, J_pre_is_matceed, J_pre_is_mffd; 15506f41313SJames Wright 156f0b65372SJed Brown PetscFunctionBeginUser; 15704855949SJed Brown PetscCall(PetscObjectTypeCompare((PetscObject)J, MATMFFD, &J_is_mffd)); 158ebfabadfSJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)J, MATCEED, &J_is_matceed)); 159ebfabadfSJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATMFFD, &J_pre_is_mffd)); 160ebfabadfSJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATCEED, &J_pre_is_matceed)); 161ebfabadfSJames Wright 16251bb547fSJames Wright PetscCall(MatCeedSetContextReal(user->mat_ijacobian, "ijacobian time shift", shift)); 163ebfabadfSJames Wright 164ebfabadfSJames Wright if (J_is_matceed || J_is_mffd) { 16504855949SJed Brown PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 16604855949SJed Brown PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 167ebfabadfSJames Wright } else PetscCall(MatCeedAssembleCOO(user->mat_ijacobian, J)); 168ebfabadfSJames Wright 169ebfabadfSJames Wright if (J_pre_is_matceed && J != J_pre) { 170ebfabadfSJames Wright PetscCall(MatAssemblyBegin(J_pre, MAT_FINAL_ASSEMBLY)); 171ebfabadfSJames Wright PetscCall(MatAssemblyEnd(J_pre, MAT_FINAL_ASSEMBLY)); 172ebfabadfSJames Wright } else if (!J_pre_is_matceed && !J_pre_is_mffd && J != J_pre) { 173ebfabadfSJames Wright PetscCall(MatCeedAssembleCOO(user->mat_ijacobian, J_pre)); 17404855949SJed Brown } 175d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 176f0b65372SJed Brown } 177f0b65372SJed Brown 1782b916ea7SJeremy L Thompson PetscErrorCode WriteOutput(User user, Vec Q, PetscInt step_no, PetscScalar time) { 179a515125bSLeila Ghaffari Vec Q_loc; 180a515125bSLeila Ghaffari char file_path[PETSC_MAX_PATH_LEN]; 181a515125bSLeila Ghaffari PetscViewer viewer; 182a515125bSLeila Ghaffari 18306f41313SJames Wright PetscFunctionBeginUser; 184852e5969SJed Brown if (user->app_ctx->checkpoint_vtk) { 185a515125bSLeila Ghaffari // Set up output 1867538d537SJames Wright PetscCall(DMGetLocalVector(user->dm, &Q_loc)); 1877538d537SJames Wright PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec")); 1887538d537SJames Wright PetscCall(VecZeroEntries(Q_loc)); 1897538d537SJames Wright PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc)); 190a515125bSLeila Ghaffari 191a515125bSLeila Ghaffari // Output 192852e5969SJed Brown PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no)); 1937538d537SJames Wright 1942b916ea7SJeremy L Thompson PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path, FILE_MODE_WRITE, &viewer)); 1957538d537SJames Wright PetscCall(VecView(Q_loc, viewer)); 1967538d537SJames Wright PetscCall(PetscViewerDestroy(&viewer)); 197a515125bSLeila Ghaffari if (user->dm_viz) { 198a515125bSLeila Ghaffari Vec Q_refined, Q_refined_loc; 199a515125bSLeila Ghaffari char file_path_refined[PETSC_MAX_PATH_LEN]; 200a515125bSLeila Ghaffari PetscViewer viewer_refined; 201a515125bSLeila Ghaffari 2027538d537SJames Wright PetscCall(DMGetGlobalVector(user->dm_viz, &Q_refined)); 2037538d537SJames Wright PetscCall(DMGetLocalVector(user->dm_viz, &Q_refined_loc)); 2047538d537SJames Wright PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined")); 2057538d537SJames Wright 2067538d537SJames Wright PetscCall(MatInterpolate(user->interp_viz, Q, Q_refined)); 2077538d537SJames Wright PetscCall(VecZeroEntries(Q_refined_loc)); 2082b916ea7SJeremy L Thompson PetscCall(DMGlobalToLocal(user->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc)); 2097538d537SJames Wright 210852e5969SJed Brown PetscCall( 211852e5969SJed Brown PetscSNPrintf(file_path_refined, sizeof file_path_refined, "%s/nsrefined-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no)); 2127538d537SJames Wright 2132b916ea7SJeremy L Thompson PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined), file_path_refined, FILE_MODE_WRITE, &viewer_refined)); 2147538d537SJames Wright PetscCall(VecView(Q_refined_loc, viewer_refined)); 2157538d537SJames Wright PetscCall(DMRestoreLocalVector(user->dm_viz, &Q_refined_loc)); 2167538d537SJames Wright PetscCall(DMRestoreGlobalVector(user->dm_viz, &Q_refined)); 2177538d537SJames Wright PetscCall(PetscViewerDestroy(&viewer_refined)); 218a515125bSLeila Ghaffari } 2197538d537SJames Wright PetscCall(DMRestoreLocalVector(user->dm, &Q_loc)); 220852e5969SJed Brown } 221a515125bSLeila Ghaffari 222a515125bSLeila Ghaffari // Save data in a binary file for continuation of simulations 22391a36801SJames Wright if (user->app_ctx->add_stepnum2bin) { 224852e5969SJed Brown PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", user->app_ctx->output_dir, step_no)); 22591a36801SJames Wright } else { 2262b916ea7SJeremy L Thompson PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", user->app_ctx->output_dir)); 22791a36801SJames Wright } 2282b916ea7SJeremy L Thompson PetscCall(PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, &viewer)); 2297538d537SJames Wright 230e1233009SJames Wright PetscInt32 token = PetscDefined(USE_64BIT_INDICES) ? FLUIDS_FILE_TOKEN_64 : FLUIDS_FILE_TOKEN_32; 231e1233009SJames Wright PetscCall(PetscViewerBinaryWrite(viewer, &token, 1, PETSC_INT32)); 2329293eaa1SJed Brown PetscCall(PetscViewerBinaryWrite(viewer, &step_no, 1, PETSC_INT)); 2339293eaa1SJed Brown time /= user->units->second; // Dimensionalize time back 2349293eaa1SJed Brown PetscCall(PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL)); 2357538d537SJames Wright PetscCall(VecView(Q, viewer)); 2367538d537SJames Wright PetscCall(PetscViewerDestroy(&viewer)); 237d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 2387538d537SJames Wright } 2397538d537SJames Wright 240c5e9980aSAdeleke O. Bankole // CSV Monitor 241c5e9980aSAdeleke O. Bankole PetscErrorCode TSMonitor_WallForce(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) { 242c5e9980aSAdeleke O. Bankole User user = ctx; 243c5e9980aSAdeleke O. Bankole Vec G_loc; 244c5e9980aSAdeleke O. Bankole PetscInt num_wall = user->app_ctx->wall_forces.num_wall, dim = 3; 245c5e9980aSAdeleke O. Bankole const PetscInt *walls = user->app_ctx->wall_forces.walls; 246c5e9980aSAdeleke O. Bankole PetscViewer viewer = user->app_ctx->wall_forces.viewer; 247c5e9980aSAdeleke O. Bankole PetscViewerFormat format = user->app_ctx->wall_forces.viewer_format; 248c5e9980aSAdeleke O. Bankole PetscScalar *reaction_force; 249c5e9980aSAdeleke O. Bankole PetscBool iascii; 250c5e9980aSAdeleke O. Bankole 251c5e9980aSAdeleke O. Bankole PetscFunctionBeginUser; 252d949ddfcSJames Wright if (!viewer) PetscFunctionReturn(PETSC_SUCCESS); 253c5e9980aSAdeleke O. Bankole PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); 254c5e9980aSAdeleke O. Bankole PetscCall(PetscMalloc1(num_wall * dim, &reaction_force)); 255c5e9980aSAdeleke O. Bankole PetscCall(Surface_Forces_NS(user->dm, G_loc, num_wall, walls, reaction_force)); 256c5e9980aSAdeleke O. Bankole PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); 257c5e9980aSAdeleke O. Bankole 258c5e9980aSAdeleke O. Bankole PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 259c5e9980aSAdeleke O. Bankole 260c5e9980aSAdeleke O. Bankole if (iascii) { 261c5e9980aSAdeleke O. Bankole if (format == PETSC_VIEWER_ASCII_CSV && !user->app_ctx->wall_forces.header_written) { 262c5e9980aSAdeleke O. Bankole PetscCall(PetscViewerASCIIPrintf(viewer, "Step,Time,Wall,ForceX,ForceY,ForceZ\n")); 263c5e9980aSAdeleke O. Bankole user->app_ctx->wall_forces.header_written = PETSC_TRUE; 264c5e9980aSAdeleke O. Bankole } 265c5e9980aSAdeleke O. Bankole for (PetscInt w = 0; w < num_wall; w++) { 266c5e9980aSAdeleke O. Bankole PetscInt wall = walls[w]; 267c5e9980aSAdeleke O. Bankole if (format == PETSC_VIEWER_ASCII_CSV) { 268c5e9980aSAdeleke O. Bankole PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT ",%g,%" PetscInt_FMT ",%g,%g,%g\n", step_no, time, wall, 269c5e9980aSAdeleke O. Bankole reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2])); 270c5e9980aSAdeleke O. Bankole 271c5e9980aSAdeleke O. Bankole } else { 272c5e9980aSAdeleke O. Bankole PetscCall(PetscViewerASCIIPrintf(viewer, "Wall %" PetscInt_FMT " Forces: Force_x = %12g, Force_y = %12g, Force_z = %12g\n", wall, 273c5e9980aSAdeleke O. Bankole reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2])); 274c5e9980aSAdeleke O. Bankole } 275c5e9980aSAdeleke O. Bankole } 276c5e9980aSAdeleke O. Bankole } 277c5e9980aSAdeleke O. Bankole PetscCall(PetscFree(reaction_force)); 278d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 279c5e9980aSAdeleke O. Bankole } 280c5e9980aSAdeleke O. Bankole 2817538d537SJames Wright // User provided TS Monitor 2822b916ea7SJeremy L Thompson PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) { 2837538d537SJames Wright User user = ctx; 2847538d537SJames Wright 28506f41313SJames Wright PetscFunctionBeginUser; 286852e5969SJed Brown // Print every 'checkpoint_interval' steps 287c539088bSJames Wright if (user->app_ctx->checkpoint_interval <= 0 || step_no % user->app_ctx->checkpoint_interval != 0 || 288e419654dSJeremy L Thompson (user->app_ctx->cont_steps == step_no && step_no != 0)) { 289d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 290e419654dSJeremy L Thompson } 2917538d537SJames Wright 2927538d537SJames Wright PetscCall(WriteOutput(user, Q, step_no, time)); 293d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 294a515125bSLeila Ghaffari } 295a515125bSLeila Ghaffari 296a515125bSLeila Ghaffari // TS: Create, setup, and solve 297e539bbbaSJames Wright PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, ProblemData problem, Vec *Q, PetscScalar *f_time, TS *ts) { 298a515125bSLeila Ghaffari MPI_Comm comm = user->comm; 299a515125bSLeila Ghaffari TSAdapt adapt; 300a515125bSLeila Ghaffari PetscScalar final_time; 301a515125bSLeila Ghaffari 30206f41313SJames Wright PetscFunctionBeginUser; 3032b916ea7SJeremy L Thompson PetscCall(TSCreate(comm, ts)); 3042b916ea7SJeremy L Thompson PetscCall(TSSetDM(*ts, dm)); 305632a41e1SJames Wright PetscCall(TSSetApplicationContext(*ts, user)); 306a515125bSLeila Ghaffari if (phys->implicit) { 3072b916ea7SJeremy L Thompson PetscCall(TSSetType(*ts, TSBDF)); 308a515125bSLeila Ghaffari if (user->op_ifunction) { 3092b916ea7SJeremy L Thompson PetscCall(TSSetIFunction(*ts, NULL, IFunction_NS, &user)); 310a515125bSLeila Ghaffari } else { // Implicit integrators can fall back to using an RHSFunction 3112b916ea7SJeremy L Thompson PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user)); 312a515125bSLeila Ghaffari } 313ebfabadfSJames Wright if (user->mat_ijacobian) { 3142b916ea7SJeremy L Thompson PetscCall(DMTSSetIJacobian(dm, FormIJacobian_NS, &user)); 315f0b65372SJed Brown } 316a515125bSLeila Ghaffari } else { 317da5fe0e4SJames Wright PetscCheck(user->op_rhs_ctx, comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction"); 3182b916ea7SJeremy L Thompson PetscCall(TSSetType(*ts, TSRK)); 3192b916ea7SJeremy L Thompson PetscCall(TSRKSetType(*ts, TSRK5F)); 3202b916ea7SJeremy L Thompson PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user)); 321a515125bSLeila Ghaffari } 3222b916ea7SJeremy L Thompson PetscCall(TSSetMaxTime(*ts, 500. * user->units->second)); 3232b916ea7SJeremy L Thompson PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER)); 32422387d3aSJames Wright if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(TSSetErrorIfStepFails(*ts, PETSC_FALSE)); 3252b916ea7SJeremy L Thompson PetscCall(TSSetTimeStep(*ts, 1.e-2 * user->units->second)); 3262b916ea7SJeremy L Thompson PetscCall(TSGetAdapt(*ts, &adapt)); 3272b916ea7SJeremy L Thompson PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * user->units->second, 1.e2 * user->units->second)); 3282b916ea7SJeremy L Thompson PetscCall(TSSetFromOptions(*ts)); 329ebfabadfSJames Wright if (user->mat_ijacobian) { 330ebfabadfSJames Wright if (app_ctx->amat_type && !strcmp(app_ctx->amat_type, MATSHELL)) { 331ebfabadfSJames Wright SNES snes; 332ebfabadfSJames Wright KSP ksp; 3333170c09fSJames Wright Mat Pmat, Amat; 334ebfabadfSJames Wright 335ebfabadfSJames Wright PetscCall(TSGetSNES(*ts, &snes)); 336ebfabadfSJames Wright PetscCall(SNESGetKSP(snes, &ksp)); 3373170c09fSJames Wright PetscCall(CreateSolveOperatorsFromMatCeed(ksp, user->mat_ijacobian, PETSC_FALSE, &Amat, &Pmat)); 338ebfabadfSJames Wright PetscCall(TSSetIJacobian(*ts, user->mat_ijacobian, Pmat, NULL, NULL)); 3393170c09fSJames Wright PetscCall(MatDestroy(&Amat)); 3403170c09fSJames Wright PetscCall(MatDestroy(&Pmat)); 341ebfabadfSJames Wright } 342ebfabadfSJames Wright } 34391f639d2SJames Wright user->time_bc_set = -1.0; // require all BCs be updated 344c26b555cSJames Wright if (app_ctx->cont_steps) { // continue from previous timestep data 345a515125bSLeila Ghaffari PetscInt count; 346a515125bSLeila Ghaffari PetscViewer viewer; 3472b916ea7SJeremy L Thompson 3489293eaa1SJed Brown if (app_ctx->cont_time <= 0) { // Legacy files did not include step number and time 3492b916ea7SJeremy L Thompson PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_time_file, FILE_MODE_READ, &viewer)); 3509293eaa1SJed Brown PetscCall(PetscViewerBinaryRead(viewer, &app_ctx->cont_time, 1, &count, PETSC_REAL)); 3512b916ea7SJeremy L Thompson PetscCall(PetscViewerDestroy(&viewer)); 3529293eaa1SJed Brown PetscCheck(app_ctx->cont_steps != -1, comm, PETSC_ERR_ARG_INCOMP, 3539293eaa1SJed Brown "-continue step number not specified, but checkpoint file does not contain a step number (likely written by older code version)"); 3549293eaa1SJed Brown } 3559293eaa1SJed Brown PetscCall(TSSetTime(*ts, app_ctx->cont_time * user->units->second)); 35674a6f4ddSJed Brown PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps)); 357a515125bSLeila Ghaffari } 3580e1e9333SJames Wright if (app_ctx->test_type == TESTTYPE_NONE) { 3592b916ea7SJeremy L Thompson PetscCall(TSMonitorSet(*ts, TSMonitor_NS, user, NULL)); 3600e1e9333SJames Wright } 361c5e9980aSAdeleke O. Bankole if (app_ctx->wall_forces.viewer) { 362c5e9980aSAdeleke O. Bankole PetscCall(TSMonitorSet(*ts, TSMonitor_WallForce, user, NULL)); 363c5e9980aSAdeleke O. Bankole } 364c931fa59SJames Wright if (app_ctx->turb_spanstats_enable) { 36591933550SJames Wright PetscCall(TSMonitorSet(*ts, TSMonitor_TurbulenceStatistics, user, NULL)); 366b8daee98SJames Wright CeedScalar previous_time = app_ctx->cont_time * user->units->second; 367b4c37c5cSJames Wright PetscCallCeed(user->ceed, 368b4c37c5cSJames Wright CeedOperatorSetContextDouble(user->spanstats.op_stats_collect_ctx->op, user->spanstats.previous_time_label, &previous_time)); 369b0488d1fSJames Wright } 37088b07121SJames Wright if (app_ctx->diff_filter_monitor) PetscCall(TSMonitorSet(*ts, TSMonitor_DifferentialFilter, user, NULL)); 371a515125bSLeila Ghaffari 3721c17f66aSJames Wright if (app_ctx->sgs_train_enable) { 3731c17f66aSJames Wright PetscCall(TSMonitorSet(*ts, TSMonitor_SGS_DD_Training, user, NULL)); 3741c17f66aSJames Wright PetscCall(TSSetPostStep(*ts, TSPostStep_SGS_DD_Training)); 3751c17f66aSJames Wright } 376e539bbbaSJames Wright 3773678fae3SJames Wright if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(PrintRunInfo(user, user->phys, problem, *ts)); 378a515125bSLeila Ghaffari // Solve 37974a6f4ddSJed Brown PetscReal start_time; 38074a6f4ddSJed Brown PetscInt start_step; 3812b916ea7SJeremy L Thompson PetscCall(TSGetTime(*ts, &start_time)); 38274a6f4ddSJed Brown PetscCall(TSGetStepNumber(*ts, &start_step)); 38391982731SJeremy L Thompson 384df4304b5SJed Brown PetscCall(PetscLogDefaultBegin()); // So we can use PetscLogStageGetPerfInfo without -log_view 38591982731SJeremy L Thompson PetscPreLoadBegin(PETSC_FALSE, "Fluids Solve"); 38691982731SJeremy L Thompson PetscCall(TSSetTime(*ts, start_time)); 38774a6f4ddSJed Brown PetscCall(TSSetStepNumber(*ts, start_step)); 38891982731SJeremy L Thompson if (PetscPreLoadingOn) { 38991982731SJeremy L Thompson // LCOV_EXCL_START 39091982731SJeremy L Thompson SNES snes; 39191982731SJeremy L Thompson Vec Q_preload; 39291982731SJeremy L Thompson PetscReal rtol; 39391982731SJeremy L Thompson PetscCall(VecDuplicate(*Q, &Q_preload)); 39491982731SJeremy L Thompson PetscCall(VecCopy(*Q, Q_preload)); 39591982731SJeremy L Thompson PetscCall(TSGetSNES(*ts, &snes)); 39691982731SJeremy L Thompson PetscCall(SNESGetTolerances(snes, NULL, &rtol, NULL, NULL, NULL)); 3972b916ea7SJeremy L Thompson PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, .99, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 39822c1b34eSJames Wright PetscCall(TSSetSolution(*ts, Q_preload)); 39991982731SJeremy L Thompson PetscCall(TSStep(*ts)); 4002b916ea7SJeremy L Thompson PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, rtol, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 40191982731SJeremy L Thompson PetscCall(VecDestroy(&Q_preload)); 40291982731SJeremy L Thompson // LCOV_EXCL_STOP 40391982731SJeremy L Thompson } else { 4042b916ea7SJeremy L Thompson PetscCall(PetscBarrier((PetscObject)*ts)); 4052b916ea7SJeremy L Thompson PetscCall(TSSolve(*ts, *Q)); 40691982731SJeremy L Thompson } 40791982731SJeremy L Thompson PetscPreLoadEnd(); 40891982731SJeremy L Thompson 40991982731SJeremy L Thompson PetscCall(TSGetSolveTime(*ts, &final_time)); 410a515125bSLeila Ghaffari *f_time = final_time; 41191982731SJeremy L Thompson 4120e1e9333SJames Wright if (app_ctx->test_type == TESTTYPE_NONE) { 4137538d537SJames Wright PetscInt step_no; 4147538d537SJames Wright PetscCall(TSGetStepNumber(*ts, &step_no)); 415b0488d1fSJames Wright if (user->app_ctx->checkpoint_interval > 0 || user->app_ctx->checkpoint_interval == -1) { 4167538d537SJames Wright PetscCall(WriteOutput(user, *Q, step_no, final_time)); 4177538d537SJames Wright } 4187538d537SJames Wright 4197eedc94cSJames Wright PetscLogStage stage_id; 420df4304b5SJed Brown PetscEventPerfInfo stage_perf; 42191982731SJeremy L Thompson 42291982731SJeremy L Thompson PetscCall(PetscLogStageGetId("Fluids Solve", &stage_id)); 423df4304b5SJed Brown PetscCall(PetscLogStageGetPerfInfo(stage_id, &stage_perf)); 424df4304b5SJed Brown PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time taken for solution (sec): %g\n", stage_perf.time)); 425a515125bSLeila Ghaffari } 426d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 427a515125bSLeila Ghaffari } 428