13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 377841947SLeila Ghaffari // 43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 577841947SLeila Ghaffari // 63d8e8822SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 777841947SLeila Ghaffari 877841947SLeila Ghaffari /// @file 977841947SLeila Ghaffari /// Time-stepping functions for Navier-Stokes example using PETSc 1077841947SLeila Ghaffari 1149aac155SJeremy L Thompson #include <ceed.h> 1249aac155SJeremy L Thompson #include <petscdmplex.h> 1349aac155SJeremy L Thompson #include <petscts.h> 1449aac155SJeremy L Thompson 1577841947SLeila Ghaffari #include "../navierstokes.h" 16ca69d878SAdeleke O. Bankole #include "../qfunctions/newtonian_state.h" 1777841947SLeila Ghaffari 1877841947SLeila Ghaffari // Compute mass matrix for explicit scheme 192b730f8bSJeremy L Thompson PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm, CeedData ceed_data, Vec M) { 2077841947SLeila Ghaffari CeedQFunction qf_mass; 2177841947SLeila Ghaffari CeedOperator op_mass; 223221f4d3SJames Wright OperatorApplyContext op_mass_ctx; 233221f4d3SJames Wright Vec Ones_loc; 2477841947SLeila Ghaffari CeedInt num_comp_q, q_data_size; 2577841947SLeila Ghaffari 26f17d818dSJames Wright PetscFunctionBeginUser; 27a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q)); 28a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &q_data_size)); 2977841947SLeila Ghaffari 30ef080ff9SJames Wright PetscCall(CreateMassQFunction(ceed, num_comp_q, q_data_size, &qf_mass)); 31a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass)); 32a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 33*356036faSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data)); 34a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 3577841947SLeila Ghaffari 363221f4d3SJames Wright PetscCall(OperatorApplyContextCreate(NULL, dm, ceed, op_mass, NULL, NULL, NULL, NULL, &op_mass_ctx)); 3777841947SLeila Ghaffari 383221f4d3SJames Wright PetscCall(DMGetLocalVector(dm, &Ones_loc)); 393221f4d3SJames Wright PetscCall(VecSet(Ones_loc, 1)); 403221f4d3SJames Wright PetscCall(ApplyCeedOperatorLocalToGlobal(Ones_loc, M, op_mass_ctx)); 4177841947SLeila Ghaffari 4277841947SLeila Ghaffari // Invert diagonally lumped mass vector for RHS function 432b730f8bSJeremy L Thompson PetscCall(VecReciprocal(M)); 4477841947SLeila Ghaffari 4577841947SLeila Ghaffari // Cleanup 463221f4d3SJames Wright PetscCall(OperatorApplyContextDestroy(op_mass_ctx)); 473221f4d3SJames Wright PetscCall(DMRestoreLocalVector(dm, &Ones_loc)); 48a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 49a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass)); 50ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 5177841947SLeila Ghaffari } 5277841947SLeila Ghaffari 53a0b9a424SJames Wright // Insert Boundary values if it's a new time 54a0b9a424SJames Wright PetscErrorCode UpdateBoundaryValues(User user, Vec Q_loc, PetscReal t) { 55a0b9a424SJames Wright PetscFunctionBeginUser; 56a0b9a424SJames Wright if (user->time_bc_set != t) { 57a0b9a424SJames Wright PetscCall(DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL)); 58a0b9a424SJames Wright user->time_bc_set = t; 59a0b9a424SJames Wright } 60ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 61a0b9a424SJames Wright } 62a0b9a424SJames Wright 6315637395SJames Wright // @brief Update the context label value to new value if necessary. 6415637395SJames Wright // @note This only supports labels with scalar label values (ie. not arrays) 65a424bcd0SJames Wright PetscErrorCode UpdateContextLabel(Ceed ceed, MPI_Comm comm, PetscScalar update_value, CeedOperator op, CeedContextFieldLabel label) { 6615637395SJames Wright PetscScalar label_value; 67e42c1df6SJames Wright 6815637395SJames Wright PetscFunctionBeginUser; 6915637395SJames Wright PetscCheck(label, comm, PETSC_ERR_ARG_BADPTR, "Label should be non-NULL"); 7015637395SJames Wright 7115637395SJames Wright { 7215637395SJames Wright size_t num_elements; 7315637395SJames Wright const PetscScalar *label_values; 74a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorGetContextDoubleRead(op, label, &num_elements, &label_values)); 7515637395SJames Wright PetscCheck(num_elements == 1, comm, PETSC_ERR_SUP, "%s does not support labels with more than 1 value. Label has %zu values", __func__, 7615637395SJames Wright num_elements); 7715637395SJames Wright label_value = *label_values; 78a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorRestoreContextDoubleRead(op, label, &label_values)); 79e42c1df6SJames Wright } 8015637395SJames Wright 8115637395SJames Wright if (label_value != update_value) { 82a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorSetContextDouble(op, label, &update_value)); 83e42c1df6SJames Wright } 84ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 85e42c1df6SJames Wright } 86e42c1df6SJames Wright 8777841947SLeila Ghaffari // RHS (Explicit time-stepper) function setup 8877841947SLeila Ghaffari // This is the RHS of the ODE, given as u_t = G(t,u) 8977841947SLeila Ghaffari // This function takes in a state vector Q and writes into G 9077841947SLeila Ghaffari PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) { 9177841947SLeila Ghaffari User user = *(User *)user_data; 9215637395SJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)ts); 93c798d55aSJames Wright PetscScalar dt; 949ad5e8e4SJames Wright Vec Q_loc = user->Q_loc; 9577841947SLeila Ghaffari 96f17d818dSJames Wright PetscFunctionBeginUser; 975e82a6e1SJeremy L Thompson // Update time dependent data 98a0b9a424SJames Wright PetscCall(UpdateBoundaryValues(user, Q_loc, t)); 99a424bcd0SJames Wright if (user->phys->solution_time_label) PetscCall(UpdateContextLabel(user->ceed, comm, t, user->op_rhs_ctx->op, user->phys->solution_time_label)); 1002b730f8bSJeremy L Thompson PetscCall(TSGetTimeStep(ts, &dt)); 101a424bcd0SJames Wright if (user->phys->timestep_size_label) PetscCall(UpdateContextLabel(user->ceed, comm, dt, user->op_rhs_ctx->op, user->phys->timestep_size_label)); 10277841947SLeila Ghaffari 1039ad5e8e4SJames Wright PetscCall(ApplyCeedOperatorGlobalToGlobal(Q, G, user->op_rhs_ctx)); 10477841947SLeila Ghaffari 105d6e67e47SJames Wright // Inverse of the lumped mass matrix 106186595e6SJames Wright PetscCall(VecPointwiseMult(G, G, user->M_inv)); 107ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 10877841947SLeila Ghaffari } 10977841947SLeila Ghaffari 110ca69d878SAdeleke O. Bankole // Surface forces function setup 111ca69d878SAdeleke O. Bankole static PetscErrorCode Surface_Forces_NS(DM dm, Vec G_loc, PetscInt num_walls, const PetscInt walls[], PetscScalar *reaction_force) { 112ca69d878SAdeleke O. Bankole DMLabel face_label; 113ca69d878SAdeleke O. Bankole const PetscScalar *g; 114d6734f85SAdeleke O. Bankole PetscInt dof, dim = 3; 115ca69d878SAdeleke O. Bankole MPI_Comm comm; 116d6734f85SAdeleke O. Bankole PetscSection s; 117ca69d878SAdeleke O. Bankole 118ca69d878SAdeleke O. Bankole PetscFunctionBeginUser; 119ca69d878SAdeleke O. Bankole PetscCall(PetscArrayzero(reaction_force, num_walls * dim)); 120ca69d878SAdeleke O. Bankole PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 121ca69d878SAdeleke O. Bankole PetscCall(DMGetLabel(dm, "Face Sets", &face_label)); 122ca69d878SAdeleke O. Bankole PetscCall(VecGetArrayRead(G_loc, &g)); 123ca69d878SAdeleke O. Bankole for (PetscInt w = 0; w < num_walls; w++) { 124ca69d878SAdeleke O. Bankole const PetscInt wall = walls[w]; 125ca69d878SAdeleke O. Bankole IS wall_is; 126d6734f85SAdeleke O. Bankole PetscCall(DMGetLocalSection(dm, &s)); 127ca69d878SAdeleke O. Bankole PetscCall(DMLabelGetStratumIS(face_label, wall, &wall_is)); 128ca69d878SAdeleke O. Bankole if (wall_is) { // There exist such points on this process 129ca69d878SAdeleke O. Bankole PetscInt num_points; 130d6734f85SAdeleke O. Bankole PetscInt num_comp = 0; 131ca69d878SAdeleke O. Bankole const PetscInt *points; 132d6734f85SAdeleke O. Bankole PetscCall(PetscSectionGetFieldComponents(s, 0, &num_comp)); 133ca69d878SAdeleke O. Bankole PetscCall(ISGetSize(wall_is, &num_points)); 134ca69d878SAdeleke O. Bankole PetscCall(ISGetIndices(wall_is, &points)); 135ca69d878SAdeleke O. Bankole for (PetscInt i = 0; i < num_points; i++) { 136ca69d878SAdeleke O. Bankole const PetscInt p = points[i]; 137ca69d878SAdeleke O. Bankole const StateConservative *r; 138ca69d878SAdeleke O. Bankole PetscCall(DMPlexPointLocalRead(dm, p, g, &r)); 139d6734f85SAdeleke O. Bankole PetscCall(PetscSectionGetDof(s, p, &dof)); 140d6734f85SAdeleke O. Bankole for (PetscInt node = 0; node < dof / num_comp; node++) { 141ca69d878SAdeleke O. Bankole for (PetscInt j = 0; j < 3; j++) { 142d6734f85SAdeleke O. Bankole reaction_force[w * dim + j] -= r[node].momentum[j]; 143d6734f85SAdeleke O. Bankole } 144ca69d878SAdeleke O. Bankole } 145ca69d878SAdeleke O. Bankole } 146ca69d878SAdeleke O. Bankole PetscCall(ISRestoreIndices(wall_is, &points)); 147ca69d878SAdeleke O. Bankole } 148ca69d878SAdeleke O. Bankole PetscCall(ISDestroy(&wall_is)); 149ca69d878SAdeleke O. Bankole } 150ca69d878SAdeleke O. Bankole PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, reaction_force, dim * num_walls, MPIU_SCALAR, MPI_SUM, comm)); 151ca69d878SAdeleke O. Bankole // Restore Vectors 152ca69d878SAdeleke O. Bankole PetscCall(VecRestoreArrayRead(G_loc, &g)); 153ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 154ca69d878SAdeleke O. Bankole } 155ca69d878SAdeleke O. Bankole 15677841947SLeila Ghaffari // Implicit time-stepper function setup 1572b730f8bSJeremy L Thompson PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) { 15877841947SLeila Ghaffari User user = *(User *)user_data; 15915637395SJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)ts); 160c798d55aSJames Wright PetscScalar dt; 1615e82a6e1SJeremy L Thompson Vec Q_loc = user->Q_loc, Q_dot_loc = user->Q_dot_loc, G_loc; 16277841947SLeila Ghaffari PetscMemType q_mem_type, q_dot_mem_type, g_mem_type; 16377841947SLeila Ghaffari 164f17d818dSJames Wright PetscFunctionBeginUser; 1655e82a6e1SJeremy L Thompson // Get local vectors 166ca69d878SAdeleke O. Bankole PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); 1675e82a6e1SJeremy L Thompson 1685e82a6e1SJeremy L Thompson // Update time dependent data 169a0b9a424SJames Wright PetscCall(UpdateBoundaryValues(user, Q_loc, t)); 170a424bcd0SJames Wright if (user->phys->solution_time_label) PetscCall(UpdateContextLabel(user->ceed, comm, t, user->op_ifunction, user->phys->solution_time_label)); 1712b730f8bSJeremy L Thompson PetscCall(TSGetTimeStep(ts, &dt)); 172a424bcd0SJames Wright if (user->phys->timestep_size_label) PetscCall(UpdateContextLabel(user->ceed, comm, dt, user->op_ifunction, user->phys->timestep_size_label)); 17377841947SLeila Ghaffari 17477841947SLeila Ghaffari // Global-to-local 175878eb0e7SJames Wright PetscCall(DMGlobalToLocalBegin(user->dm, Q, INSERT_VALUES, Q_loc)); 176878eb0e7SJames Wright PetscCall(DMGlobalToLocalBegin(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc)); 177878eb0e7SJames Wright PetscCall(DMGlobalToLocalEnd(user->dm, Q, INSERT_VALUES, Q_loc)); 178878eb0e7SJames Wright PetscCall(DMGlobalToLocalEnd(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc)); 17977841947SLeila Ghaffari 18077841947SLeila Ghaffari // Place PETSc vectors in CEED vectors 181c798d55aSJames Wright PetscCall(VecReadP2C(Q_loc, &q_mem_type, user->q_ceed)); 182c798d55aSJames Wright PetscCall(VecReadP2C(Q_dot_loc, &q_dot_mem_type, user->q_dot_ceed)); 183c798d55aSJames Wright PetscCall(VecP2C(G_loc, &g_mem_type, user->g_ceed)); 18477841947SLeila Ghaffari 18577841947SLeila Ghaffari // Apply CEED operator 18675d1798cSJames Wright PetscCall(PetscLogEventBegin(FLUIDS_CeedOperatorApply, Q, G, 0, 0)); 18775d1798cSJames Wright PetscCall(PetscLogGpuTimeBegin()); 188a424bcd0SJames Wright PetscCallCeed(user->ceed, CeedOperatorApply(user->op_ifunction, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE)); 18975d1798cSJames Wright PetscCall(PetscLogGpuTimeEnd()); 19075d1798cSJames Wright PetscCall(PetscLogEventEnd(FLUIDS_CeedOperatorApply, Q, G, 0, 0)); 19177841947SLeila Ghaffari 19277841947SLeila Ghaffari // Restore vectors 193c798d55aSJames Wright PetscCall(VecReadC2P(user->q_ceed, q_mem_type, Q_loc)); 194c798d55aSJames Wright PetscCall(VecReadC2P(user->q_dot_ceed, q_dot_mem_type, Q_dot_loc)); 195c798d55aSJames Wright PetscCall(VecC2P(user->g_ceed, g_mem_type, G_loc)); 19677841947SLeila Ghaffari 19719ffbc25SJames Wright if (user->app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) { 198cbef7084SJames Wright PetscCall(SgsDDModelApplyIFunction(user, Q_loc, G_loc)); 19919ffbc25SJames Wright } 200be34b3b7SJames Wright 20177841947SLeila Ghaffari // Local-to-Global 2022b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(G)); 2032b730f8bSJeremy L Thompson PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G)); 20477841947SLeila Ghaffari 20577841947SLeila Ghaffari // Restore vectors 206ca69d878SAdeleke O. Bankole PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); 207ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 20877841947SLeila Ghaffari } 20977841947SLeila Ghaffari 2102b730f8bSJeremy L Thompson static PetscErrorCode FormPreallocation(User user, PetscBool pbdiagonal, Mat J, CeedVector *coo_values) { 211544be873SJed Brown PetscCount ncoo; 212457e73b2SJames Wright PetscInt *rows_petsc, *cols_petsc; 213544be873SJed Brown 214544be873SJed Brown PetscFunctionBeginUser; 215544be873SJed Brown if (pbdiagonal) { 216544be873SJed Brown CeedSize l_size; 217a424bcd0SJames Wright PetscCallCeed(user->ceed, CeedOperatorGetActiveVectorLengths(user->op_ijacobian, &l_size, NULL)); 218544be873SJed Brown ncoo = l_size * 5; 219457e73b2SJames Wright rows_petsc = malloc(ncoo * sizeof(rows_petsc)); 220457e73b2SJames Wright cols_petsc = malloc(ncoo * sizeof(cols_petsc)); 221544be873SJed Brown for (PetscCount n = 0; n < l_size / 5; n++) { 222544be873SJed Brown for (PetscInt i = 0; i < 5; i++) { 223544be873SJed Brown for (PetscInt j = 0; j < 5; j++) { 224457e73b2SJames Wright rows_petsc[(n * 5 + i) * 5 + j] = n * 5 + i; 225457e73b2SJames Wright cols_petsc[(n * 5 + i) * 5 + j] = n * 5 + j; 226544be873SJed Brown } 227544be873SJed Brown } 228544be873SJed Brown } 229544be873SJed Brown } else { 230457e73b2SJames Wright CeedInt *rows_ceed, *cols_ceed; 231a424bcd0SJames Wright PetscCallCeed(user->ceed, CeedOperatorLinearAssembleSymbolic(user->op_ijacobian, &ncoo, &rows_ceed, &cols_ceed)); 232457e73b2SJames Wright PetscCall(IntArrayC2P(ncoo, &rows_ceed, &rows_petsc)); 233457e73b2SJames Wright PetscCall(IntArrayC2P(ncoo, &cols_ceed, &cols_petsc)); 234544be873SJed Brown } 235457e73b2SJames Wright PetscCall(MatSetPreallocationCOOLocal(J, ncoo, rows_petsc, cols_petsc)); 236457e73b2SJames Wright free(rows_petsc); 237457e73b2SJames Wright free(cols_petsc); 238a424bcd0SJames Wright PetscCallCeed(user->ceed, CeedVectorCreate(user->ceed, ncoo, coo_values)); 239ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 240544be873SJed Brown } 241544be873SJed Brown 2422b730f8bSJeremy L Thompson static PetscErrorCode FormSetValues(User user, PetscBool pbdiagonal, Mat J, CeedVector coo_values) { 243544be873SJed Brown CeedMemType mem_type = CEED_MEM_HOST; 244544be873SJed Brown const PetscScalar *values; 245544be873SJed Brown MatType mat_type; 246544be873SJed Brown 247544be873SJed Brown PetscFunctionBeginUser; 248544be873SJed Brown PetscCall(MatGetType(J, &mat_type)); 2492b730f8bSJeremy L Thompson if (strstr(mat_type, "kokkos") || strstr(mat_type, "cusparse")) mem_type = CEED_MEM_DEVICE; 250fbb546ddSJames Wright if (pbdiagonal) { 25175d1798cSJames Wright PetscCall(PetscLogEventBegin(FLUIDS_CeedOperatorAssemblePointBlockDiagonal, J, 0, 0, 0)); 25275d1798cSJames Wright PetscCall(PetscLogGpuTimeBegin()); 253a424bcd0SJames Wright PetscCallCeed(user->ceed, CeedOperatorLinearAssemblePointBlockDiagonal(user->op_ijacobian, coo_values, CEED_REQUEST_IMMEDIATE)); 25475d1798cSJames Wright PetscCall(PetscLogGpuTimeEnd()); 25575d1798cSJames Wright PetscCall(PetscLogEventEnd(FLUIDS_CeedOperatorAssemblePointBlockDiagonal, J, 0, 0, 0)); 256544be873SJed Brown } else { 25775d1798cSJames Wright PetscCall(PetscLogEventBegin(FLUIDS_CeedOperatorAssemble, J, 0, 0, 0)); 25875d1798cSJames Wright PetscCall(PetscLogGpuTimeBegin()); 259a424bcd0SJames Wright PetscCallCeed(user->ceed, CeedOperatorLinearAssemble(user->op_ijacobian, coo_values)); 26075d1798cSJames Wright PetscCall(PetscLogGpuTimeEnd()); 26175d1798cSJames Wright PetscCall(PetscLogEventEnd(FLUIDS_CeedOperatorAssemble, J, 0, 0, 0)); 262544be873SJed Brown } 263a424bcd0SJames Wright PetscCallCeed(user->ceed, CeedVectorGetArrayRead(coo_values, mem_type, &values)); 264544be873SJed Brown PetscCall(MatSetValuesCOO(J, values, INSERT_VALUES)); 265a424bcd0SJames Wright PetscCallCeed(user->ceed, CeedVectorRestoreArrayRead(coo_values, &values)); 266ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 267544be873SJed Brown } 268544be873SJed Brown 2692b730f8bSJeremy L Thompson PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, PetscReal shift, Mat J, Mat J_pre, void *user_data) { 270e334ad8fSJed Brown User user = *(User *)user_data; 271a424bcd0SJames Wright Ceed ceed = user->ceed; 272d6bf345cSJed Brown PetscBool J_is_shell, J_is_mffd, J_pre_is_shell; 273f17d818dSJames Wright 274e334ad8fSJed Brown PetscFunctionBeginUser; 275a424bcd0SJames Wright if (user->phys->ijacobian_time_shift_label) 276a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_ijacobian, user->phys->ijacobian_time_shift_label, &shift)); 277d6bf345cSJed Brown PetscCall(PetscObjectTypeCompare((PetscObject)J, MATMFFD, &J_is_mffd)); 278e334ad8fSJed Brown PetscCall(PetscObjectTypeCompare((PetscObject)J, MATSHELL, &J_is_shell)); 2792b730f8bSJeremy L Thompson PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATSHELL, &J_pre_is_shell)); 280e334ad8fSJed Brown if (!user->matrices_set_up) { 281e334ad8fSJed Brown if (J_is_shell) { 282b07979f9SJames Wright OperatorApplyContext op_ijacobian_ctx; 283b07979f9SJames Wright OperatorApplyContextCreate(user->dm, user->dm, user->ceed, user->op_ijacobian, user->q_ceed, user->g_ceed, user->Q_dot_loc, NULL, 284b07979f9SJames Wright &op_ijacobian_ctx); 285b07979f9SJames Wright PetscCall(MatShellSetContext(J, op_ijacobian_ctx)); 286b07979f9SJames Wright PetscCall(MatShellSetContextDestroy(J, (PetscErrorCode(*)(void *))OperatorApplyContextDestroy)); 287b07979f9SJames Wright PetscCall(MatShellSetOperation(J, MATOP_MULT, (void (*)(void))MatMult_Ceed)); 288b07979f9SJames Wright PetscCall(MatShellSetOperation(J, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiag_Ceed)); 289e334ad8fSJed Brown PetscCall(MatSetUp(J)); 290e334ad8fSJed Brown } 291e334ad8fSJed Brown if (!J_pre_is_shell) { 2922b730f8bSJeremy L Thompson PetscCall(FormPreallocation(user, user->app_ctx->pmat_pbdiagonal, J_pre, &user->coo_values_pmat)); 293544be873SJed Brown } 294d6bf345cSJed Brown if (J != J_pre && !J_is_shell && !J_is_mffd) { 295544be873SJed Brown PetscCall(FormPreallocation(user, PETSC_FALSE, J, &user->coo_values_amat)); 296544be873SJed Brown } 297e334ad8fSJed Brown user->matrices_set_up = true; 298e334ad8fSJed Brown } 299e334ad8fSJed Brown if (!J_pre_is_shell) { 3002b730f8bSJeremy L Thompson PetscCall(FormSetValues(user, user->app_ctx->pmat_pbdiagonal, J_pre, user->coo_values_pmat)); 301e334ad8fSJed Brown } 302d6bf345cSJed Brown if (user->coo_values_amat) { 303d6bf345cSJed Brown PetscCall(FormSetValues(user, PETSC_FALSE, J, user->coo_values_amat)); 304d6bf345cSJed Brown } else if (J_is_mffd) { 305d6bf345cSJed Brown PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 306d6bf345cSJed Brown PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 307d6bf345cSJed Brown } 308ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 309e334ad8fSJed Brown } 310e334ad8fSJed Brown 3112b730f8bSJeremy L Thompson PetscErrorCode WriteOutput(User user, Vec Q, PetscInt step_no, PetscScalar time) { 31277841947SLeila Ghaffari Vec Q_loc; 31377841947SLeila Ghaffari char file_path[PETSC_MAX_PATH_LEN]; 31477841947SLeila Ghaffari PetscViewer viewer; 31577841947SLeila Ghaffari 316f17d818dSJames Wright PetscFunctionBeginUser; 31737cbb16aSJed Brown if (user->app_ctx->checkpoint_vtk) { 31877841947SLeila Ghaffari // Set up output 319f14660a4SJames Wright PetscCall(DMGetLocalVector(user->dm, &Q_loc)); 320f14660a4SJames Wright PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec")); 321f14660a4SJames Wright PetscCall(VecZeroEntries(Q_loc)); 322f14660a4SJames Wright PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc)); 32377841947SLeila Ghaffari 32477841947SLeila Ghaffari // Output 32537cbb16aSJed Brown PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no)); 326f14660a4SJames Wright 3272b730f8bSJeremy L Thompson PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path, FILE_MODE_WRITE, &viewer)); 328f14660a4SJames Wright PetscCall(VecView(Q_loc, viewer)); 329f14660a4SJames Wright PetscCall(PetscViewerDestroy(&viewer)); 33077841947SLeila Ghaffari if (user->dm_viz) { 33177841947SLeila Ghaffari Vec Q_refined, Q_refined_loc; 33277841947SLeila Ghaffari char file_path_refined[PETSC_MAX_PATH_LEN]; 33377841947SLeila Ghaffari PetscViewer viewer_refined; 33477841947SLeila Ghaffari 335f14660a4SJames Wright PetscCall(DMGetGlobalVector(user->dm_viz, &Q_refined)); 336f14660a4SJames Wright PetscCall(DMGetLocalVector(user->dm_viz, &Q_refined_loc)); 337f14660a4SJames Wright PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined")); 338f14660a4SJames Wright 339f14660a4SJames Wright PetscCall(MatInterpolate(user->interp_viz, Q, Q_refined)); 340f14660a4SJames Wright PetscCall(VecZeroEntries(Q_refined_loc)); 3412b730f8bSJeremy L Thompson PetscCall(DMGlobalToLocal(user->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc)); 342f14660a4SJames Wright 34337cbb16aSJed Brown PetscCall( 34437cbb16aSJed Brown PetscSNPrintf(file_path_refined, sizeof file_path_refined, "%s/nsrefined-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no)); 345f14660a4SJames Wright 3462b730f8bSJeremy L Thompson PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined), file_path_refined, FILE_MODE_WRITE, &viewer_refined)); 347f14660a4SJames Wright PetscCall(VecView(Q_refined_loc, viewer_refined)); 348f14660a4SJames Wright PetscCall(DMRestoreLocalVector(user->dm_viz, &Q_refined_loc)); 349f14660a4SJames Wright PetscCall(DMRestoreGlobalVector(user->dm_viz, &Q_refined)); 350f14660a4SJames Wright PetscCall(PetscViewerDestroy(&viewer_refined)); 35177841947SLeila Ghaffari } 352f14660a4SJames Wright PetscCall(DMRestoreLocalVector(user->dm, &Q_loc)); 35337cbb16aSJed Brown } 35477841947SLeila Ghaffari 35577841947SLeila Ghaffari // Save data in a binary file for continuation of simulations 35669293791SJames Wright if (user->app_ctx->add_stepnum2bin) { 35737cbb16aSJed Brown PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", user->app_ctx->output_dir, step_no)); 35869293791SJames Wright } else { 3592b730f8bSJeremy L Thompson PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", user->app_ctx->output_dir)); 36069293791SJames Wright } 3612b730f8bSJeremy L Thompson PetscCall(PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, &viewer)); 362f14660a4SJames Wright 3630de6a49fSJames Wright PetscInt32 token = PetscDefined(USE_64BIT_INDICES) ? FLUIDS_FILE_TOKEN_64 : FLUIDS_FILE_TOKEN_32; 3640de6a49fSJames Wright PetscCall(PetscViewerBinaryWrite(viewer, &token, 1, PETSC_INT32)); 3654de8550aSJed Brown PetscCall(PetscViewerBinaryWrite(viewer, &step_no, 1, PETSC_INT)); 3664de8550aSJed Brown time /= user->units->second; // Dimensionalize time back 3674de8550aSJed Brown PetscCall(PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL)); 368f14660a4SJames Wright PetscCall(VecView(Q, viewer)); 369f14660a4SJames Wright PetscCall(PetscViewerDestroy(&viewer)); 370ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 371f14660a4SJames Wright } 372f14660a4SJames Wright 373ca69d878SAdeleke O. Bankole // CSV Monitor 374ca69d878SAdeleke O. Bankole PetscErrorCode TSMonitor_WallForce(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) { 375ca69d878SAdeleke O. Bankole User user = ctx; 376ca69d878SAdeleke O. Bankole Vec G_loc; 377ca69d878SAdeleke O. Bankole PetscInt num_wall = user->app_ctx->wall_forces.num_wall, dim = 3; 378ca69d878SAdeleke O. Bankole const PetscInt *walls = user->app_ctx->wall_forces.walls; 379ca69d878SAdeleke O. Bankole PetscViewer viewer = user->app_ctx->wall_forces.viewer; 380ca69d878SAdeleke O. Bankole PetscViewerFormat format = user->app_ctx->wall_forces.viewer_format; 381ca69d878SAdeleke O. Bankole PetscScalar *reaction_force; 382ca69d878SAdeleke O. Bankole PetscBool iascii; 383ca69d878SAdeleke O. Bankole 384ca69d878SAdeleke O. Bankole PetscFunctionBeginUser; 385ee4ca9cbSJames Wright if (!viewer) PetscFunctionReturn(PETSC_SUCCESS); 386ca69d878SAdeleke O. Bankole PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); 387ca69d878SAdeleke O. Bankole PetscCall(PetscMalloc1(num_wall * dim, &reaction_force)); 388ca69d878SAdeleke O. Bankole PetscCall(Surface_Forces_NS(user->dm, G_loc, num_wall, walls, reaction_force)); 389ca69d878SAdeleke O. Bankole PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); 390ca69d878SAdeleke O. Bankole 391ca69d878SAdeleke O. Bankole PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 392ca69d878SAdeleke O. Bankole 393ca69d878SAdeleke O. Bankole if (iascii) { 394ca69d878SAdeleke O. Bankole if (format == PETSC_VIEWER_ASCII_CSV && !user->app_ctx->wall_forces.header_written) { 395ca69d878SAdeleke O. Bankole PetscCall(PetscViewerASCIIPrintf(viewer, "Step,Time,Wall,ForceX,ForceY,ForceZ\n")); 396ca69d878SAdeleke O. Bankole user->app_ctx->wall_forces.header_written = PETSC_TRUE; 397ca69d878SAdeleke O. Bankole } 398ca69d878SAdeleke O. Bankole for (PetscInt w = 0; w < num_wall; w++) { 399ca69d878SAdeleke O. Bankole PetscInt wall = walls[w]; 400ca69d878SAdeleke O. Bankole if (format == PETSC_VIEWER_ASCII_CSV) { 401ca69d878SAdeleke O. Bankole PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT ",%g,%" PetscInt_FMT ",%g,%g,%g\n", step_no, time, wall, 402ca69d878SAdeleke O. Bankole reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2])); 403ca69d878SAdeleke O. Bankole 404ca69d878SAdeleke O. Bankole } else { 405ca69d878SAdeleke O. Bankole PetscCall(PetscViewerASCIIPrintf(viewer, "Wall %" PetscInt_FMT " Forces: Force_x = %12g, Force_y = %12g, Force_z = %12g\n", wall, 406ca69d878SAdeleke O. Bankole reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2])); 407ca69d878SAdeleke O. Bankole } 408ca69d878SAdeleke O. Bankole } 409ca69d878SAdeleke O. Bankole } 410ca69d878SAdeleke O. Bankole PetscCall(PetscFree(reaction_force)); 411ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 412ca69d878SAdeleke O. Bankole } 413ca69d878SAdeleke O. Bankole 414f14660a4SJames Wright // User provided TS Monitor 4152b730f8bSJeremy L Thompson PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) { 416f14660a4SJames Wright User user = ctx; 417f14660a4SJames Wright 418f17d818dSJames Wright PetscFunctionBeginUser; 41937cbb16aSJed Brown // Print every 'checkpoint_interval' steps 420894de27cSJames Wright if (user->app_ctx->checkpoint_interval <= 0 || step_no % user->app_ctx->checkpoint_interval != 0 || 42149aac155SJeremy L Thompson (user->app_ctx->cont_steps == step_no && step_no != 0)) { 422ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 42349aac155SJeremy L Thompson } 424f14660a4SJames Wright 425f14660a4SJames Wright PetscCall(WriteOutput(user, Q, step_no, time)); 426ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 42777841947SLeila Ghaffari } 42877841947SLeila Ghaffari 42977841947SLeila Ghaffari // TS: Create, setup, and solve 4302b730f8bSJeremy L Thompson PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, Vec *Q, PetscScalar *f_time, TS *ts) { 43177841947SLeila Ghaffari MPI_Comm comm = user->comm; 43277841947SLeila Ghaffari TSAdapt adapt; 43377841947SLeila Ghaffari PetscScalar final_time; 43477841947SLeila Ghaffari 435f17d818dSJames Wright PetscFunctionBeginUser; 4362b730f8bSJeremy L Thompson PetscCall(TSCreate(comm, ts)); 4372b730f8bSJeremy L Thompson PetscCall(TSSetDM(*ts, dm)); 43877841947SLeila Ghaffari if (phys->implicit) { 4392b730f8bSJeremy L Thompson PetscCall(TSSetType(*ts, TSBDF)); 44077841947SLeila Ghaffari if (user->op_ifunction) { 4412b730f8bSJeremy L Thompson PetscCall(TSSetIFunction(*ts, NULL, IFunction_NS, &user)); 44277841947SLeila Ghaffari } else { // Implicit integrators can fall back to using an RHSFunction 4432b730f8bSJeremy L Thompson PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user)); 44477841947SLeila Ghaffari } 445e334ad8fSJed Brown if (user->op_ijacobian) { 4462b730f8bSJeremy L Thompson PetscCall(DMTSSetIJacobian(dm, FormIJacobian_NS, &user)); 447544be873SJed Brown if (app_ctx->amat_type) { 448544be873SJed Brown Mat Pmat, Amat; 4492b730f8bSJeremy L Thompson PetscCall(DMCreateMatrix(dm, &Pmat)); 4502b730f8bSJeremy L Thompson PetscCall(DMSetMatType(dm, app_ctx->amat_type)); 4512b730f8bSJeremy L Thompson PetscCall(DMCreateMatrix(dm, &Amat)); 4522b730f8bSJeremy L Thompson PetscCall(TSSetIJacobian(*ts, Amat, Pmat, NULL, NULL)); 4532b730f8bSJeremy L Thompson PetscCall(MatDestroy(&Amat)); 4542b730f8bSJeremy L Thompson PetscCall(MatDestroy(&Pmat)); 455544be873SJed Brown } 456e334ad8fSJed Brown } 45777841947SLeila Ghaffari } else { 4589ad5e8e4SJames Wright PetscCheck(user->op_rhs_ctx, comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction"); 4592b730f8bSJeremy L Thompson PetscCall(TSSetType(*ts, TSRK)); 4602b730f8bSJeremy L Thompson PetscCall(TSRKSetType(*ts, TSRK5F)); 4612b730f8bSJeremy L Thompson PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user)); 46277841947SLeila Ghaffari } 4632b730f8bSJeremy L Thompson PetscCall(TSSetMaxTime(*ts, 500. * user->units->second)); 4642b730f8bSJeremy L Thompson PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER)); 46551b00d91SJames Wright if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(TSSetErrorIfStepFails(*ts, PETSC_FALSE)); 4662b730f8bSJeremy L Thompson PetscCall(TSSetTimeStep(*ts, 1.e-2 * user->units->second)); 4678fb33541SJames Wright if (app_ctx->test_type != TESTTYPE_NONE) { 4682b730f8bSJeremy L Thompson PetscCall(TSSetMaxSteps(*ts, 10)); 4692b730f8bSJeremy L Thompson } 4702b730f8bSJeremy L Thompson PetscCall(TSGetAdapt(*ts, &adapt)); 4712b730f8bSJeremy L Thompson PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * user->units->second, 1.e2 * user->units->second)); 4722b730f8bSJeremy L Thompson PetscCall(TSSetFromOptions(*ts)); 47315637395SJames Wright user->time_bc_set = -1.0; // require all BCs be updated 474d55646a4SJames Wright if (app_ctx->cont_steps) { // continue from previous timestep data 47577841947SLeila Ghaffari PetscInt count; 47677841947SLeila Ghaffari PetscViewer viewer; 4772b730f8bSJeremy L Thompson 4784de8550aSJed Brown if (app_ctx->cont_time <= 0) { // Legacy files did not include step number and time 4792b730f8bSJeremy L Thompson PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_time_file, FILE_MODE_READ, &viewer)); 4804de8550aSJed Brown PetscCall(PetscViewerBinaryRead(viewer, &app_ctx->cont_time, 1, &count, PETSC_REAL)); 4812b730f8bSJeremy L Thompson PetscCall(PetscViewerDestroy(&viewer)); 4824de8550aSJed Brown PetscCheck(app_ctx->cont_steps != -1, comm, PETSC_ERR_ARG_INCOMP, 4834de8550aSJed Brown "-continue step number not specified, but checkpoint file does not contain a step number (likely written by older code version)"); 4844de8550aSJed Brown } 4854de8550aSJed Brown PetscCall(TSSetTime(*ts, app_ctx->cont_time * user->units->second)); 486d8e0aecdSJed Brown PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps)); 48777841947SLeila Ghaffari } 4888fb33541SJames Wright if (app_ctx->test_type == TESTTYPE_NONE) { 4892b730f8bSJeremy L Thompson PetscCall(TSMonitorSet(*ts, TSMonitor_NS, user, NULL)); 4908fb33541SJames Wright } 491ca69d878SAdeleke O. Bankole if (app_ctx->wall_forces.viewer) { 492ca69d878SAdeleke O. Bankole PetscCall(TSMonitorSet(*ts, TSMonitor_WallForce, user, NULL)); 493ca69d878SAdeleke O. Bankole } 494b7d66439SJames Wright if (app_ctx->turb_spanstats_enable) { 495f5452247SJames Wright PetscCall(TSMonitorSet(*ts, TSMonitor_TurbulenceStatistics, user, NULL)); 496495b9769SJames Wright CeedScalar previous_time = app_ctx->cont_time * user->units->second; 497a424bcd0SJames Wright PetscCallCeed(user->ceed, 498a424bcd0SJames Wright CeedOperatorSetContextDouble(user->spanstats.op_stats_collect_ctx->op, user->spanstats.previous_time_label, &previous_time)); 499a175e481SJames Wright } 5004e9802d1SJames Wright if (app_ctx->diff_filter_monitor) PetscCall(TSMonitorSet(*ts, TSMonitor_DifferentialFilter, user, NULL)); 50177841947SLeila Ghaffari 50277841947SLeila Ghaffari // Solve 503d8e0aecdSJed Brown PetscReal start_time; 504d8e0aecdSJed Brown PetscInt start_step; 5052b730f8bSJeremy L Thompson PetscCall(TSGetTime(*ts, &start_time)); 506d8e0aecdSJed Brown PetscCall(TSGetStepNumber(*ts, &start_step)); 5077b39487dSJeremy L Thompson 508711a423aSJed Brown PetscCall(PetscLogDefaultBegin()); // So we can use PetscLogStageGetPerfInfo without -log_view 5097b39487dSJeremy L Thompson PetscPreLoadBegin(PETSC_FALSE, "Fluids Solve"); 5107b39487dSJeremy L Thompson PetscCall(TSSetTime(*ts, start_time)); 511d8e0aecdSJed Brown PetscCall(TSSetStepNumber(*ts, start_step)); 5127b39487dSJeremy L Thompson if (PetscPreLoadingOn) { 5137b39487dSJeremy L Thompson // LCOV_EXCL_START 5147b39487dSJeremy L Thompson SNES snes; 5157b39487dSJeremy L Thompson Vec Q_preload; 5167b39487dSJeremy L Thompson PetscReal rtol; 5177b39487dSJeremy L Thompson PetscCall(VecDuplicate(*Q, &Q_preload)); 5187b39487dSJeremy L Thompson PetscCall(VecCopy(*Q, Q_preload)); 5197b39487dSJeremy L Thompson PetscCall(TSGetSNES(*ts, &snes)); 5207b39487dSJeremy L Thompson PetscCall(SNESGetTolerances(snes, NULL, &rtol, NULL, NULL, NULL)); 5212b730f8bSJeremy L Thompson PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, .99, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 522fc4c3d69SJames Wright PetscCall(TSSetSolution(*ts, Q_preload)); 5237b39487dSJeremy L Thompson PetscCall(TSStep(*ts)); 5242b730f8bSJeremy L Thompson PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, rtol, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 5257b39487dSJeremy L Thompson PetscCall(VecDestroy(&Q_preload)); 5267b39487dSJeremy L Thompson // LCOV_EXCL_STOP 5277b39487dSJeremy L Thompson } else { 5282b730f8bSJeremy L Thompson PetscCall(PetscBarrier((PetscObject)*ts)); 5292b730f8bSJeremy L Thompson PetscCall(TSSolve(*ts, *Q)); 5307b39487dSJeremy L Thompson } 5317b39487dSJeremy L Thompson PetscPreLoadEnd(); 5327b39487dSJeremy L Thompson 5337b39487dSJeremy L Thompson PetscCall(TSGetSolveTime(*ts, &final_time)); 53477841947SLeila Ghaffari *f_time = final_time; 5357b39487dSJeremy L Thompson 5368fb33541SJames Wright if (app_ctx->test_type == TESTTYPE_NONE) { 537f14660a4SJames Wright PetscInt step_no; 538f14660a4SJames Wright PetscCall(TSGetStepNumber(*ts, &step_no)); 539a175e481SJames Wright if (user->app_ctx->checkpoint_interval > 0 || user->app_ctx->checkpoint_interval == -1) { 540f14660a4SJames Wright PetscCall(WriteOutput(user, *Q, step_no, final_time)); 541f14660a4SJames Wright } 542f14660a4SJames Wright 54375d1798cSJames Wright PetscLogStage stage_id; 544711a423aSJed Brown PetscEventPerfInfo stage_perf; 5457b39487dSJeremy L Thompson 5467b39487dSJeremy L Thompson PetscCall(PetscLogStageGetId("Fluids Solve", &stage_id)); 547711a423aSJed Brown PetscCall(PetscLogStageGetPerfInfo(stage_id, &stage_perf)); 548711a423aSJed Brown PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time taken for solution (sec): %g\n", stage_perf.time)); 54977841947SLeila Ghaffari } 550ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 55177841947SLeila Ghaffari } 552