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 Navier-Stokes example using PETSc 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(User user, Vec Q_loc, PetscReal t) { 16 PetscFunctionBeginUser; 17 if (user->time_bc_set != t) { 18 PetscCall(DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL)); 19 user->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 User user = *(User *)user_data; 29 Ceed ceed = user->ceed; 30 PetscScalar dt; 31 Vec Q_loc = user->Q_loc; 32 PetscMemType q_mem_type; 33 34 PetscFunctionBeginUser; 35 // Update time dependent data 36 PetscCall(UpdateBoundaryValues(user, Q_loc, t)); 37 if (user->phys->solution_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_rhs_ctx->op, user->phys->solution_time_label, &t)); 38 PetscCall(TSGetTimeStep(ts, &dt)); 39 if (user->phys->timestep_size_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_rhs_ctx->op, user->phys->timestep_size_label, &dt)); 40 41 PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc)); 42 if (user->app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) PetscCall(DivDiffFluxProjectionApply(user->diff_flux_proj, Q_loc)); 43 PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, G, user->op_rhs_ctx)); 44 45 PetscCall(VecReadPetscToCeed(Q_loc, &q_mem_type, user->q_ceed)); 46 47 // Inverse of the mass matrix 48 PetscCall(KSPSolve(user->mass_ksp, G, G)); 49 50 PetscCall(VecReadCeedToPetsc(user->q_ceed, q_mem_type, Q_loc)); 51 PetscFunctionReturn(PETSC_SUCCESS); 52 } 53 54 // Surface forces function setup 55 static PetscErrorCode Surface_Forces_NS(DM dm, Vec G_loc, PetscInt num_walls, const PetscInt walls[], PetscScalar *reaction_force) { 56 DMLabel face_label; 57 const PetscScalar *g_array; 58 PetscInt dim = 3; 59 MPI_Comm comm = PetscObjectComm((PetscObject)dm); 60 PetscSection section; 61 62 PetscFunctionBeginUser; 63 PetscCall(DMGetLabel(dm, "Face Sets", &face_label)); 64 PetscCall(VecGetArrayRead(G_loc, &g_array)); 65 for (PetscInt w = 0; w < num_walls; w++) { 66 const PetscInt wall = walls[w], *points; 67 IS wall_is; 68 PetscInt num_points, num_comp = 0; 69 70 PetscCall(DMLabelGetStratumIS(face_label, wall, &wall_is)); 71 if (!wall_is) continue; // No wall points on this process, skip 72 73 PetscCall(DMGetLocalSection(dm, §ion)); 74 PetscCall(PetscSectionGetFieldComponents(section, 0, &num_comp)); 75 PetscCall(ISGetSize(wall_is, &num_points)); 76 PetscCall(ISGetIndices(wall_is, &points)); 77 for (PetscInt i = 0; i < num_points; i++) { 78 const PetscInt p = points[i]; 79 const StateConservative *r; 80 PetscInt dof; 81 82 PetscCall(DMPlexPointLocalRead(dm, p, g_array, &r)); 83 PetscCall(PetscSectionGetDof(section, p, &dof)); 84 for (PetscInt node = 0; node < dof / num_comp; node++) { 85 for (PetscInt j = 0; j < dim; j++) { 86 reaction_force[w * dim + j] -= r[node].momentum[j]; 87 } 88 } 89 } 90 PetscCall(ISRestoreIndices(wall_is, &points)); 91 PetscCall(ISDestroy(&wall_is)); 92 } 93 PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, reaction_force, dim * num_walls, MPIU_SCALAR, MPI_SUM, comm)); 94 PetscCall(VecRestoreArrayRead(G_loc, &g_array)); 95 PetscFunctionReturn(PETSC_SUCCESS); 96 } 97 98 // Implicit time-stepper function setup 99 PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) { 100 User user = *(User *)user_data; 101 Ceed ceed = user->ceed; 102 PetscScalar dt; 103 Vec Q_loc = user->Q_loc, Q_dot_loc = user->Q_dot_loc, G_loc; 104 PetscMemType q_mem_type, q_dot_mem_type, g_mem_type; 105 106 PetscFunctionBeginUser; 107 // Get local vectors 108 PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); 109 110 // Update time dependent data 111 PetscCall(UpdateBoundaryValues(user, Q_loc, t)); 112 if (user->phys->solution_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_ifunction, user->phys->solution_time_label, &t)); 113 PetscCall(TSGetTimeStep(ts, &dt)); 114 if (user->phys->timestep_size_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_ifunction, user->phys->timestep_size_label, &dt)); 115 116 // Global-to-local 117 PetscCall(DMGlobalToLocalBegin(user->dm, Q, INSERT_VALUES, Q_loc)); 118 PetscCall(DMGlobalToLocalBegin(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc)); 119 PetscCall(DMGlobalToLocalEnd(user->dm, Q, INSERT_VALUES, Q_loc)); 120 if (user->app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) PetscCall(DivDiffFluxProjectionApply(user->diff_flux_proj, Q_loc)); 121 PetscCall(DMGlobalToLocalEnd(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc)); 122 123 // Place PETSc vectors in CEED vectors 124 PetscCall(VecReadPetscToCeed(Q_loc, &q_mem_type, user->q_ceed)); 125 PetscCall(VecReadPetscToCeed(Q_dot_loc, &q_dot_mem_type, user->q_dot_ceed)); 126 PetscCall(VecPetscToCeed(G_loc, &g_mem_type, user->g_ceed)); 127 128 // Apply CEED operator 129 PetscCall(PetscLogEventBegin(FLUIDS_CeedOperatorApply, Q, G, 0, 0)); 130 PetscCall(PetscLogGpuTimeBegin()); 131 PetscCallCeed(user->ceed, CeedOperatorApply(user->op_ifunction, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE)); 132 PetscCall(PetscLogGpuTimeEnd()); 133 PetscCall(PetscLogEventEnd(FLUIDS_CeedOperatorApply, Q, G, 0, 0)); 134 135 // Restore vectors 136 PetscCall(VecReadCeedToPetsc(user->q_ceed, q_mem_type, Q_loc)); 137 PetscCall(VecReadCeedToPetsc(user->q_dot_ceed, q_dot_mem_type, Q_dot_loc)); 138 PetscCall(VecCeedToPetsc(user->g_ceed, g_mem_type, G_loc)); 139 140 if (user->app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) { 141 PetscCall(SgsDDApplyIFunction(user, Q_loc, G_loc)); 142 } 143 144 // Local-to-Global 145 PetscCall(VecZeroEntries(G)); 146 PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G)); 147 148 // Restore vectors 149 PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); 150 PetscFunctionReturn(PETSC_SUCCESS); 151 } 152 153 PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, PetscReal shift, Mat J, Mat J_pre, void *user_data) { 154 User user = *(User *)user_data; 155 PetscBool J_is_matceed, J_is_mffd, J_pre_is_matceed, J_pre_is_mffd; 156 157 PetscFunctionBeginUser; 158 PetscCall(PetscObjectTypeCompare((PetscObject)J, MATMFFD, &J_is_mffd)); 159 PetscCall(PetscObjectTypeCompare((PetscObject)J, MATCEED, &J_is_matceed)); 160 PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATMFFD, &J_pre_is_mffd)); 161 PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATCEED, &J_pre_is_matceed)); 162 163 PetscCall(MatCeedSetContextReal(user->mat_ijacobian, "ijacobian time shift", shift)); 164 165 if (J_is_matceed || J_is_mffd) { 166 PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 167 PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 168 } else PetscCall(MatCeedAssembleCOO(user->mat_ijacobian, J)); 169 170 if (J_pre_is_matceed && J != J_pre) { 171 PetscCall(MatAssemblyBegin(J_pre, MAT_FINAL_ASSEMBLY)); 172 PetscCall(MatAssemblyEnd(J_pre, MAT_FINAL_ASSEMBLY)); 173 } else if (!J_pre_is_matceed && !J_pre_is_mffd && J != J_pre) { 174 PetscCall(MatCeedAssembleCOO(user->mat_ijacobian, J_pre)); 175 } 176 PetscFunctionReturn(PETSC_SUCCESS); 177 } 178 179 PetscErrorCode WriteOutput(User user, Vec Q, PetscInt step_no, PetscScalar time) { 180 Vec Q_loc; 181 char file_path[PETSC_MAX_PATH_LEN]; 182 PetscViewer viewer; 183 184 PetscFunctionBeginUser; 185 if (user->app_ctx->checkpoint_vtk) { 186 // Set up output 187 PetscCall(DMGetLocalVector(user->dm, &Q_loc)); 188 PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec")); 189 PetscCall(VecZeroEntries(Q_loc)); 190 PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc)); 191 192 // Output 193 PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no)); 194 195 PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path, FILE_MODE_WRITE, &viewer)); 196 PetscCall(VecView(Q_loc, viewer)); 197 PetscCall(PetscViewerDestroy(&viewer)); 198 if (user->dm_viz) { 199 Vec Q_refined, Q_refined_loc; 200 char file_path_refined[PETSC_MAX_PATH_LEN]; 201 PetscViewer viewer_refined; 202 203 PetscCall(DMGetGlobalVector(user->dm_viz, &Q_refined)); 204 PetscCall(DMGetLocalVector(user->dm_viz, &Q_refined_loc)); 205 PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined")); 206 207 PetscCall(MatInterpolate(user->interp_viz, Q, Q_refined)); 208 PetscCall(VecZeroEntries(Q_refined_loc)); 209 PetscCall(DMGlobalToLocal(user->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc)); 210 211 PetscCall( 212 PetscSNPrintf(file_path_refined, sizeof file_path_refined, "%s/nsrefined-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no)); 213 214 PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined), file_path_refined, FILE_MODE_WRITE, &viewer_refined)); 215 PetscCall(VecView(Q_refined_loc, viewer_refined)); 216 PetscCall(DMRestoreLocalVector(user->dm_viz, &Q_refined_loc)); 217 PetscCall(DMRestoreGlobalVector(user->dm_viz, &Q_refined)); 218 PetscCall(PetscViewerDestroy(&viewer_refined)); 219 } 220 PetscCall(DMRestoreLocalVector(user->dm, &Q_loc)); 221 } 222 223 // Save data in a binary file for continuation of simulations 224 if (user->app_ctx->add_stepnum2bin) { 225 PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", user->app_ctx->output_dir, step_no)); 226 } else { 227 PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", user->app_ctx->output_dir)); 228 } 229 PetscCall(PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, &viewer)); 230 231 time /= user->units->second; // Dimensionalize time back 232 PetscCall(HoneeWriteBinaryVec(viewer, Q, time, step_no)); 233 PetscCall(PetscViewerDestroy(&viewer)); 234 PetscFunctionReturn(PETSC_SUCCESS); 235 } 236 237 // CSV Monitor 238 PetscErrorCode TSMonitor_WallForce(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) { 239 User user = ctx; 240 Vec G_loc; 241 PetscInt num_wall = user->app_ctx->wall_forces.num_wall, dim = 3; 242 const PetscInt *walls = user->app_ctx->wall_forces.walls; 243 PetscViewer viewer = user->app_ctx->wall_forces.viewer; 244 PetscViewerFormat format = user->app_ctx->wall_forces.viewer_format; 245 PetscScalar *reaction_force; 246 PetscBool is_ascii; 247 248 PetscFunctionBeginUser; 249 if (!viewer) PetscFunctionReturn(PETSC_SUCCESS); 250 PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); 251 PetscCall(PetscCalloc1(num_wall * dim, &reaction_force)); 252 PetscCall(Surface_Forces_NS(user->dm, G_loc, num_wall, walls, reaction_force)); 253 PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); 254 255 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &is_ascii)); 256 257 if (is_ascii) { 258 if (format == PETSC_VIEWER_ASCII_CSV && !user->app_ctx->wall_forces.header_written) { 259 PetscCall(PetscViewerASCIIPrintf(viewer, "Step,Time,Wall,ForceX,ForceY,ForceZ\n")); 260 user->app_ctx->wall_forces.header_written = PETSC_TRUE; 261 } 262 for (PetscInt w = 0; w < num_wall; w++) { 263 PetscInt wall = walls[w]; 264 if (format == PETSC_VIEWER_ASCII_CSV) { 265 PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT ",%g,%" PetscInt_FMT ",%g,%g,%g\n", step_no, time, wall, 266 reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2])); 267 268 } else { 269 PetscCall(PetscViewerASCIIPrintf(viewer, "Wall %" PetscInt_FMT " Forces: Force_x = %12g, Force_y = %12g, Force_z = %12g\n", wall, 270 reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2])); 271 } 272 } 273 } 274 PetscCall(PetscFree(reaction_force)); 275 PetscFunctionReturn(PETSC_SUCCESS); 276 } 277 278 // User provided TS Monitor 279 PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) { 280 User user = ctx; 281 282 PetscFunctionBeginUser; 283 // Print every 'checkpoint_interval' steps 284 if (user->app_ctx->checkpoint_interval <= 0 || step_no % user->app_ctx->checkpoint_interval != 0 || 285 (user->app_ctx->cont_steps == step_no && step_no != 0)) { 286 PetscFunctionReturn(PETSC_SUCCESS); 287 } 288 289 PetscCall(WriteOutput(user, Q, step_no, time)); 290 PetscFunctionReturn(PETSC_SUCCESS); 291 } 292 293 // TS: Create, setup, and solve 294 PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, ProblemData problem, Vec *Q, PetscScalar *f_time, TS *ts) { 295 MPI_Comm comm = user->comm; 296 TSAdapt adapt; 297 PetscScalar final_time; 298 299 PetscFunctionBeginUser; 300 PetscCall(TSCreate(comm, ts)); 301 PetscCall(TSSetDM(*ts, dm)); 302 PetscCall(TSSetApplicationContext(*ts, user)); 303 if (phys->implicit) { 304 PetscCall(TSSetType(*ts, TSBDF)); 305 if (user->op_ifunction) { 306 PetscCall(TSSetIFunction(*ts, NULL, IFunction_NS, &user)); 307 } else { // Implicit integrators can fall back to using an RHSFunction 308 PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user)); 309 } 310 if (user->mat_ijacobian) { 311 PetscCall(DMTSSetIJacobian(dm, FormIJacobian_NS, &user)); 312 } 313 } else { 314 PetscCheck(user->op_rhs_ctx, comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction"); 315 PetscCall(TSSetType(*ts, TSRK)); 316 PetscCall(TSRKSetType(*ts, TSRK5F)); 317 PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user)); 318 } 319 PetscCall(TSSetMaxTime(*ts, 500. * user->units->second)); 320 PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER)); 321 if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(TSSetErrorIfStepFails(*ts, PETSC_FALSE)); 322 PetscCall(TSSetTimeStep(*ts, 1.e-2 * user->units->second)); 323 PetscCall(TSGetAdapt(*ts, &adapt)); 324 PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * user->units->second, 1.e2 * user->units->second)); 325 PetscCall(TSSetFromOptions(*ts)); 326 if (user->mat_ijacobian) { 327 if (app_ctx->amat_type && !strcmp(app_ctx->amat_type, MATSHELL)) { 328 SNES snes; 329 KSP ksp; 330 Mat Pmat, Amat; 331 332 PetscCall(TSGetSNES(*ts, &snes)); 333 PetscCall(SNESGetKSP(snes, &ksp)); 334 PetscCall(CreateSolveOperatorsFromMatCeed(ksp, user->mat_ijacobian, PETSC_FALSE, &Amat, &Pmat)); 335 PetscCall(TSSetIJacobian(*ts, user->mat_ijacobian, Pmat, NULL, NULL)); 336 PetscCall(MatDestroy(&Amat)); 337 PetscCall(MatDestroy(&Pmat)); 338 } 339 } 340 user->time_bc_set = -1.0; // require all BCs be updated 341 if (app_ctx->cont_steps) { // continue from previous timestep data 342 PetscInt count; 343 PetscViewer viewer; 344 345 if (app_ctx->cont_time <= 0) { // Legacy files did not include step number and time 346 PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_time_file, FILE_MODE_READ, &viewer)); 347 PetscCall(PetscViewerBinaryRead(viewer, &app_ctx->cont_time, 1, &count, PETSC_REAL)); 348 PetscCall(PetscViewerDestroy(&viewer)); 349 PetscCheck(app_ctx->cont_steps != -1, comm, PETSC_ERR_ARG_INCOMP, 350 "-continue step number not specified, but checkpoint file does not contain a step number (likely written by older code version)"); 351 } 352 PetscCall(TSSetTime(*ts, app_ctx->cont_time * user->units->second)); 353 PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps)); 354 } 355 if (app_ctx->test_type == TESTTYPE_NONE) { 356 PetscCall(TSMonitorSet(*ts, TSMonitor_NS, user, NULL)); 357 } 358 if (app_ctx->wall_forces.viewer) { 359 PetscCall(TSMonitorSet(*ts, TSMonitor_WallForce, user, NULL)); 360 } 361 if (app_ctx->turb_spanstats_enable) { 362 PetscCall(TSMonitorSet(*ts, TSMonitor_TurbulenceStatistics, user, NULL)); 363 CeedScalar previous_time = app_ctx->cont_time * user->units->second; 364 PetscCallCeed(user->ceed, 365 CeedOperatorSetContextDouble(user->spanstats.op_stats_collect_ctx->op, user->spanstats.previous_time_label, &previous_time)); 366 } 367 if (app_ctx->diff_filter_monitor) PetscCall(TSMonitorSet(*ts, TSMonitor_DifferentialFilter, user, NULL)); 368 369 if (app_ctx->sgs_train_enable) { 370 PetscCall(TSMonitorSet(*ts, TSMonitor_SGS_DD_Training, user, NULL)); 371 PetscCall(TSSetPostStep(*ts, TSPostStep_SGS_DD_Training)); 372 } 373 374 if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(PrintRunInfo(user, user->phys, problem, *ts)); 375 // Solve 376 PetscReal start_time; 377 PetscInt start_step; 378 PetscCall(TSGetTime(*ts, &start_time)); 379 PetscCall(TSGetStepNumber(*ts, &start_step)); 380 381 PetscCall(PetscLogDefaultBegin()); // So we can use PetscLogStageGetPerfInfo without -log_view 382 PetscPreLoadBegin(PETSC_FALSE, "Fluids Solve"); 383 PetscCall(TSSetTime(*ts, start_time)); 384 PetscCall(TSSetStepNumber(*ts, start_step)); 385 if (PetscPreLoadingOn) { 386 // LCOV_EXCL_START 387 SNES snes; 388 Vec Q_preload; 389 PetscReal rtol; 390 PetscCall(VecDuplicate(*Q, &Q_preload)); 391 PetscCall(VecCopy(*Q, Q_preload)); 392 PetscCall(TSGetSNES(*ts, &snes)); 393 PetscCall(SNESGetTolerances(snes, NULL, &rtol, NULL, NULL, NULL)); 394 PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, .99, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 395 PetscCall(TSSetSolution(*ts, Q_preload)); 396 PetscCall(TSStep(*ts)); 397 PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, rtol, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 398 PetscCall(VecDestroy(&Q_preload)); 399 // LCOV_EXCL_STOP 400 } else { 401 PetscCall(PetscBarrier((PetscObject)*ts)); 402 PetscCall(TSSolve(*ts, *Q)); 403 } 404 PetscPreLoadEnd(); 405 406 PetscCall(TSGetSolveTime(*ts, &final_time)); 407 *f_time = final_time; 408 409 if (app_ctx->test_type == TESTTYPE_NONE) { 410 PetscInt step_no; 411 PetscCall(TSGetStepNumber(*ts, &step_no)); 412 if (user->app_ctx->checkpoint_interval > 0 || user->app_ctx->checkpoint_interval == -1) { 413 PetscCall(WriteOutput(user, *Q, step_no, final_time)); 414 } 415 416 PetscLogStage stage_id; 417 PetscEventPerfInfo stage_perf; 418 419 PetscCall(PetscLogStageGetId("Fluids Solve", &stage_id)); 420 PetscCall(PetscLogStageGetPerfInfo(stage_id, &stage_perf)); 421 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time taken for solution (sec): %g\n", stage_perf.time)); 422 } 423 PetscFunctionReturn(PETSC_SUCCESS); 424 } 425