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