1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 3 4 /// @file 5 /// Time-stepping functions for HONEE 6 7 #include <ceed.h> 8 #include <petscdmplex.h> 9 #include <petscts.h> 10 11 #include <navierstokes.h> 12 #include "../qfunctions/newtonian_state.h" 13 14 // @brief Insert Boundary values if it's a new time 15 PetscErrorCode UpdateBoundaryValues(Honee honee, Vec Q_loc, PetscReal t) { 16 PetscFunctionBeginUser; 17 if (honee->time_bc_set != t) { 18 PetscCall(DMPlexInsertBoundaryValues(honee->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL)); 19 honee->time_bc_set = t; 20 } 21 PetscFunctionReturn(PETSC_SUCCESS); 22 } 23 24 // RHS (Explicit time-stepper) function setup 25 // This is the RHS of the ODE, given as u_t = G(t,u) 26 // This function takes in a state vector Q and writes into G 27 PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) { 28 Honee honee = *(Honee *)user_data; 29 Ceed ceed = honee->ceed; 30 PetscScalar dt; 31 Vec Q_loc = honee->Q_loc, R; 32 PetscMemType q_mem_type; 33 34 PetscFunctionBeginUser; 35 // Update time dependent data 36 PetscCall(UpdateBoundaryValues(honee, Q_loc, t)); 37 if (honee->phys->solution_time_label) 38 PetscCallCeed(ceed, CeedOperatorSetContextDouble(honee->op_rhs_ctx->op, honee->phys->solution_time_label, &t)); 39 PetscCall(TSGetTimeStep(ts, &dt)); 40 if (honee->phys->timestep_size_label) 41 PetscCallCeed(ceed, CeedOperatorSetContextDouble(honee->op_rhs_ctx->op, honee->phys->timestep_size_label, &dt)); 42 43 PetscCall(DMGetNamedGlobalVector(honee->dm, "RHS Residual", &R)); 44 PetscCall(DMGlobalToLocal(honee->dm, Q, INSERT_VALUES, Q_loc)); 45 if (honee->app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) PetscCall(DivDiffFluxProjectionApply(honee->diff_flux_proj, Q_loc)); 46 PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, R, honee->op_rhs_ctx)); 47 48 // Inverse of the mass matrix 49 PetscCall(VecReadPetscToCeed(Q_loc, &q_mem_type, honee->q_ceed)); 50 PetscCall(KSPSolve(honee->mass_ksp, R, G)); 51 PetscCall(VecReadCeedToPetsc(honee->q_ceed, q_mem_type, Q_loc)); 52 53 PetscCall(DMRestoreNamedGlobalVector(honee->dm, "RHS Residual", &R)); 54 PetscFunctionReturn(PETSC_SUCCESS); 55 } 56 57 // Surface forces function setup 58 static PetscErrorCode Surface_Forces_NS(DM dm, Vec G_loc, PetscInt num_walls, const PetscInt walls[], PetscScalar *reaction_force) { 59 DMLabel face_label; 60 const PetscScalar *g_array; 61 PetscInt dim = 3; 62 MPI_Comm comm = PetscObjectComm((PetscObject)dm); 63 PetscSection section; 64 65 PetscFunctionBeginUser; 66 PetscCall(DMGetLabel(dm, "Face Sets", &face_label)); 67 PetscCall(VecGetArrayRead(G_loc, &g_array)); 68 for (PetscInt w = 0; w < num_walls; w++) { 69 const PetscInt wall = walls[w], *points; 70 IS wall_is; 71 PetscInt num_points, num_comp = 0; 72 73 PetscCall(DMLabelGetStratumIS(face_label, wall, &wall_is)); 74 if (!wall_is) continue; // No wall points on this process, skip 75 76 PetscCall(DMGetLocalSection(dm, §ion)); 77 PetscCall(PetscSectionGetFieldComponents(section, 0, &num_comp)); 78 PetscCall(ISGetSize(wall_is, &num_points)); 79 PetscCall(ISGetIndices(wall_is, &points)); 80 for (PetscInt i = 0; i < num_points; i++) { 81 const PetscInt p = points[i]; 82 const StateConservative *r; 83 PetscInt dof; 84 85 PetscCall(DMPlexPointLocalRead(dm, p, g_array, &r)); 86 PetscCall(PetscSectionGetDof(section, p, &dof)); 87 for (PetscInt node = 0; node < dof / num_comp; node++) { 88 for (PetscInt j = 0; j < dim; j++) { 89 reaction_force[w * dim + j] -= r[node].momentum[j]; 90 } 91 } 92 } 93 PetscCall(ISRestoreIndices(wall_is, &points)); 94 PetscCall(ISDestroy(&wall_is)); 95 } 96 PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, reaction_force, dim * num_walls, MPIU_SCALAR, MPI_SUM, comm)); 97 PetscCall(VecRestoreArrayRead(G_loc, &g_array)); 98 PetscFunctionReturn(PETSC_SUCCESS); 99 } 100 101 // Implicit time-stepper function setup 102 PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) { 103 Honee honee = *(Honee *)user_data; 104 Ceed ceed = honee->ceed; 105 PetscScalar dt; 106 Vec Q_loc = honee->Q_loc, Q_dot_loc = honee->Q_dot_loc, G_loc; 107 PetscMemType q_mem_type, q_dot_mem_type, g_mem_type; 108 109 PetscFunctionBeginUser; 110 // Get local vectors 111 PetscCall(DMGetNamedLocalVector(honee->dm, "ResidualLocal", &G_loc)); 112 113 // Update time dependent data 114 PetscCall(UpdateBoundaryValues(honee, Q_loc, t)); 115 if (honee->phys->solution_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(honee->op_ifunction, honee->phys->solution_time_label, &t)); 116 PetscCall(TSGetTimeStep(ts, &dt)); 117 if (honee->phys->timestep_size_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(honee->op_ifunction, honee->phys->timestep_size_label, &dt)); 118 119 // Global-to-local 120 PetscCall(DMGlobalToLocalBegin(honee->dm, Q, INSERT_VALUES, Q_loc)); 121 PetscCall(DMGlobalToLocalBegin(honee->dm, Q_dot, INSERT_VALUES, Q_dot_loc)); 122 PetscCall(DMGlobalToLocalEnd(honee->dm, Q, INSERT_VALUES, Q_loc)); 123 if (honee->app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) PetscCall(DivDiffFluxProjectionApply(honee->diff_flux_proj, Q_loc)); 124 PetscCall(DMGlobalToLocalEnd(honee->dm, Q_dot, INSERT_VALUES, Q_dot_loc)); 125 126 // Place PETSc vectors in CEED vectors 127 PetscCall(VecReadPetscToCeed(Q_loc, &q_mem_type, honee->q_ceed)); 128 PetscCall(VecReadPetscToCeed(Q_dot_loc, &q_dot_mem_type, honee->q_dot_ceed)); 129 PetscCall(VecPetscToCeed(G_loc, &g_mem_type, honee->g_ceed)); 130 131 // Apply CEED operator 132 PetscCall(PetscLogEventBegin(HONEE_CeedOperatorApply, Q, G, 0, 0)); 133 PetscCall(PetscLogGpuTimeBegin()); 134 PetscCallCeed(honee->ceed, CeedOperatorApply(honee->op_ifunction, honee->q_ceed, honee->g_ceed, CEED_REQUEST_IMMEDIATE)); 135 PetscCall(PetscLogGpuTimeEnd()); 136 PetscCall(PetscLogEventEnd(HONEE_CeedOperatorApply, Q, G, 0, 0)); 137 138 // Restore vectors 139 PetscCall(VecReadCeedToPetsc(honee->q_ceed, q_mem_type, Q_loc)); 140 PetscCall(VecReadCeedToPetsc(honee->q_dot_ceed, q_dot_mem_type, Q_dot_loc)); 141 PetscCall(VecCeedToPetsc(honee->g_ceed, g_mem_type, G_loc)); 142 143 if (honee->app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) { 144 PetscCall(SgsDDApplyIFunction(honee, Q_loc, G_loc)); 145 } 146 147 // Local-to-Global 148 PetscCall(VecZeroEntries(G)); 149 PetscCall(DMLocalToGlobal(honee->dm, G_loc, ADD_VALUES, G)); 150 151 // Restore vectors 152 PetscCall(DMRestoreNamedLocalVector(honee->dm, "ResidualLocal", &G_loc)); 153 PetscFunctionReturn(PETSC_SUCCESS); 154 } 155 156 PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, PetscReal shift, Mat J, Mat J_pre, void *user_data) { 157 Honee honee = *(Honee *)user_data; 158 PetscBool J_is_matceed, J_is_mffd, J_pre_is_matceed, J_pre_is_mffd; 159 160 PetscFunctionBeginUser; 161 PetscCall(PetscObjectTypeCompare((PetscObject)J, MATMFFD, &J_is_mffd)); 162 PetscCall(PetscObjectTypeCompare((PetscObject)J, MATCEED, &J_is_matceed)); 163 PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATMFFD, &J_pre_is_mffd)); 164 PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATCEED, &J_pre_is_matceed)); 165 166 PetscCall(MatCeedSetContextReal(honee->mat_ijacobian, "ijacobian time shift", shift)); 167 168 if (J_is_matceed || J_is_mffd) { 169 PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 170 PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 171 } else PetscCall(MatCeedAssembleCOO(honee->mat_ijacobian, J)); 172 173 if (J_pre_is_matceed && J != J_pre) { 174 PetscCall(MatAssemblyBegin(J_pre, MAT_FINAL_ASSEMBLY)); 175 PetscCall(MatAssemblyEnd(J_pre, MAT_FINAL_ASSEMBLY)); 176 } else if (!J_pre_is_matceed && !J_pre_is_mffd && J != J_pre) { 177 PetscCall(MatCeedAssembleCOO(honee->mat_ijacobian, J_pre)); 178 } 179 PetscFunctionReturn(PETSC_SUCCESS); 180 } 181 182 PetscErrorCode WriteOutput(Honee honee, Vec Q, PetscInt step_no, PetscScalar time) { 183 Vec Q_loc; 184 char file_path[PETSC_MAX_PATH_LEN]; 185 PetscViewer viewer; 186 187 PetscFunctionBeginUser; 188 if (honee->app_ctx->checkpoint_vtk) { 189 // Set up output 190 PetscCall(DMGetLocalVector(honee->dm, &Q_loc)); 191 PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec")); 192 PetscCall(VecZeroEntries(Q_loc)); 193 PetscCall(DMGlobalToLocal(honee->dm, Q, INSERT_VALUES, Q_loc)); 194 195 // Output 196 PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-%03" PetscInt_FMT ".vtu", honee->app_ctx->output_dir, step_no)); 197 198 PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path, FILE_MODE_WRITE, &viewer)); 199 PetscCall(VecView(Q_loc, viewer)); 200 PetscCall(PetscViewerDestroy(&viewer)); 201 if (honee->dm_viz) { 202 Vec Q_refined, Q_refined_loc; 203 char file_path_refined[PETSC_MAX_PATH_LEN]; 204 PetscViewer viewer_refined; 205 206 PetscCall(DMGetGlobalVector(honee->dm_viz, &Q_refined)); 207 PetscCall(DMGetLocalVector(honee->dm_viz, &Q_refined_loc)); 208 PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined")); 209 210 PetscCall(MatInterpolate(honee->interp_viz, Q, Q_refined)); 211 PetscCall(VecZeroEntries(Q_refined_loc)); 212 PetscCall(DMGlobalToLocal(honee->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc)); 213 214 PetscCall( 215 PetscSNPrintf(file_path_refined, sizeof file_path_refined, "%s/nsrefined-%03" PetscInt_FMT ".vtu", honee->app_ctx->output_dir, step_no)); 216 217 PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined), file_path_refined, FILE_MODE_WRITE, &viewer_refined)); 218 PetscCall(VecView(Q_refined_loc, viewer_refined)); 219 PetscCall(DMRestoreLocalVector(honee->dm_viz, &Q_refined_loc)); 220 PetscCall(DMRestoreGlobalVector(honee->dm_viz, &Q_refined)); 221 PetscCall(PetscViewerDestroy(&viewer_refined)); 222 } 223 PetscCall(DMRestoreLocalVector(honee->dm, &Q_loc)); 224 } 225 226 // Save data in a binary file for continuation of simulations 227 if (honee->app_ctx->add_stepnum2bin) { 228 PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", honee->app_ctx->output_dir, step_no)); 229 } else { 230 PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", honee->app_ctx->output_dir)); 231 } 232 PetscCall(PetscViewerBinaryOpen(honee->comm, file_path, FILE_MODE_WRITE, &viewer)); 233 234 time /= honee->units->second; // Dimensionalize time back 235 PetscCall(HoneeWriteBinaryVec(viewer, Q, time, step_no)); 236 PetscCall(PetscViewerDestroy(&viewer)); 237 PetscFunctionReturn(PETSC_SUCCESS); 238 } 239 240 // CSV Monitor 241 PetscErrorCode TSMonitor_WallForce(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) { 242 Honee honee = ctx; 243 Vec G_loc; 244 PetscInt num_wall = honee->app_ctx->wall_forces.num_wall, dim = 3; 245 const PetscInt *walls = honee->app_ctx->wall_forces.walls; 246 PetscViewer viewer = honee->app_ctx->wall_forces.viewer; 247 PetscViewerFormat format = honee->app_ctx->wall_forces.viewer_format; 248 PetscScalar *reaction_force; 249 PetscBool is_ascii; 250 251 PetscFunctionBeginUser; 252 if (!viewer) PetscFunctionReturn(PETSC_SUCCESS); 253 PetscCall(DMGetNamedLocalVector(honee->dm, "ResidualLocal", &G_loc)); 254 PetscCall(PetscCalloc1(num_wall * dim, &reaction_force)); 255 PetscCall(Surface_Forces_NS(honee->dm, G_loc, num_wall, walls, reaction_force)); 256 PetscCall(DMRestoreNamedLocalVector(honee->dm, "ResidualLocal", &G_loc)); 257 258 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &is_ascii)); 259 260 if (is_ascii) { 261 if (format == PETSC_VIEWER_ASCII_CSV && !honee->app_ctx->wall_forces.header_written) { 262 PetscCall(PetscViewerASCIIPrintf(viewer, "Step,Time,Wall,ForceX,ForceY,ForceZ\n")); 263 honee->app_ctx->wall_forces.header_written = PETSC_TRUE; 264 } 265 for (PetscInt w = 0; w < num_wall; w++) { 266 PetscInt wall = walls[w]; 267 if (format == PETSC_VIEWER_ASCII_CSV) { 268 PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT ",%g,%" PetscInt_FMT ",%g,%g,%g\n", step_no, time, wall, 269 reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2])); 270 271 } else { 272 PetscCall(PetscViewerASCIIPrintf(viewer, "Wall %" PetscInt_FMT " Forces: Force_x = %12g, Force_y = %12g, Force_z = %12g\n", wall, 273 reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2])); 274 } 275 } 276 } 277 PetscCall(PetscFree(reaction_force)); 278 PetscFunctionReturn(PETSC_SUCCESS); 279 } 280 281 // User provided TS Monitor 282 PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) { 283 Honee honee = ctx; 284 285 PetscFunctionBeginUser; 286 // Print every 'checkpoint_interval' steps 287 if (honee->app_ctx->checkpoint_interval <= 0 || step_no % honee->app_ctx->checkpoint_interval != 0 || 288 (honee->app_ctx->cont_steps == step_no && step_no != 0)) { 289 PetscFunctionReturn(PETSC_SUCCESS); 290 } 291 292 PetscCall(WriteOutput(honee, Q, step_no, time)); 293 PetscFunctionReturn(PETSC_SUCCESS); 294 } 295 296 PetscErrorCode TSPostStep_CheckStep(TS ts) { 297 Honee honee; 298 PetscReal norm; 299 PetscInt step; 300 Vec Q; 301 302 PetscFunctionBeginUser; 303 PetscCall(TSGetApplicationContext(ts, &honee)); 304 PetscCall(TSGetStepNumber(ts, &step)); 305 PetscCall(TSGetSolution(ts, &Q)); 306 if (step % honee->app_ctx->check_step_interval) PetscFunctionReturn(PETSC_SUCCESS); 307 PetscCall(VecNorm(Q, NORM_1, &norm)); 308 if (PetscIsInfOrNanReal(norm)) { 309 PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "Solution diverged: Nans found in solution\n")); 310 PetscCall(TSSetConvergedReason(ts, TS_DIVERGED_NONLINEAR_SOLVE)); 311 } 312 PetscFunctionReturn(PETSC_SUCCESS); 313 } 314 315 PetscErrorCode TSPostStep_MaxWallTime(TS ts) { 316 Honee honee; 317 PetscInt step; 318 PetscMPIInt rank; 319 MPI_Comm comm; 320 PetscBool is_wall_time_exceeded = PETSC_FALSE; 321 322 PetscFunctionBeginUser; 323 PetscCall(TSGetApplicationContext(ts, &honee)); 324 PetscCall(TSGetStepNumber(ts, &step)); 325 if (step % honee->max_wall_time_interval) PetscFunctionReturn(PETSC_SUCCESS); 326 PetscCall(PetscObjectGetComm((PetscObject)ts, &comm)); 327 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 328 if (rank == 0) is_wall_time_exceeded = time(NULL) > honee->max_wall_time ? PETSC_TRUE : PETSC_FALSE; 329 // Must broadcast to avoid race condition 330 PetscCallMPI(MPI_Bcast(&is_wall_time_exceeded, 1, MPIU_BOOL, 0, comm)); 331 if (PetscUnlikely(is_wall_time_exceeded)) { 332 PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "Stopping TSSolve: Set max wall time exceeded\n")); 333 PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_USER)); 334 } 335 PetscFunctionReturn(PETSC_SUCCESS); 336 } 337 338 /** 339 @brief TSPostStep for HONEE 340 341 `TSSetPostStep()` only accepts a single function argument, so this function groups together all post-step 342 functionality needed for HONEE features 343 344 @param[in] ts TS object 345 **/ 346 PetscErrorCode TSPostStep_Honee(TS ts) { 347 Honee honee; 348 349 PetscFunctionBeginUser; 350 PetscCall(TSGetApplicationContext(ts, &honee)); 351 if (honee->max_wall_time > 0) PetscCall(TSPostStep_MaxWallTime(ts)); 352 if (honee->app_ctx->sgs_train_enable) PetscCall(TSPostStep_SGS_DD_Training(ts)); 353 if (honee->app_ctx->check_step_interval > 0) PetscCall(TSPostStep_CheckStep(ts)); 354 PetscFunctionReturn(PETSC_SUCCESS); 355 } 356 357 // TS: Create, setup, and solve 358 PetscErrorCode TSSolve_NS(DM dm, Honee honee, AppCtx app_ctx, Physics phys, ProblemData problem, Vec *Q, PetscScalar *f_time, TS *ts) { 359 MPI_Comm comm = honee->comm; 360 TSAdapt adapt; 361 PetscScalar final_time; 362 363 PetscFunctionBeginUser; 364 PetscCall(TSCreate(comm, ts)); 365 PetscCall(TSSetDM(*ts, dm)); 366 PetscCall(TSSetApplicationContext(*ts, honee)); 367 if (phys->implicit) { 368 PetscCall(TSSetType(*ts, TSBDF)); 369 if (honee->op_ifunction) { 370 PetscCall(TSSetIFunction(*ts, NULL, IFunction_NS, &honee)); 371 } else { // Implicit integrators can fall back to using an RHSFunction 372 PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &honee)); 373 } 374 if (honee->mat_ijacobian) { 375 PetscCall(DMTSSetIJacobian(dm, FormIJacobian_NS, &honee)); 376 } 377 } else { 378 PetscCheck(honee->op_rhs_ctx, comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction"); 379 PetscCall(TSSetType(*ts, TSRK)); 380 PetscCall(TSRKSetType(*ts, TSRK5F)); 381 PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &honee)); 382 } 383 PetscCall(TSSetMaxTime(*ts, 500. * honee->units->second)); 384 PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER)); 385 if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(TSSetErrorIfStepFails(*ts, PETSC_FALSE)); 386 PetscCall(TSSetTimeStep(*ts, 1.e-2 * honee->units->second)); 387 PetscCall(TSGetAdapt(*ts, &adapt)); 388 PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * honee->units->second, 1.e2 * honee->units->second)); 389 PetscCall(TSSetFromOptions(*ts)); 390 if (honee->mat_ijacobian) { 391 if (app_ctx->amat_type && !strcmp(app_ctx->amat_type, MATSHELL)) { 392 SNES snes; 393 KSP ksp; 394 Mat Pmat, Amat; 395 396 PetscCall(TSGetSNES(*ts, &snes)); 397 PetscCall(SNESGetKSP(snes, &ksp)); 398 PetscCall(CreateSolveOperatorsFromMatCeed(ksp, honee->mat_ijacobian, PETSC_FALSE, &Amat, &Pmat)); 399 PetscCall(TSSetIJacobian(*ts, honee->mat_ijacobian, Pmat, NULL, NULL)); 400 PetscCall(MatDestroy(&Amat)); 401 PetscCall(MatDestroy(&Pmat)); 402 } 403 } 404 honee->time_bc_set = -1.0; // require all BCs be updated 405 if (app_ctx->cont_steps) { // continue from previous timestep data 406 PetscCall(TSSetTime(*ts, app_ctx->cont_time * honee->units->second)); 407 PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps)); 408 } 409 if (honee->set_poststep) PetscCall(TSSetPostStep(*ts, TSPostStep_Honee)); 410 if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(TSMonitorSet(*ts, TSMonitor_NS, honee, NULL)); 411 if (app_ctx->wall_forces.viewer) PetscCall(TSMonitorSet(*ts, TSMonitor_WallForce, honee, NULL)); 412 if (app_ctx->turb_spanstats_enable) { 413 PetscCall(TSMonitorSet(*ts, TSMonitor_TurbulenceStatistics, honee, NULL)); 414 CeedScalar previous_time = app_ctx->cont_time * honee->units->second; 415 PetscCallCeed(honee->ceed, 416 CeedOperatorSetContextDouble(honee->spanstats.op_stats_collect_ctx->op, honee->spanstats.previous_time_label, &previous_time)); 417 } 418 PetscCall(TSMonitorSetFromOptions(*ts, "-ts_monitor_total_kinetic_energy", "Monitor total kinetic energy balance terms in the domain", NULL, 419 TSMonitor_TotalKineticEnergy, SetupMontiorTotalKineticEnergy)); 420 PetscCall(TSMonitorSetFromOptions(*ts, "-ts_monitor_cfl", "Monitor element CFL statistics", NULL, TSMonitor_Cfl, SetupMontiorCfl)); 421 if (app_ctx->diff_filter_monitor) PetscCall(TSMonitorSet(*ts, TSMonitor_DifferentialFilter, honee, NULL)); 422 if (app_ctx->sgs_train_enable) PetscCall(TSMonitorSet(*ts, TSMonitor_SGS_DD_Training, honee, NULL)); 423 424 if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(PrintRunInfo(honee, honee->phys, problem, *ts)); 425 // Solve 426 PetscReal start_time; 427 PetscInt start_step; 428 PetscCall(TSGetTime(*ts, &start_time)); 429 PetscCall(TSGetStepNumber(*ts, &start_step)); 430 431 PetscCall(PetscLogDefaultBegin()); // So we can use PetscLogStageGetPerfInfo without -log_view 432 PetscPreLoadBegin(PETSC_FALSE, "HONEE Solve"); 433 PetscCall(TSSetTime(*ts, start_time)); 434 PetscCall(TSSetStepNumber(*ts, start_step)); 435 if (PetscPreLoadingOn) { 436 // LCOV_EXCL_START 437 SNES snes; 438 KSP ksp; 439 Vec Q_preload; 440 PetscReal rtol_snes, rtol_ksp; 441 PetscCall(VecDuplicate(*Q, &Q_preload)); 442 PetscCall(VecCopy(*Q, Q_preload)); 443 PetscCall(TSGetSNES(*ts, &snes)); 444 PetscCall(SNESGetTolerances(snes, NULL, &rtol_snes, NULL, NULL, NULL)); 445 PetscCall(SNESSetTolerances(snes, PETSC_CURRENT, .99, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT)); 446 PetscCall(SNESGetKSP(snes, &ksp)); 447 PetscCall(KSPGetTolerances(ksp, &rtol_ksp, NULL, NULL, NULL)); 448 PetscCall(KSPSetTolerances(ksp, .99, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT)); 449 PetscCall(TSSetSolution(*ts, Q_preload)); 450 PetscCall(TSStep(*ts)); 451 PetscCall(SNESSetTolerances(snes, PETSC_CURRENT, rtol_snes, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT)); 452 PetscCall(KSPSetTolerances(ksp, rtol_ksp, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT)); 453 PetscCall(VecDestroy(&Q_preload)); 454 // LCOV_EXCL_STOP 455 } else { 456 PetscCall(PetscBarrier((PetscObject)*ts)); 457 PetscCall(TSSolve(*ts, *Q)); 458 } 459 PetscPreLoadEnd(); 460 461 PetscCall(TSGetSolveTime(*ts, &final_time)); 462 *f_time = final_time; 463 464 if (app_ctx->test_type == TESTTYPE_NONE) { 465 PetscInt step_no; 466 PetscCall(TSGetStepNumber(*ts, &step_no)); 467 if (honee->app_ctx->checkpoint_interval > 0 || honee->app_ctx->checkpoint_interval == -1) { 468 PetscCall(WriteOutput(honee, *Q, step_no, final_time)); 469 } 470 471 PetscLogStage stage_id; 472 PetscEventPerfInfo stage_perf; 473 474 PetscCall(PetscLogStageGetId("HONEE Solve", &stage_id)); 475 PetscCall(PetscLogStageGetPerfInfo(stage_id, &stage_perf)); 476 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time taken for solution (sec): %g\n", stage_perf.time)); 477 } 478 PetscFunctionReturn(PETSC_SUCCESS); 479 } 480