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 118d78d7c8SJames Wright #include <differential_filter.h> 12149fb536SJames Wright #include <navierstokes.h> 139ae013d6SJames Wright #include <smartsim.h> 14f0d618eaSJames Wright #include <spanstats.h> 15c5e9980aSAdeleke O. Bankole #include "../qfunctions/newtonian_state.h" 16a515125bSLeila Ghaffari 174cb09fddSJames Wright // @brief Insert Boundary values if it's a new time 180c373b74SJames Wright PetscErrorCode UpdateBoundaryValues(Honee honee, Vec Q_loc, PetscReal t) { 19c996854bSJames Wright PetscFunctionBeginUser; 200c373b74SJames Wright if (honee->time_bc_set != t) { 210c373b74SJames Wright PetscCall(DMPlexInsertBoundaryValues(honee->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL)); 220c373b74SJames Wright honee->time_bc_set = t; 23c996854bSJames Wright } 24d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 25c996854bSJames Wright } 26c996854bSJames Wright 27a515125bSLeila Ghaffari // RHS (Explicit time-stepper) function setup 28a515125bSLeila Ghaffari // This is the RHS of the ODE, given as u_t = G(t,u) 29a515125bSLeila Ghaffari // This function takes in a state vector Q and writes into G 30a515125bSLeila Ghaffari PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) { 310c373b74SJames Wright Honee honee = *(Honee *)user_data; 320c373b74SJames Wright Ceed ceed = honee->ceed; 33fd969b44SJames Wright PetscScalar dt; 34dc3c760aSJames Wright Vec Q_loc = honee->Q_loc, R; 35a78efa86SJames Wright PetscMemType q_mem_type; 36a515125bSLeila Ghaffari 3706f41313SJames Wright PetscFunctionBeginUser; 38e2f84137SJeremy L Thompson // Update time dependent data 390c373b74SJames Wright PetscCall(UpdateBoundaryValues(honee, Q_loc, t)); 400c373b74SJames Wright if (honee->phys->solution_time_label) 410c373b74SJames Wright PetscCallCeed(ceed, CeedOperatorSetContextDouble(honee->op_rhs_ctx->op, honee->phys->solution_time_label, &t)); 422b916ea7SJeremy L Thompson PetscCall(TSGetTimeStep(ts, &dt)); 430c373b74SJames Wright if (honee->phys->timestep_size_label) 440c373b74SJames Wright PetscCallCeed(ceed, CeedOperatorSetContextDouble(honee->op_rhs_ctx->op, honee->phys->timestep_size_label, &dt)); 45a515125bSLeila Ghaffari 46dc3c760aSJames Wright PetscCall(DMGetNamedGlobalVector(honee->dm, "RHS Residual", &R)); 470c373b74SJames Wright PetscCall(DMGlobalToLocal(honee->dm, Q, INSERT_VALUES, Q_loc)); 480c373b74SJames Wright if (honee->app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) PetscCall(DivDiffFluxProjectionApply(honee->diff_flux_proj, Q_loc)); 49dc3c760aSJames Wright PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, R, honee->op_rhs_ctx)); 50a78efa86SJames Wright 51a78efa86SJames Wright // Inverse of the mass matrix 52dc3c760aSJames Wright PetscCall(VecReadPetscToCeed(Q_loc, &q_mem_type, honee->q_ceed)); 53bc5573aaSJames Wright { 54bc5573aaSJames Wright // Run PCApply manually if using ksp_type preonly -pc_type jacobi 55bc5573aaSJames Wright // This is to avoid an AllReduce call in KSPSolve_Preonly, which causes significant slowdowns for lumped mass matrix solves. 56bc5573aaSJames Wright // See https://gitlab.com/petsc/petsc/-/merge_requests/8048 for more details and a possible fix 57bc5573aaSJames Wright PC pc; 58bc5573aaSJames Wright PetscBool ispreonly, isjacobi; 59bc5573aaSJames Wright PetscCall(KSPGetPC(honee->mass_ksp, &pc)); 60bc5573aaSJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)honee->mass_ksp, KSPPREONLY, &ispreonly)); 61bc5573aaSJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &isjacobi)); 62bc5573aaSJames Wright if (ispreonly && isjacobi) PetscCall(PCApply(pc, R, G)); 63bc5573aaSJames Wright else PetscCall(KSPSolve(honee->mass_ksp, R, G)); 64bc5573aaSJames Wright } 650c373b74SJames Wright PetscCall(VecReadCeedToPetsc(honee->q_ceed, q_mem_type, Q_loc)); 66dc3c760aSJames Wright 67dc3c760aSJames Wright PetscCall(DMRestoreNamedGlobalVector(honee->dm, "RHS Residual", &R)); 68d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 69a515125bSLeila Ghaffari } 70a515125bSLeila Ghaffari 71c5e9980aSAdeleke O. Bankole // Surface forces function setup 72c5e9980aSAdeleke O. Bankole static PetscErrorCode Surface_Forces_NS(DM dm, Vec G_loc, PetscInt num_walls, const PetscInt walls[], PetscScalar *reaction_force) { 73c5e9980aSAdeleke O. Bankole DMLabel face_label; 7462ab3c0bSJames Wright const PetscScalar *g_array; 7562ab3c0bSJames Wright PetscInt dim = 3; 7662ab3c0bSJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)dm); 7762ab3c0bSJames Wright PetscSection section; 78c5e9980aSAdeleke O. Bankole 79c5e9980aSAdeleke O. Bankole PetscFunctionBeginUser; 80c5e9980aSAdeleke O. Bankole PetscCall(DMGetLabel(dm, "Face Sets", &face_label)); 8162ab3c0bSJames Wright PetscCall(VecGetArrayRead(G_loc, &g_array)); 82c5e9980aSAdeleke O. Bankole for (PetscInt w = 0; w < num_walls; w++) { 8362ab3c0bSJames Wright const PetscInt wall = walls[w], *points; 84c5e9980aSAdeleke O. Bankole IS wall_is; 8562ab3c0bSJames Wright PetscInt num_points, num_comp = 0; 8662ab3c0bSJames Wright 87c5e9980aSAdeleke O. Bankole PetscCall(DMLabelGetStratumIS(face_label, wall, &wall_is)); 8862ab3c0bSJames Wright if (!wall_is) continue; // No wall points on this process, skip 8962ab3c0bSJames Wright 9062ab3c0bSJames Wright PetscCall(DMGetLocalSection(dm, §ion)); 9162ab3c0bSJames Wright PetscCall(PetscSectionGetFieldComponents(section, 0, &num_comp)); 92c5e9980aSAdeleke O. Bankole PetscCall(ISGetSize(wall_is, &num_points)); 93c5e9980aSAdeleke O. Bankole PetscCall(ISGetIndices(wall_is, &points)); 94c5e9980aSAdeleke O. Bankole for (PetscInt i = 0; i < num_points; i++) { 95c5e9980aSAdeleke O. Bankole const PetscInt p = points[i]; 96c5e9980aSAdeleke O. Bankole const StateConservative *r; 9762ab3c0bSJames Wright PetscInt dof; 9862ab3c0bSJames Wright 9962ab3c0bSJames Wright PetscCall(DMPlexPointLocalRead(dm, p, g_array, &r)); 10062ab3c0bSJames Wright PetscCall(PetscSectionGetDof(section, p, &dof)); 1012004e3acSAdeleke O. Bankole for (PetscInt node = 0; node < dof / num_comp; node++) { 10262ab3c0bSJames Wright for (PetscInt j = 0; j < dim; j++) { 1032004e3acSAdeleke O. Bankole reaction_force[w * dim + j] -= r[node].momentum[j]; 1042004e3acSAdeleke O. Bankole } 105c5e9980aSAdeleke O. Bankole } 106c5e9980aSAdeleke O. Bankole } 107c5e9980aSAdeleke O. Bankole PetscCall(ISRestoreIndices(wall_is, &points)); 108c5e9980aSAdeleke O. Bankole PetscCall(ISDestroy(&wall_is)); 109c5e9980aSAdeleke O. Bankole } 110*14bd2a07SJames Wright PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, reaction_force, dim * num_walls, MPIU_SCALAR, MPI_SUM, comm)); 11162ab3c0bSJames Wright PetscCall(VecRestoreArrayRead(G_loc, &g_array)); 112d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 113c5e9980aSAdeleke O. Bankole } 114c5e9980aSAdeleke O. Bankole 115a515125bSLeila Ghaffari // Implicit time-stepper function setup 1162b916ea7SJeremy L Thompson PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) { 1170c373b74SJames Wright Honee honee = *(Honee *)user_data; 1180c373b74SJames Wright Ceed ceed = honee->ceed; 119fd969b44SJames Wright PetscScalar dt; 1200c373b74SJames Wright Vec Q_loc = honee->Q_loc, Q_dot_loc = honee->Q_dot_loc, G_loc; 121a515125bSLeila Ghaffari PetscMemType q_mem_type, q_dot_mem_type, g_mem_type; 122a515125bSLeila Ghaffari 12306f41313SJames Wright PetscFunctionBeginUser; 124dabd2275SJames Wright PetscCall(DMGlobalToLocalBegin(honee->dm, Q_dot, INSERT_VALUES, Q_dot_loc)); 1250c373b74SJames Wright PetscCall(DMGetNamedLocalVector(honee->dm, "ResidualLocal", &G_loc)); 126e2f84137SJeremy L Thompson 127e2f84137SJeremy L Thompson // Update time dependent data 1280c373b74SJames Wright PetscCall(UpdateBoundaryValues(honee, Q_loc, t)); 1290c373b74SJames Wright if (honee->phys->solution_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(honee->op_ifunction, honee->phys->solution_time_label, &t)); 1302b916ea7SJeremy L Thompson PetscCall(TSGetTimeStep(ts, &dt)); 1310c373b74SJames Wright if (honee->phys->timestep_size_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(honee->op_ifunction, honee->phys->timestep_size_label, &dt)); 132a515125bSLeila Ghaffari 133a515125bSLeila Ghaffari // Global-to-local 1340c373b74SJames Wright PetscCall(DMGlobalToLocalBegin(honee->dm, Q, INSERT_VALUES, Q_loc)); 1350c373b74SJames Wright PetscCall(DMGlobalToLocalEnd(honee->dm, Q, INSERT_VALUES, Q_loc)); 1360c373b74SJames Wright if (honee->app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) PetscCall(DivDiffFluxProjectionApply(honee->diff_flux_proj, Q_loc)); 1370c373b74SJames Wright PetscCall(DMGlobalToLocalEnd(honee->dm, Q_dot, INSERT_VALUES, Q_dot_loc)); 138a515125bSLeila Ghaffari 139a515125bSLeila Ghaffari // Place PETSc vectors in CEED vectors 1400c373b74SJames Wright PetscCall(VecReadPetscToCeed(Q_loc, &q_mem_type, honee->q_ceed)); 1410c373b74SJames Wright PetscCall(VecReadPetscToCeed(Q_dot_loc, &q_dot_mem_type, honee->q_dot_ceed)); 1420c373b74SJames Wright PetscCall(VecPetscToCeed(G_loc, &g_mem_type, honee->g_ceed)); 143a515125bSLeila Ghaffari 144a515125bSLeila Ghaffari // Apply CEED operator 145ea615d4cSJames Wright PetscCall(PetscLogEventBegin(HONEE_CeedOperatorApply, Q, G, 0, 0)); 1467eedc94cSJames Wright PetscCall(PetscLogGpuTimeBegin()); 1470c373b74SJames Wright PetscCallCeed(honee->ceed, CeedOperatorApply(honee->op_ifunction, honee->q_ceed, honee->g_ceed, CEED_REQUEST_IMMEDIATE)); 1487eedc94cSJames Wright PetscCall(PetscLogGpuTimeEnd()); 149ea615d4cSJames Wright PetscCall(PetscLogEventEnd(HONEE_CeedOperatorApply, Q, G, 0, 0)); 150a515125bSLeila Ghaffari 151a515125bSLeila Ghaffari // Restore vectors 1520c373b74SJames Wright PetscCall(VecReadCeedToPetsc(honee->q_ceed, q_mem_type, Q_loc)); 1530c373b74SJames Wright PetscCall(VecReadCeedToPetsc(honee->q_dot_ceed, q_dot_mem_type, Q_dot_loc)); 1540c373b74SJames Wright PetscCall(VecCeedToPetsc(honee->g_ceed, g_mem_type, G_loc)); 155a515125bSLeila Ghaffari 1560c373b74SJames Wright if (honee->app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) { 1570c373b74SJames Wright PetscCall(SgsDDApplyIFunction(honee, Q_loc, G_loc)); 15801ab89c1SJames Wright } 1599c678832SJames Wright 160a515125bSLeila Ghaffari // Local-to-Global 1612b916ea7SJeremy L Thompson PetscCall(VecZeroEntries(G)); 1620c373b74SJames Wright PetscCall(DMLocalToGlobal(honee->dm, G_loc, ADD_VALUES, G)); 163a515125bSLeila Ghaffari 164a515125bSLeila Ghaffari // Restore vectors 1650c373b74SJames Wright PetscCall(DMRestoreNamedLocalVector(honee->dm, "ResidualLocal", &G_loc)); 166d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 167a515125bSLeila Ghaffari } 168a515125bSLeila Ghaffari 1692b916ea7SJeremy L Thompson PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, PetscReal shift, Mat J, Mat J_pre, void *user_data) { 1700c373b74SJames Wright Honee honee = *(Honee *)user_data; 171ebfabadfSJames Wright PetscBool J_is_matceed, J_is_mffd, J_pre_is_matceed, J_pre_is_mffd; 17206f41313SJames Wright 173f0b65372SJed Brown PetscFunctionBeginUser; 17404855949SJed Brown PetscCall(PetscObjectTypeCompare((PetscObject)J, MATMFFD, &J_is_mffd)); 175ebfabadfSJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)J, MATCEED, &J_is_matceed)); 176ebfabadfSJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATMFFD, &J_pre_is_mffd)); 177ebfabadfSJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATCEED, &J_pre_is_matceed)); 178ebfabadfSJames Wright 1790c373b74SJames Wright PetscCall(MatCeedSetContextReal(honee->mat_ijacobian, "ijacobian time shift", shift)); 180ebfabadfSJames Wright 181ebfabadfSJames Wright if (J_is_matceed || J_is_mffd) { 18204855949SJed Brown PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 18304855949SJed Brown PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 1840c373b74SJames Wright } else PetscCall(MatCeedAssembleCOO(honee->mat_ijacobian, J)); 185ebfabadfSJames Wright 186ebfabadfSJames Wright if (J_pre_is_matceed && J != J_pre) { 187ebfabadfSJames Wright PetscCall(MatAssemblyBegin(J_pre, MAT_FINAL_ASSEMBLY)); 188ebfabadfSJames Wright PetscCall(MatAssemblyEnd(J_pre, MAT_FINAL_ASSEMBLY)); 189ebfabadfSJames Wright } else if (!J_pre_is_matceed && !J_pre_is_mffd && J != J_pre) { 1900c373b74SJames Wright PetscCall(MatCeedAssembleCOO(honee->mat_ijacobian, J_pre)); 19104855949SJed Brown } 192d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 193f0b65372SJed Brown } 194f0b65372SJed Brown 1950c373b74SJames Wright PetscErrorCode WriteOutput(Honee honee, Vec Q, PetscInt step_no, PetscScalar time) { 196a515125bSLeila Ghaffari Vec Q_loc; 197a515125bSLeila Ghaffari char file_path[PETSC_MAX_PATH_LEN]; 198a515125bSLeila Ghaffari PetscViewer viewer; 199a515125bSLeila Ghaffari 20006f41313SJames Wright PetscFunctionBeginUser; 2010c373b74SJames Wright if (honee->app_ctx->checkpoint_vtk) { 202a515125bSLeila Ghaffari // Set up output 2030c373b74SJames Wright PetscCall(DMGetLocalVector(honee->dm, &Q_loc)); 2047538d537SJames Wright PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec")); 2057538d537SJames Wright PetscCall(VecZeroEntries(Q_loc)); 2060c373b74SJames Wright PetscCall(DMGlobalToLocal(honee->dm, Q, INSERT_VALUES, Q_loc)); 207a515125bSLeila Ghaffari 208a515125bSLeila Ghaffari // Output 2090c373b74SJames Wright PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-%03" PetscInt_FMT ".vtu", honee->app_ctx->output_dir, step_no)); 2107538d537SJames Wright 2112b916ea7SJeremy L Thompson PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path, FILE_MODE_WRITE, &viewer)); 2127538d537SJames Wright PetscCall(VecView(Q_loc, viewer)); 2137538d537SJames Wright PetscCall(PetscViewerDestroy(&viewer)); 2140c373b74SJames Wright if (honee->dm_viz) { 215a515125bSLeila Ghaffari Vec Q_refined, Q_refined_loc; 216a515125bSLeila Ghaffari char file_path_refined[PETSC_MAX_PATH_LEN]; 217a515125bSLeila Ghaffari PetscViewer viewer_refined; 218a515125bSLeila Ghaffari 2190c373b74SJames Wright PetscCall(DMGetGlobalVector(honee->dm_viz, &Q_refined)); 2200c373b74SJames Wright PetscCall(DMGetLocalVector(honee->dm_viz, &Q_refined_loc)); 2217538d537SJames Wright PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined")); 2227538d537SJames Wright 2230c373b74SJames Wright PetscCall(MatInterpolate(honee->interp_viz, Q, Q_refined)); 2247538d537SJames Wright PetscCall(VecZeroEntries(Q_refined_loc)); 2250c373b74SJames Wright PetscCall(DMGlobalToLocal(honee->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc)); 2267538d537SJames Wright 227852e5969SJed Brown PetscCall( 2280c373b74SJames Wright PetscSNPrintf(file_path_refined, sizeof file_path_refined, "%s/nsrefined-%03" PetscInt_FMT ".vtu", honee->app_ctx->output_dir, step_no)); 2297538d537SJames Wright 2302b916ea7SJeremy L Thompson PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined), file_path_refined, FILE_MODE_WRITE, &viewer_refined)); 2317538d537SJames Wright PetscCall(VecView(Q_refined_loc, viewer_refined)); 2320c373b74SJames Wright PetscCall(DMRestoreLocalVector(honee->dm_viz, &Q_refined_loc)); 2330c373b74SJames Wright PetscCall(DMRestoreGlobalVector(honee->dm_viz, &Q_refined)); 2347538d537SJames Wright PetscCall(PetscViewerDestroy(&viewer_refined)); 235a515125bSLeila Ghaffari } 2360c373b74SJames Wright PetscCall(DMRestoreLocalVector(honee->dm, &Q_loc)); 237852e5969SJed Brown } 238a515125bSLeila Ghaffari 239a515125bSLeila Ghaffari // Save data in a binary file for continuation of simulations 2400c373b74SJames Wright if (honee->app_ctx->add_stepnum2bin) { 2410c373b74SJames Wright PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", honee->app_ctx->output_dir, step_no)); 24291a36801SJames Wright } else { 2430c373b74SJames Wright PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", honee->app_ctx->output_dir)); 24491a36801SJames Wright } 2450c373b74SJames Wright PetscCall(PetscViewerBinaryOpen(honee->comm, file_path, FILE_MODE_WRITE, &viewer)); 2467538d537SJames Wright 2470c373b74SJames Wright time /= honee->units->second; // Dimensionalize time back 248b237916aSJames Wright PetscCall(HoneeWriteBinaryVec(viewer, Q, time, step_no)); 2497538d537SJames Wright PetscCall(PetscViewerDestroy(&viewer)); 250d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 2517538d537SJames Wright } 2527538d537SJames Wright 253c5e9980aSAdeleke O. Bankole // CSV Monitor 254c5e9980aSAdeleke O. Bankole PetscErrorCode TSMonitor_WallForce(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) { 2550c373b74SJames Wright Honee honee = ctx; 256c5e9980aSAdeleke O. Bankole Vec G_loc; 2570c373b74SJames Wright PetscInt num_wall = honee->app_ctx->wall_forces.num_wall, dim = 3; 2580c373b74SJames Wright const PetscInt *walls = honee->app_ctx->wall_forces.walls; 2590c373b74SJames Wright PetscViewer viewer = honee->app_ctx->wall_forces.viewer; 2600c373b74SJames Wright PetscViewerFormat format = honee->app_ctx->wall_forces.viewer_format; 261c5e9980aSAdeleke O. Bankole PetscScalar *reaction_force; 26262ab3c0bSJames Wright PetscBool is_ascii; 263c5e9980aSAdeleke O. Bankole 264c5e9980aSAdeleke O. Bankole PetscFunctionBeginUser; 265d949ddfcSJames Wright if (!viewer) PetscFunctionReturn(PETSC_SUCCESS); 2660c373b74SJames Wright PetscCall(DMGetNamedLocalVector(honee->dm, "ResidualLocal", &G_loc)); 26762ab3c0bSJames Wright PetscCall(PetscCalloc1(num_wall * dim, &reaction_force)); 2680c373b74SJames Wright PetscCall(Surface_Forces_NS(honee->dm, G_loc, num_wall, walls, reaction_force)); 2690c373b74SJames Wright PetscCall(DMRestoreNamedLocalVector(honee->dm, "ResidualLocal", &G_loc)); 270c5e9980aSAdeleke O. Bankole 27162ab3c0bSJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &is_ascii)); 272c5e9980aSAdeleke O. Bankole 27362ab3c0bSJames Wright if (is_ascii) { 2740c373b74SJames Wright if (format == PETSC_VIEWER_ASCII_CSV && !honee->app_ctx->wall_forces.header_written) { 275c5e9980aSAdeleke O. Bankole PetscCall(PetscViewerASCIIPrintf(viewer, "Step,Time,Wall,ForceX,ForceY,ForceZ\n")); 2760c373b74SJames Wright honee->app_ctx->wall_forces.header_written = PETSC_TRUE; 277c5e9980aSAdeleke O. Bankole } 278c5e9980aSAdeleke O. Bankole for (PetscInt w = 0; w < num_wall; w++) { 279c5e9980aSAdeleke O. Bankole PetscInt wall = walls[w]; 280c5e9980aSAdeleke O. Bankole if (format == PETSC_VIEWER_ASCII_CSV) { 281c5e9980aSAdeleke O. Bankole PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT ",%g,%" PetscInt_FMT ",%g,%g,%g\n", step_no, time, wall, 282c5e9980aSAdeleke O. Bankole reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2])); 283c5e9980aSAdeleke O. Bankole 284c5e9980aSAdeleke O. Bankole } else { 285c5e9980aSAdeleke O. Bankole PetscCall(PetscViewerASCIIPrintf(viewer, "Wall %" PetscInt_FMT " Forces: Force_x = %12g, Force_y = %12g, Force_z = %12g\n", wall, 286c5e9980aSAdeleke O. Bankole reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2])); 287c5e9980aSAdeleke O. Bankole } 288c5e9980aSAdeleke O. Bankole } 289c5e9980aSAdeleke O. Bankole } 290c5e9980aSAdeleke O. Bankole PetscCall(PetscFree(reaction_force)); 291d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 292c5e9980aSAdeleke O. Bankole } 293c5e9980aSAdeleke O. Bankole 2947538d537SJames Wright // User provided TS Monitor 2952b916ea7SJeremy L Thompson PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) { 2960c373b74SJames Wright Honee honee = ctx; 2977538d537SJames Wright 29806f41313SJames Wright PetscFunctionBeginUser; 299852e5969SJed Brown // Print every 'checkpoint_interval' steps 3000c373b74SJames Wright if (honee->app_ctx->checkpoint_interval <= 0 || step_no % honee->app_ctx->checkpoint_interval != 0 || 3010c373b74SJames Wright (honee->app_ctx->cont_steps == step_no && step_no != 0)) { 302d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 303e419654dSJeremy L Thompson } 3047538d537SJames Wright 3050c373b74SJames Wright PetscCall(WriteOutput(honee, Q, step_no, time)); 306d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 307a515125bSLeila Ghaffari } 308a515125bSLeila Ghaffari 3098b774af8SJames Wright PetscErrorCode TSPostStep_CheckStep(TS ts) { 3108b774af8SJames Wright Honee honee; 3118b774af8SJames Wright PetscReal norm; 3128b774af8SJames Wright PetscInt step; 3138b774af8SJames Wright Vec Q; 3148b774af8SJames Wright 3158b774af8SJames Wright PetscFunctionBeginUser; 3168b774af8SJames Wright PetscCall(TSGetApplicationContext(ts, &honee)); 3178b774af8SJames Wright PetscCall(TSGetStepNumber(ts, &step)); 3188b774af8SJames Wright PetscCall(TSGetSolution(ts, &Q)); 3198b774af8SJames Wright if (step % honee->app_ctx->check_step_interval) PetscFunctionReturn(PETSC_SUCCESS); 3208b774af8SJames Wright PetscCall(VecNorm(Q, NORM_1, &norm)); 321827052a7SJames Wright if (PetscIsInfOrNanReal(norm)) { 3228b774af8SJames Wright PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "Solution diverged: Nans found in solution\n")); 3238b774af8SJames Wright PetscCall(TSSetConvergedReason(ts, TS_DIVERGED_NONLINEAR_SOLVE)); 3248b774af8SJames Wright } 3258b774af8SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3268b774af8SJames Wright } 3278b774af8SJames Wright 328354560d1SJames Wright PetscErrorCode TSPostStep_MaxWallTime(TS ts) { 329354560d1SJames Wright Honee honee; 330354560d1SJames Wright PetscInt step; 331354560d1SJames Wright PetscMPIInt rank; 332354560d1SJames Wright MPI_Comm comm; 333354560d1SJames Wright PetscBool is_wall_time_exceeded = PETSC_FALSE; 334354560d1SJames Wright 335354560d1SJames Wright PetscFunctionBeginUser; 336354560d1SJames Wright PetscCall(TSGetApplicationContext(ts, &honee)); 337354560d1SJames Wright PetscCall(TSGetStepNumber(ts, &step)); 338354560d1SJames Wright if (step % honee->max_wall_time_interval) PetscFunctionReturn(PETSC_SUCCESS); 339354560d1SJames Wright PetscCall(PetscObjectGetComm((PetscObject)ts, &comm)); 340354560d1SJames Wright PetscCallMPI(MPI_Comm_rank(comm, &rank)); 341354560d1SJames Wright if (rank == 0) is_wall_time_exceeded = time(NULL) > honee->max_wall_time ? PETSC_TRUE : PETSC_FALSE; 342354560d1SJames Wright // Must broadcast to avoid race condition 343354560d1SJames Wright PetscCallMPI(MPI_Bcast(&is_wall_time_exceeded, 1, MPIU_BOOL, 0, comm)); 344354560d1SJames Wright if (PetscUnlikely(is_wall_time_exceeded)) { 345354560d1SJames Wright PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "Stopping TSSolve: Set max wall time exceeded\n")); 346354560d1SJames Wright PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_USER)); 347354560d1SJames Wright } 348354560d1SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 349354560d1SJames Wright } 350354560d1SJames Wright 3514c6ae86eSJames Wright /** 3524c6ae86eSJames Wright @brief TSPostStep for HONEE 3534c6ae86eSJames Wright 3544c6ae86eSJames Wright `TSSetPostStep()` only accepts a single function argument, so this function groups together all post-step 3554c6ae86eSJames Wright functionality needed for HONEE features 3564c6ae86eSJames Wright 3574c6ae86eSJames Wright @param[in] ts TS object 3584c6ae86eSJames Wright **/ 3594c6ae86eSJames Wright PetscErrorCode TSPostStep_Honee(TS ts) { 3604c6ae86eSJames Wright Honee honee; 3614c6ae86eSJames Wright 3624c6ae86eSJames Wright PetscFunctionBeginUser; 3634c6ae86eSJames Wright PetscCall(TSGetApplicationContext(ts, &honee)); 3644c6ae86eSJames Wright if (honee->max_wall_time > 0) PetscCall(TSPostStep_MaxWallTime(ts)); 3654c6ae86eSJames Wright if (honee->app_ctx->sgs_train_enable) PetscCall(TSPostStep_SGS_DD_Training(ts)); 3664c6ae86eSJames Wright if (honee->app_ctx->check_step_interval > 0) PetscCall(TSPostStep_CheckStep(ts)); 3674c6ae86eSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3684c6ae86eSJames Wright } 3694c6ae86eSJames Wright 3704e1dcb99SJames Wright /** 3714e1dcb99SJames Wright @brief Save `TSEvaluationSolutions()` to file defined by `-ts_eval_solutions_view` 3724e1dcb99SJames Wright 3734e1dcb99SJames Wright Intended to be used with `TSSetEvaluationTimes()`/`-ts_eval_times` 3744e1dcb99SJames Wright 3754e1dcb99SJames Wright @param[in] ts TS object that has the evaluation solutions 3764e1dcb99SJames Wright **/ 3774e1dcb99SJames Wright static PetscErrorCode HoneeTSEvaluationSolutions(TS ts) { 3784e1dcb99SJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)ts); 3794e1dcb99SJames Wright const PetscReal *sol_times; 3804e1dcb99SJames Wright PetscReal orig_time; 3814e1dcb99SJames Wright PetscInt num_sols, orig_step; 3824e1dcb99SJames Wright Vec *sol_vecs; 3834e1dcb99SJames Wright DM dm; 3844e1dcb99SJames Wright PetscBool is_viewer_set; 3854e1dcb99SJames Wright PetscViewer viewer; 3864e1dcb99SJames Wright PetscViewerFormat format; 3874e1dcb99SJames Wright 3884e1dcb99SJames Wright PetscFunctionBeginUser; 3894e1dcb99SJames Wright PetscCall(TSGetEvaluationSolutions(ts, &num_sols, &sol_times, &sol_vecs)); 3904e1dcb99SJames Wright if (num_sols == 0) PetscFunctionReturn(PETSC_SUCCESS); 3914e1dcb99SJames Wright PetscCall(TSGetDM(ts, &dm)); 3924e1dcb99SJames Wright PetscCall(DMGetOutputSequenceNumber(dm, &orig_step, &orig_time)); 3934e1dcb99SJames Wright 3944e1dcb99SJames Wright const char *option = "-ts_eval_solutions_view"; 3954e1dcb99SJames Wright PetscCall(PetscOptionsCreateViewer(comm, NULL, NULL, option, &viewer, &format, &is_viewer_set)); 3964e1dcb99SJames Wright if (!is_viewer_set) { 3974e1dcb99SJames Wright PetscCall(PetscPrintf(comm, "\n\nWARNING: Viewer not set for TSEvaluationSolutions. Set %s to save this data. Now throwing it away!\n", option)); 3984e1dcb99SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3994e1dcb99SJames Wright } 4004e1dcb99SJames Wright for (PetscInt i = 0; i < num_sols; i++) { 4014e1dcb99SJames Wright PetscCall(DMSetOutputSequenceNumber(dm, i, sol_times[i])); 4024e1dcb99SJames Wright PetscCall(VecView(sol_vecs[i], viewer)); 4034e1dcb99SJames Wright } 4044e1dcb99SJames Wright PetscCall(PetscViewerPopFormat(viewer)); 4054e1dcb99SJames Wright PetscCall(PetscViewerDestroy(&viewer)); 4064e1dcb99SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 4074e1dcb99SJames Wright } 4084e1dcb99SJames Wright 409a515125bSLeila Ghaffari // TS: Create, setup, and solve 4102a9a4b51SJames Wright PetscErrorCode TSSolve_NS(DM dm, Honee honee, AppCtx app_ctx, Physics phys, ProblemData problem, Vec Q, PetscScalar *f_time, TS *ts) { 4110c373b74SJames Wright MPI_Comm comm = honee->comm; 412a515125bSLeila Ghaffari TSAdapt adapt; 413a515125bSLeila Ghaffari PetscScalar final_time; 414a515125bSLeila Ghaffari 41506f41313SJames Wright PetscFunctionBeginUser; 4162b916ea7SJeremy L Thompson PetscCall(TSCreate(comm, ts)); 4172b916ea7SJeremy L Thompson PetscCall(TSSetDM(*ts, dm)); 4180c373b74SJames Wright PetscCall(TSSetApplicationContext(*ts, honee)); 419a515125bSLeila Ghaffari if (phys->implicit) { 4202b916ea7SJeremy L Thompson PetscCall(TSSetType(*ts, TSBDF)); 4210c373b74SJames Wright if (honee->op_ifunction) { 4220c373b74SJames Wright PetscCall(TSSetIFunction(*ts, NULL, IFunction_NS, &honee)); 423a515125bSLeila Ghaffari } else { // Implicit integrators can fall back to using an RHSFunction 4240c373b74SJames Wright PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &honee)); 425a515125bSLeila Ghaffari } 4260c373b74SJames Wright if (honee->mat_ijacobian) { 4270c373b74SJames Wright PetscCall(DMTSSetIJacobian(dm, FormIJacobian_NS, &honee)); 428f0b65372SJed Brown } 429a515125bSLeila Ghaffari } else { 4300c373b74SJames Wright PetscCheck(honee->op_rhs_ctx, comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction"); 4312b916ea7SJeremy L Thompson PetscCall(TSSetType(*ts, TSRK)); 4322b916ea7SJeremy L Thompson PetscCall(TSRKSetType(*ts, TSRK5F)); 4330c373b74SJames Wright PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &honee)); 434a515125bSLeila Ghaffari } 4350c373b74SJames Wright PetscCall(TSSetMaxTime(*ts, 500. * honee->units->second)); 4362b916ea7SJeremy L Thompson PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER)); 43722387d3aSJames Wright if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(TSSetErrorIfStepFails(*ts, PETSC_FALSE)); 4380c373b74SJames Wright PetscCall(TSSetTimeStep(*ts, 1.e-2 * honee->units->second)); 4392b916ea7SJeremy L Thompson PetscCall(TSGetAdapt(*ts, &adapt)); 4400c373b74SJames Wright PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * honee->units->second, 1.e2 * honee->units->second)); 4412b916ea7SJeremy L Thompson PetscCall(TSSetFromOptions(*ts)); 4420c373b74SJames Wright if (honee->mat_ijacobian) { 443ebfabadfSJames Wright if (app_ctx->amat_type && !strcmp(app_ctx->amat_type, MATSHELL)) { 444ebfabadfSJames Wright SNES snes; 445ebfabadfSJames Wright KSP ksp; 4463170c09fSJames Wright Mat Pmat, Amat; 447ebfabadfSJames Wright 448ebfabadfSJames Wright PetscCall(TSGetSNES(*ts, &snes)); 449ebfabadfSJames Wright PetscCall(SNESGetKSP(snes, &ksp)); 4500c373b74SJames Wright PetscCall(CreateSolveOperatorsFromMatCeed(ksp, honee->mat_ijacobian, PETSC_FALSE, &Amat, &Pmat)); 4510c373b74SJames Wright PetscCall(TSSetIJacobian(*ts, honee->mat_ijacobian, Pmat, NULL, NULL)); 4523170c09fSJames Wright PetscCall(MatDestroy(&Amat)); 4533170c09fSJames Wright PetscCall(MatDestroy(&Pmat)); 454ebfabadfSJames Wright } 455ebfabadfSJames Wright } 4560c373b74SJames Wright honee->time_bc_set = -1.0; // require all BCs be updated 457c26b555cSJames Wright if (app_ctx->cont_steps) { // continue from previous timestep data 4580c373b74SJames Wright PetscCall(TSSetTime(*ts, app_ctx->cont_time * honee->units->second)); 45974a6f4ddSJed Brown PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps)); 460a515125bSLeila Ghaffari } 46116cb6b6bSJames Wright 46216cb6b6bSJames Wright PetscBool add_ksp_postsolve_residual = PETSC_FALSE; 46316cb6b6bSJames Wright PetscCall(PetscOptionsGetBool(NULL, NULL, "-ksp_post_solve_residual", &add_ksp_postsolve_residual, NULL)); 46416cb6b6bSJames Wright if (add_ksp_postsolve_residual) { 46516cb6b6bSJames Wright SNES snes; 46616cb6b6bSJames Wright KSP ksp; 46716cb6b6bSJames Wright 46816cb6b6bSJames Wright PetscCall(TSGetSNES(*ts, &snes)); 46916cb6b6bSJames Wright PetscCall(SNESGetKSP(snes, &ksp)); 47016cb6b6bSJames Wright PetscCall(KSPSetResidualHistory(ksp, NULL, PETSC_DECIDE, PETSC_TRUE)); 47116cb6b6bSJames Wright PetscCall(KSPSetPostSolve(ksp, KSPPostSolve_Honee, NULL)); 47216cb6b6bSJames Wright } 4734c6ae86eSJames Wright if (honee->set_poststep) PetscCall(TSSetPostStep(*ts, TSPostStep_Honee)); 4744c6ae86eSJames Wright if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(TSMonitorSet(*ts, TSMonitor_NS, honee, NULL)); 4754c6ae86eSJames Wright if (app_ctx->wall_forces.viewer) PetscCall(TSMonitorSet(*ts, TSMonitor_WallForce, honee, NULL)); 4768fbf1d3fSJames Wright 47778c5b8e5SJames Wright PetscOptionsBegin(comm, NULL, "HONEE TS Monitor Options", NULL); 47825125139SJames Wright PetscCall(TSMonitorSetFromOptions(*ts, "-ts_monitor_total_kinetic_energy", "Monitor total kinetic energy balance terms in the domain", NULL, 47925125139SJames Wright TSMonitor_TotalKineticEnergy, SetupMontiorTotalKineticEnergy)); 48087fd7f33SJames Wright PetscCall(TSMonitorSetFromOptions(*ts, "-ts_monitor_cfl", "Monitor element CFL statistics", NULL, TSMonitor_Cfl, SetupMontiorCfl)); 48178c5b8e5SJames Wright 48278c5b8e5SJames Wright PetscCall( 48378c5b8e5SJames Wright PetscOptionsDeprecated("-ts_monitor_turbulence_spanstats_viewer_interval", "-ts_monitor_spanstats_turbulence_interval", "HONEE 0.0", NULL)); 48478c5b8e5SJames Wright PetscCall(PetscOptionsDeprecated("-ts_monitor_turbulence_spanstats_viewer", "-ts_monitor_spanstats_turbulence", "HONEE 0.0", NULL)); 48578c5b8e5SJames Wright PetscCall(TSMonitorSetFromOptions(*ts, "-ts_monitor_spanstats_turbulence", "Setup spanwise statistics collection", NULL, 48617dc1449SJames Wright TSMonitor_SpanwiseStatisticsTurbulence, SpanwiseStatisticsSetup_Turbulence)); 487b30619f6SJames Wright PetscCall(TSMonitorSetFromOptions(*ts, "-ts_monitor_spanstats_cflpe", "Setup spanwise statistics collection of CFL and Pe", NULL, 488b30619f6SJames Wright TSMonitor_SpanwiseStatisticsCflPe, SpanwiseStatisticsSetup_CflPe)); 489cb8a476cSJames Wright PetscCall(TSMonitorSetFromOptions(*ts, "-diff_filter_monitor", "Run differential filtering for every timestep (used for testing)", NULL, 490cb8a476cSJames Wright TSMonitor_DifferentialFilter, TSMonitor_DifferentialFilterSetup)); 4919ae013d6SJames Wright PetscCall(TSMonitorSetFromOptions(*ts, "-ts_monitor_smartsim_solution", "Write solution to SmartSim database", NULL, TSMonitor_SmartSimSolution, 4929ae013d6SJames Wright TSMonitor_SmartSimSolutionSetup)); 49378c5b8e5SJames Wright PetscOptionsEnd(); 49478c5b8e5SJames Wright 4954c6ae86eSJames Wright if (app_ctx->sgs_train_enable) PetscCall(TSMonitorSet(*ts, TSMonitor_SGS_DD_Training, honee, NULL)); 496e539bbbaSJames Wright 4970c373b74SJames Wright if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(PrintRunInfo(honee, honee->phys, problem, *ts)); 49871f2ed29SJames Wright if (honee->mass_ksp) PetscCall(KSPViewFromOptions(honee->mass_ksp, NULL, "-ksp_view_pre_ts_solve")); 49971f2ed29SJames Wright 500a515125bSLeila Ghaffari // Solve 50174a6f4ddSJed Brown PetscReal start_time; 50274a6f4ddSJed Brown PetscInt start_step; 5032b916ea7SJeremy L Thompson PetscCall(TSGetTime(*ts, &start_time)); 50474a6f4ddSJed Brown PetscCall(TSGetStepNumber(*ts, &start_step)); 50591982731SJeremy L Thompson 506df4304b5SJed Brown PetscCall(PetscLogDefaultBegin()); // So we can use PetscLogStageGetPerfInfo without -log_view 507ea615d4cSJames Wright PetscPreLoadBegin(PETSC_FALSE, "HONEE Solve"); 50891982731SJeremy L Thompson PetscCall(TSSetTime(*ts, start_time)); 50974a6f4ddSJed Brown PetscCall(TSSetStepNumber(*ts, start_step)); 51091982731SJeremy L Thompson if (PetscPreLoadingOn) { 51191982731SJeremy L Thompson // LCOV_EXCL_START 51291982731SJeremy L Thompson SNES snes; 5134f25acd9SJames Wright KSP ksp; 51491982731SJeremy L Thompson Vec Q_preload; 5154f25acd9SJames Wright PetscReal rtol_snes, rtol_ksp; 5162a9a4b51SJames Wright PetscCall(VecDuplicate(Q, &Q_preload)); 5172a9a4b51SJames Wright PetscCall(VecCopy(Q, Q_preload)); 51891982731SJeremy L Thompson PetscCall(TSGetSNES(*ts, &snes)); 5194f25acd9SJames Wright PetscCall(SNESGetTolerances(snes, NULL, &rtol_snes, NULL, NULL, NULL)); 5204f25acd9SJames Wright PetscCall(SNESSetTolerances(snes, PETSC_CURRENT, .99, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT)); 5214f25acd9SJames Wright PetscCall(SNESGetKSP(snes, &ksp)); 5224f25acd9SJames Wright PetscCall(KSPGetTolerances(ksp, &rtol_ksp, NULL, NULL, NULL)); 5234f25acd9SJames Wright PetscCall(KSPSetTolerances(ksp, .99, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT)); 52422c1b34eSJames Wright PetscCall(TSSetSolution(*ts, Q_preload)); 52591982731SJeremy L Thompson PetscCall(TSStep(*ts)); 5264f25acd9SJames Wright PetscCall(SNESSetTolerances(snes, PETSC_CURRENT, rtol_snes, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT)); 5274f25acd9SJames Wright PetscCall(KSPSetTolerances(ksp, rtol_ksp, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT)); 52891982731SJeremy L Thompson PetscCall(VecDestroy(&Q_preload)); 52991982731SJeremy L Thompson // LCOV_EXCL_STOP 53091982731SJeremy L Thompson } else { 5312b916ea7SJeremy L Thompson PetscCall(PetscBarrier((PetscObject)*ts)); 5322a9a4b51SJames Wright PetscCall(TSSolve(*ts, Q)); 53391982731SJeremy L Thompson } 53491982731SJeremy L Thompson PetscPreLoadEnd(); 53591982731SJeremy L Thompson 53691982731SJeremy L Thompson PetscCall(TSGetSolveTime(*ts, &final_time)); 537a515125bSLeila Ghaffari *f_time = final_time; 53891982731SJeremy L Thompson 5390e1e9333SJames Wright if (app_ctx->test_type == TESTTYPE_NONE) { 5407538d537SJames Wright PetscInt step_no; 5414e1dcb99SJames Wright PetscLogStage stage_id; 5424e1dcb99SJames Wright PetscEventPerfInfo stage_perf; 5434e1dcb99SJames Wright 5447538d537SJames Wright PetscCall(TSGetStepNumber(*ts, &step_no)); 5450c373b74SJames Wright if (honee->app_ctx->checkpoint_interval > 0 || honee->app_ctx->checkpoint_interval == -1) { 5462a9a4b51SJames Wright PetscCall(WriteOutput(honee, Q, step_no, final_time)); 5477538d537SJames Wright } 5484e1dcb99SJames Wright PetscCall(HoneeTSEvaluationSolutions(*ts)); 54991982731SJeremy L Thompson 550ea615d4cSJames Wright PetscCall(PetscLogStageGetId("HONEE Solve", &stage_id)); 551df4304b5SJed Brown PetscCall(PetscLogStageGetPerfInfo(stage_id, &stage_perf)); 552df4304b5SJed Brown PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time taken for solution (sec): %g\n", stage_perf.time)); 553a515125bSLeila Ghaffari } 554d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 555a515125bSLeila Ghaffari } 556