1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3 // 4 // SPDX-License-Identifier: BSD-2-Clause 5 // 6 // This file is part of CEED: http://github.com/ceed 7 8 /// @file 9 /// Time-stepping functions for Navier-Stokes example using PETSc 10 11 #include <ceed.h> 12 #include <petscdmplex.h> 13 #include <petscts.h> 14 15 #include "../navierstokes.h" 16 #include "../qfunctions/newtonian_state.h" 17 18 // Compute mass matrix for explicit scheme 19 PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm, CeedData ceed_data, Vec M) { 20 CeedQFunction qf_mass; 21 CeedOperator op_mass; 22 OperatorApplyContext op_mass_ctx; 23 Vec Ones_loc; 24 CeedInt num_comp_q, q_data_size; 25 PetscFunctionBeginUser; 26 27 // CEED Restriction 28 CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q); 29 CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &q_data_size); 30 31 // CEED QFunction 32 PetscCall(CreateMassQFunction(ceed, num_comp_q, q_data_size, &qf_mass)); 33 34 // CEED Operator 35 CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass); 36 CeedOperatorSetField(op_mass, "u", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 37 CeedOperatorSetField(op_mass, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data); 38 CeedOperatorSetField(op_mass, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 39 40 PetscCall(OperatorApplyContextCreate(NULL, dm, ceed, op_mass, NULL, NULL, NULL, NULL, &op_mass_ctx)); 41 42 PetscCall(DMGetLocalVector(dm, &Ones_loc)); 43 PetscCall(VecSet(Ones_loc, 1)); 44 PetscCall(ApplyCeedOperatorLocalToGlobal(Ones_loc, M, op_mass_ctx)); 45 46 // Invert diagonally lumped mass vector for RHS function 47 PetscCall(VecReciprocal(M)); 48 49 // Cleanup 50 PetscCall(OperatorApplyContextDestroy(op_mass_ctx)); 51 PetscCall(DMRestoreLocalVector(dm, &Ones_loc)); 52 CeedQFunctionDestroy(&qf_mass); 53 CeedOperatorDestroy(&op_mass); 54 55 PetscFunctionReturn(0); 56 } 57 58 // Insert Boundary values if it's a new time 59 PetscErrorCode UpdateBoundaryValues(User user, Vec Q_loc, PetscReal t) { 60 PetscFunctionBeginUser; 61 if (user->time_bc_set != t) { 62 PetscCall(DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL)); 63 user->time_bc_set = t; 64 } 65 PetscFunctionReturn(0); 66 } 67 68 // @brief Update the context label value to new value if necessary. 69 // @note This only supports labels with scalar label values (ie. not arrays) 70 PetscErrorCode UpdateContextLabel(MPI_Comm comm, PetscScalar update_value, CeedOperator op, CeedContextFieldLabel label) { 71 PetscScalar label_value; 72 73 PetscFunctionBeginUser; 74 PetscCheck(label, comm, PETSC_ERR_ARG_BADPTR, "Label should be non-NULL"); 75 76 { 77 size_t num_elements; 78 const PetscScalar *label_values; 79 CeedOperatorGetContextDoubleRead(op, label, &num_elements, &label_values); 80 PetscCheck(num_elements == 1, comm, PETSC_ERR_SUP, "%s does not support labels with more than 1 value. Label has %zu values", __func__, 81 num_elements); 82 label_value = *label_values; 83 CeedOperatorRestoreContextDoubleRead(op, label, &label_values); 84 } 85 86 if (label_value != update_value) { 87 CeedOperatorSetContextDouble(op, label, &update_value); 88 } 89 PetscFunctionReturn(0); 90 } 91 92 // RHS (Explicit time-stepper) function setup 93 // This is the RHS of the ODE, given as u_t = G(t,u) 94 // This function takes in a state vector Q and writes into G 95 PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) { 96 User user = *(User *)user_data; 97 MPI_Comm comm = PetscObjectComm((PetscObject)ts); 98 PetscScalar dt; 99 Vec Q_loc = user->Q_loc; 100 PetscFunctionBeginUser; 101 102 // Update time dependent data 103 PetscCall(UpdateBoundaryValues(user, Q_loc, t)); 104 if (user->phys->solution_time_label) PetscCall(UpdateContextLabel(comm, t, user->op_rhs_ctx->op, user->phys->solution_time_label)); 105 PetscCall(TSGetTimeStep(ts, &dt)); 106 if (user->phys->timestep_size_label) PetscCall(UpdateContextLabel(comm, dt, user->op_rhs_ctx->op, user->phys->timestep_size_label)); 107 108 PetscCall(ApplyCeedOperatorGlobalToGlobal(Q, G, user->op_rhs_ctx)); 109 110 // Inverse of the lumped mass matrix (M is Minv) 111 PetscCall(VecPointwiseMult(G, G, user->M_inv)); 112 113 PetscFunctionReturn(0); 114 } 115 116 // Surface forces function setup 117 static PetscErrorCode Surface_Forces_NS(DM dm, Vec G_loc, PetscInt num_walls, const PetscInt walls[], PetscScalar *reaction_force) { 118 DMLabel face_label; 119 const PetscScalar *g; 120 PetscInt dof, dim = 3; 121 MPI_Comm comm; 122 PetscSection s; 123 124 PetscFunctionBeginUser; 125 PetscCall(PetscArrayzero(reaction_force, num_walls * dim)); 126 PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 127 PetscCall(DMGetLabel(dm, "Face Sets", &face_label)); 128 PetscCall(VecGetArrayRead(G_loc, &g)); 129 for (PetscInt w = 0; w < num_walls; w++) { 130 const PetscInt wall = walls[w]; 131 IS wall_is; 132 PetscCall(DMGetLocalSection(dm, &s)); 133 PetscCall(DMLabelGetStratumIS(face_label, wall, &wall_is)); 134 if (wall_is) { // There exist such points on this process 135 PetscInt num_points; 136 PetscInt num_comp = 0; 137 const PetscInt *points; 138 PetscCall(PetscSectionGetFieldComponents(s, 0, &num_comp)); 139 PetscCall(ISGetSize(wall_is, &num_points)); 140 PetscCall(ISGetIndices(wall_is, &points)); 141 for (PetscInt i = 0; i < num_points; i++) { 142 const PetscInt p = points[i]; 143 const StateConservative *r; 144 PetscCall(DMPlexPointLocalRead(dm, p, g, &r)); 145 PetscCall(PetscSectionGetDof(s, p, &dof)); 146 for (PetscInt node = 0; node < dof / num_comp; node++) { 147 for (PetscInt j = 0; j < 3; j++) { 148 reaction_force[w * dim + j] -= r[node].momentum[j]; 149 } 150 } 151 } 152 PetscCall(ISRestoreIndices(wall_is, &points)); 153 } 154 PetscCall(ISDestroy(&wall_is)); 155 } 156 PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, reaction_force, dim * num_walls, MPIU_SCALAR, MPI_SUM, comm)); 157 // Restore Vectors 158 PetscCall(VecRestoreArrayRead(G_loc, &g)); 159 160 PetscFunctionReturn(0); 161 } 162 163 // Implicit time-stepper function setup 164 PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) { 165 User user = *(User *)user_data; 166 MPI_Comm comm = PetscObjectComm((PetscObject)ts); 167 PetscScalar dt; 168 Vec Q_loc = user->Q_loc, Q_dot_loc = user->Q_dot_loc, G_loc; 169 PetscMemType q_mem_type, q_dot_mem_type, g_mem_type; 170 PetscFunctionBeginUser; 171 172 // Get local vectors 173 PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); 174 175 // Update time dependent data 176 PetscCall(UpdateBoundaryValues(user, Q_loc, t)); 177 if (user->phys->solution_time_label) PetscCall(UpdateContextLabel(comm, t, user->op_ifunction, user->phys->solution_time_label)); 178 PetscCall(TSGetTimeStep(ts, &dt)); 179 if (user->phys->timestep_size_label) PetscCall(UpdateContextLabel(comm, dt, user->op_ifunction, user->phys->timestep_size_label)); 180 181 // Global-to-local 182 PetscCall(DMGlobalToLocalBegin(user->dm, Q, INSERT_VALUES, Q_loc)); 183 PetscCall(DMGlobalToLocalBegin(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc)); 184 PetscCall(DMGlobalToLocalEnd(user->dm, Q, INSERT_VALUES, Q_loc)); 185 PetscCall(DMGlobalToLocalEnd(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc)); 186 187 // Place PETSc vectors in CEED vectors 188 PetscCall(VecReadP2C(Q_loc, &q_mem_type, user->q_ceed)); 189 PetscCall(VecReadP2C(Q_dot_loc, &q_dot_mem_type, user->q_dot_ceed)); 190 PetscCall(VecP2C(G_loc, &g_mem_type, user->g_ceed)); 191 192 // Apply CEED operator 193 CeedOperatorApply(user->op_ifunction, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE); 194 195 // Restore vectors 196 PetscCall(VecReadC2P(user->q_ceed, q_mem_type, Q_loc)); 197 PetscCall(VecReadC2P(user->q_dot_ceed, q_dot_mem_type, Q_dot_loc)); 198 PetscCall(VecC2P(user->g_ceed, g_mem_type, G_loc)); 199 200 if (user->app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) { 201 PetscCall(SGS_DD_ModelApplyIFunction(user, Q_loc, G_loc)); 202 } 203 204 // Local-to-Global 205 PetscCall(VecZeroEntries(G)); 206 PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G)); 207 208 // Restore vectors 209 PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); 210 211 PetscFunctionReturn(0); 212 } 213 214 static PetscErrorCode FormPreallocation(User user, PetscBool pbdiagonal, Mat J, CeedVector *coo_values) { 215 PetscCount ncoo; 216 PetscInt *rows, *cols; 217 218 PetscFunctionBeginUser; 219 if (pbdiagonal) { 220 CeedSize l_size; 221 CeedOperatorGetActiveVectorLengths(user->op_ijacobian, &l_size, NULL); 222 ncoo = l_size * 5; 223 rows = malloc(ncoo * sizeof(rows[0])); 224 cols = malloc(ncoo * sizeof(cols[0])); 225 for (PetscCount n = 0; n < l_size / 5; n++) { 226 for (PetscInt i = 0; i < 5; i++) { 227 for (PetscInt j = 0; j < 5; j++) { 228 rows[(n * 5 + i) * 5 + j] = n * 5 + i; 229 cols[(n * 5 + i) * 5 + j] = n * 5 + j; 230 } 231 } 232 } 233 } else { 234 PetscCall(CeedOperatorLinearAssembleSymbolic(user->op_ijacobian, &ncoo, &rows, &cols)); 235 } 236 PetscCall(MatSetPreallocationCOOLocal(J, ncoo, rows, cols)); 237 free(rows); 238 free(cols); 239 CeedVectorCreate(user->ceed, ncoo, coo_values); 240 PetscFunctionReturn(0); 241 } 242 243 static PetscErrorCode FormSetValues(User user, PetscBool pbdiagonal, Mat J, CeedVector coo_values) { 244 CeedMemType mem_type = CEED_MEM_HOST; 245 const PetscScalar *values; 246 MatType mat_type; 247 248 PetscFunctionBeginUser; 249 PetscCall(MatGetType(J, &mat_type)); 250 if (strstr(mat_type, "kokkos") || strstr(mat_type, "cusparse")) mem_type = CEED_MEM_DEVICE; 251 if (pbdiagonal) { 252 CeedOperatorLinearAssemblePointBlockDiagonal(user->op_ijacobian, coo_values, CEED_REQUEST_IMMEDIATE); 253 } else { 254 CeedOperatorLinearAssemble(user->op_ijacobian, coo_values); 255 } 256 CeedVectorGetArrayRead(coo_values, mem_type, &values); 257 PetscCall(MatSetValuesCOO(J, values, INSERT_VALUES)); 258 CeedVectorRestoreArrayRead(coo_values, &values); 259 PetscFunctionReturn(0); 260 } 261 262 PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, PetscReal shift, Mat J, Mat J_pre, void *user_data) { 263 User user = *(User *)user_data; 264 PetscBool J_is_shell, J_is_mffd, J_pre_is_shell; 265 PetscFunctionBeginUser; 266 if (user->phys->ijacobian_time_shift_label) CeedOperatorSetContextDouble(user->op_ijacobian, user->phys->ijacobian_time_shift_label, &shift); 267 PetscCall(PetscObjectTypeCompare((PetscObject)J, MATMFFD, &J_is_mffd)); 268 PetscCall(PetscObjectTypeCompare((PetscObject)J, MATSHELL, &J_is_shell)); 269 PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATSHELL, &J_pre_is_shell)); 270 if (!user->matrices_set_up) { 271 if (J_is_shell) { 272 OperatorApplyContext op_ijacobian_ctx; 273 OperatorApplyContextCreate(user->dm, user->dm, user->ceed, user->op_ijacobian, user->q_ceed, user->g_ceed, user->Q_dot_loc, NULL, 274 &op_ijacobian_ctx); 275 PetscCall(MatShellSetContext(J, op_ijacobian_ctx)); 276 PetscCall(MatShellSetContextDestroy(J, (PetscErrorCode(*)(void *))OperatorApplyContextDestroy)); 277 PetscCall(MatShellSetOperation(J, MATOP_MULT, (void (*)(void))MatMult_Ceed)); 278 PetscCall(MatShellSetOperation(J, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiag_Ceed)); 279 PetscCall(MatSetUp(J)); 280 } 281 if (!J_pre_is_shell) { 282 PetscCall(FormPreallocation(user, user->app_ctx->pmat_pbdiagonal, J_pre, &user->coo_values_pmat)); 283 } 284 if (J != J_pre && !J_is_shell && !J_is_mffd) { 285 PetscCall(FormPreallocation(user, PETSC_FALSE, J, &user->coo_values_amat)); 286 } 287 user->matrices_set_up = true; 288 } 289 if (!J_pre_is_shell) { 290 PetscCall(FormSetValues(user, user->app_ctx->pmat_pbdiagonal, J_pre, user->coo_values_pmat)); 291 } 292 if (user->coo_values_amat) { 293 PetscCall(FormSetValues(user, PETSC_FALSE, J, user->coo_values_amat)); 294 } else if (J_is_mffd) { 295 PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 296 PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 297 } 298 PetscFunctionReturn(0); 299 } 300 301 PetscErrorCode WriteOutput(User user, Vec Q, PetscInt step_no, PetscScalar time) { 302 Vec Q_loc; 303 char file_path[PETSC_MAX_PATH_LEN]; 304 PetscViewer viewer; 305 PetscFunctionBeginUser; 306 307 if (user->app_ctx->checkpoint_vtk) { 308 // Set up output 309 PetscCall(DMGetLocalVector(user->dm, &Q_loc)); 310 PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec")); 311 PetscCall(VecZeroEntries(Q_loc)); 312 PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc)); 313 314 // Output 315 PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no)); 316 317 PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path, FILE_MODE_WRITE, &viewer)); 318 PetscCall(VecView(Q_loc, viewer)); 319 PetscCall(PetscViewerDestroy(&viewer)); 320 if (user->dm_viz) { 321 Vec Q_refined, Q_refined_loc; 322 char file_path_refined[PETSC_MAX_PATH_LEN]; 323 PetscViewer viewer_refined; 324 325 PetscCall(DMGetGlobalVector(user->dm_viz, &Q_refined)); 326 PetscCall(DMGetLocalVector(user->dm_viz, &Q_refined_loc)); 327 PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined")); 328 329 PetscCall(MatInterpolate(user->interp_viz, Q, Q_refined)); 330 PetscCall(VecZeroEntries(Q_refined_loc)); 331 PetscCall(DMGlobalToLocal(user->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc)); 332 333 PetscCall( 334 PetscSNPrintf(file_path_refined, sizeof file_path_refined, "%s/nsrefined-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no)); 335 336 PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined), file_path_refined, FILE_MODE_WRITE, &viewer_refined)); 337 PetscCall(VecView(Q_refined_loc, viewer_refined)); 338 PetscCall(DMRestoreLocalVector(user->dm_viz, &Q_refined_loc)); 339 PetscCall(DMRestoreGlobalVector(user->dm_viz, &Q_refined)); 340 PetscCall(PetscViewerDestroy(&viewer_refined)); 341 } 342 PetscCall(DMRestoreLocalVector(user->dm, &Q_loc)); 343 } 344 345 // Save data in a binary file for continuation of simulations 346 if (user->app_ctx->add_stepnum2bin) { 347 PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", user->app_ctx->output_dir, step_no)); 348 } else { 349 PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", user->app_ctx->output_dir)); 350 } 351 PetscCall(PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, &viewer)); 352 353 PetscInt token = FLUIDS_FILE_TOKEN; 354 PetscCall(PetscViewerBinaryWrite(viewer, &token, 1, PETSC_INT)); 355 PetscCall(PetscViewerBinaryWrite(viewer, &step_no, 1, PETSC_INT)); 356 time /= user->units->second; // Dimensionalize time back 357 PetscCall(PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL)); 358 PetscCall(VecView(Q, viewer)); 359 PetscCall(PetscViewerDestroy(&viewer)); 360 PetscFunctionReturn(0); 361 } 362 363 // CSV Monitor 364 PetscErrorCode TSMonitor_WallForce(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) { 365 User user = ctx; 366 Vec G_loc; 367 PetscInt num_wall = user->app_ctx->wall_forces.num_wall, dim = 3; 368 const PetscInt *walls = user->app_ctx->wall_forces.walls; 369 PetscViewer viewer = user->app_ctx->wall_forces.viewer; 370 PetscViewerFormat format = user->app_ctx->wall_forces.viewer_format; 371 PetscScalar *reaction_force; 372 PetscBool iascii; 373 374 PetscFunctionBeginUser; 375 if (!viewer) PetscFunctionReturn(0); 376 PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); 377 PetscCall(PetscMalloc1(num_wall * dim, &reaction_force)); 378 PetscCall(Surface_Forces_NS(user->dm, G_loc, num_wall, walls, reaction_force)); 379 PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); 380 381 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 382 383 if (iascii) { 384 if (format == PETSC_VIEWER_ASCII_CSV && !user->app_ctx->wall_forces.header_written) { 385 PetscCall(PetscViewerASCIIPrintf(viewer, "Step,Time,Wall,ForceX,ForceY,ForceZ\n")); 386 user->app_ctx->wall_forces.header_written = PETSC_TRUE; 387 } 388 for (PetscInt w = 0; w < num_wall; w++) { 389 PetscInt wall = walls[w]; 390 if (format == PETSC_VIEWER_ASCII_CSV) { 391 PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT ",%g,%" PetscInt_FMT ",%g,%g,%g\n", step_no, time, wall, 392 reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2])); 393 394 } else { 395 PetscCall(PetscViewerASCIIPrintf(viewer, "Wall %" PetscInt_FMT " Forces: Force_x = %12g, Force_y = %12g, Force_z = %12g\n", wall, 396 reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2])); 397 } 398 } 399 } 400 PetscCall(PetscFree(reaction_force)); 401 PetscFunctionReturn(0); 402 } 403 404 // User provided TS Monitor 405 PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) { 406 User user = ctx; 407 PetscFunctionBeginUser; 408 409 // Print every 'checkpoint_interval' steps 410 if (user->app_ctx->checkpoint_interval <= 0 || step_no % user->app_ctx->checkpoint_interval != 0 || 411 (user->app_ctx->cont_steps == step_no && step_no != 0)) { 412 PetscFunctionReturn(0); 413 } 414 415 PetscCall(WriteOutput(user, Q, step_no, time)); 416 417 PetscFunctionReturn(0); 418 } 419 420 // TS: Create, setup, and solve 421 PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, Vec *Q, PetscScalar *f_time, TS *ts) { 422 MPI_Comm comm = user->comm; 423 TSAdapt adapt; 424 PetscScalar final_time; 425 PetscFunctionBeginUser; 426 427 PetscCall(TSCreate(comm, ts)); 428 PetscCall(TSSetDM(*ts, dm)); 429 if (phys->implicit) { 430 PetscCall(TSSetType(*ts, TSBDF)); 431 if (user->op_ifunction) { 432 PetscCall(TSSetIFunction(*ts, NULL, IFunction_NS, &user)); 433 } else { // Implicit integrators can fall back to using an RHSFunction 434 PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user)); 435 } 436 if (user->op_ijacobian) { 437 PetscCall(DMTSSetIJacobian(dm, FormIJacobian_NS, &user)); 438 if (app_ctx->amat_type) { 439 Mat Pmat, Amat; 440 PetscCall(DMCreateMatrix(dm, &Pmat)); 441 PetscCall(DMSetMatType(dm, app_ctx->amat_type)); 442 PetscCall(DMCreateMatrix(dm, &Amat)); 443 PetscCall(TSSetIJacobian(*ts, Amat, Pmat, NULL, NULL)); 444 PetscCall(MatDestroy(&Amat)); 445 PetscCall(MatDestroy(&Pmat)); 446 } 447 } 448 } else { 449 PetscCheck(user->op_rhs_ctx, comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction"); 450 PetscCall(TSSetType(*ts, TSRK)); 451 PetscCall(TSRKSetType(*ts, TSRK5F)); 452 PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user)); 453 } 454 PetscCall(TSSetMaxTime(*ts, 500. * user->units->second)); 455 PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER)); 456 if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(TSSetErrorIfStepFails(*ts, PETSC_FALSE)); 457 PetscCall(TSSetTimeStep(*ts, 1.e-2 * user->units->second)); 458 if (app_ctx->test_type != TESTTYPE_NONE) { 459 PetscCall(TSSetMaxSteps(*ts, 10)); 460 } 461 PetscCall(TSGetAdapt(*ts, &adapt)); 462 PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * user->units->second, 1.e2 * user->units->second)); 463 PetscCall(TSSetFromOptions(*ts)); 464 user->time_bc_set = -1.0; // require all BCs be updated 465 if (!app_ctx->cont_steps) { // print initial condition 466 if (app_ctx->test_type == TESTTYPE_NONE) { 467 PetscCall(TSMonitor_NS(*ts, 0, 0., *Q, user)); 468 } 469 } else { // continue from time of last output 470 PetscInt count; 471 PetscViewer viewer; 472 473 if (app_ctx->cont_time <= 0) { // Legacy files did not include step number and time 474 PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_time_file, FILE_MODE_READ, &viewer)); 475 PetscCall(PetscViewerBinaryRead(viewer, &app_ctx->cont_time, 1, &count, PETSC_REAL)); 476 PetscCall(PetscViewerDestroy(&viewer)); 477 PetscCheck(app_ctx->cont_steps != -1, comm, PETSC_ERR_ARG_INCOMP, 478 "-continue step number not specified, but checkpoint file does not contain a step number (likely written by older code version)"); 479 } 480 PetscCall(TSSetTime(*ts, app_ctx->cont_time * user->units->second)); 481 PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps)); 482 } 483 if (app_ctx->test_type == TESTTYPE_NONE) { 484 PetscCall(TSMonitorSet(*ts, TSMonitor_NS, user, NULL)); 485 } 486 if (app_ctx->wall_forces.viewer) { 487 PetscCall(TSMonitorSet(*ts, TSMonitor_WallForce, user, NULL)); 488 } 489 if (app_ctx->turb_spanstats_enable) { 490 PetscCall(TSMonitorSet(*ts, TSMonitor_TurbulenceStatistics, user, NULL)); 491 CeedScalar previous_time = app_ctx->cont_time * user->units->second; 492 CeedOperatorSetContextDouble(user->spanstats.op_stats_collect_ctx->op, user->spanstats.previous_time_label, &previous_time); 493 } 494 if (app_ctx->diff_filter_monitor) PetscCall(TSMonitorSet(*ts, TSMonitor_DifferentialFilter, user, NULL)); 495 496 // Solve 497 PetscReal start_time; 498 PetscInt start_step; 499 PetscCall(TSGetTime(*ts, &start_time)); 500 PetscCall(TSGetStepNumber(*ts, &start_step)); 501 502 PetscPreLoadBegin(PETSC_FALSE, "Fluids Solve"); 503 PetscCall(TSSetTime(*ts, start_time)); 504 PetscCall(TSSetStepNumber(*ts, start_step)); 505 if (PetscPreLoadingOn) { 506 // LCOV_EXCL_START 507 SNES snes; 508 Vec Q_preload; 509 PetscReal rtol; 510 PetscCall(VecDuplicate(*Q, &Q_preload)); 511 PetscCall(VecCopy(*Q, Q_preload)); 512 PetscCall(TSGetSNES(*ts, &snes)); 513 PetscCall(SNESGetTolerances(snes, NULL, &rtol, NULL, NULL, NULL)); 514 PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, .99, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 515 PetscCall(TSSetSolution(*ts, *Q)); 516 PetscCall(TSStep(*ts)); 517 PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, rtol, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 518 PetscCall(VecDestroy(&Q_preload)); 519 // LCOV_EXCL_STOP 520 } else { 521 PetscCall(PetscBarrier((PetscObject)*ts)); 522 PetscCall(TSSolve(*ts, *Q)); 523 } 524 PetscPreLoadEnd(); 525 526 PetscCall(TSGetSolveTime(*ts, &final_time)); 527 *f_time = final_time; 528 529 if (app_ctx->test_type == TESTTYPE_NONE) { 530 PetscInt step_no; 531 PetscCall(TSGetStepNumber(*ts, &step_no)); 532 if (user->app_ctx->checkpoint_interval > 0 || user->app_ctx->checkpoint_interval == -1) { 533 PetscCall(WriteOutput(user, *Q, step_no, final_time)); 534 } 535 536 PetscLogEvent stage_id; 537 PetscStageLog stage_log; 538 539 PetscCall(PetscLogStageGetId("Fluids Solve", &stage_id)); 540 PetscCall(PetscLogGetStageLog(&stage_log)); 541 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time taken for solution (sec): %g\n", stage_log->stageInfo[stage_id].perfInfo.time)); 542 } 543 PetscFunctionReturn(0); 544 } 545