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