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