1727da7e7SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 2727da7e7SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3a515125bSLeila Ghaffari // 4727da7e7SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 5a515125bSLeila Ghaffari // 6727da7e7SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 7a515125bSLeila Ghaffari 8a515125bSLeila Ghaffari /// @file 9a515125bSLeila Ghaffari /// Time-stepping functions for Navier-Stokes example using PETSc 10a515125bSLeila Ghaffari 11e419654dSJeremy L Thompson #include <ceed.h> 12e419654dSJeremy L Thompson #include <petscdmplex.h> 13e419654dSJeremy L Thompson #include <petscts.h> 14e419654dSJeremy L Thompson 15a515125bSLeila Ghaffari #include "../navierstokes.h" 16c5e9980aSAdeleke O. Bankole #include "../qfunctions/newtonian_state.h" 17a515125bSLeila Ghaffari 18a515125bSLeila Ghaffari // Compute mass matrix for explicit scheme 192b916ea7SJeremy L Thompson PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm, CeedData ceed_data, Vec M) { 20a515125bSLeila Ghaffari CeedQFunction qf_mass; 21a515125bSLeila Ghaffari CeedOperator op_mass; 220143e3daSJames Wright OperatorApplyContext op_mass_ctx; 230143e3daSJames Wright Vec Ones_loc; 24a515125bSLeila Ghaffari CeedInt num_comp_q, q_data_size; 25a515125bSLeila Ghaffari 26*06f41313SJames Wright PetscFunctionBeginUser; 27a515125bSLeila Ghaffari // CEED Restriction 28b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q)); 29b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &q_data_size)); 30a515125bSLeila Ghaffari 31a515125bSLeila Ghaffari // CEED QFunction 329f59f36eSJames Wright PetscCall(CreateMassQFunction(ceed, num_comp_q, q_data_size, &qf_mass)); 33a515125bSLeila Ghaffari 34a515125bSLeila Ghaffari // CEED Operator 35b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass)); 36b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 37b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data)); 38b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 39a515125bSLeila Ghaffari 400143e3daSJames Wright PetscCall(OperatorApplyContextCreate(NULL, dm, ceed, op_mass, NULL, NULL, NULL, NULL, &op_mass_ctx)); 41a515125bSLeila Ghaffari 420143e3daSJames Wright PetscCall(DMGetLocalVector(dm, &Ones_loc)); 430143e3daSJames Wright PetscCall(VecSet(Ones_loc, 1)); 440143e3daSJames Wright PetscCall(ApplyCeedOperatorLocalToGlobal(Ones_loc, M, op_mass_ctx)); 45a515125bSLeila Ghaffari 46a515125bSLeila Ghaffari // Invert diagonally lumped mass vector for RHS function 472b916ea7SJeremy L Thompson PetscCall(VecReciprocal(M)); 48a515125bSLeila Ghaffari 49a515125bSLeila Ghaffari // Cleanup 500143e3daSJames Wright PetscCall(OperatorApplyContextDestroy(op_mass_ctx)); 510143e3daSJames Wright PetscCall(DMRestoreLocalVector(dm, &Ones_loc)); 52b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 53b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass)); 54d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 55a515125bSLeila Ghaffari } 56a515125bSLeila Ghaffari 57c996854bSJames Wright // Insert Boundary values if it's a new time 58c996854bSJames Wright PetscErrorCode UpdateBoundaryValues(User user, Vec Q_loc, PetscReal t) { 59c996854bSJames Wright PetscFunctionBeginUser; 60c996854bSJames Wright if (user->time_bc_set != t) { 61c996854bSJames Wright PetscCall(DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL)); 62c996854bSJames Wright user->time_bc_set = t; 63c996854bSJames Wright } 64d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 65c996854bSJames Wright } 66c996854bSJames Wright 6791f639d2SJames Wright // @brief Update the context label value to new value if necessary. 6891f639d2SJames Wright // @note This only supports labels with scalar label values (ie. not arrays) 69b4c37c5cSJames Wright PetscErrorCode UpdateContextLabel(Ceed ceed, MPI_Comm comm, PetscScalar update_value, CeedOperator op, CeedContextFieldLabel label) { 7091f639d2SJames Wright PetscScalar label_value; 71c2cb7fc8SJames Wright 7291f639d2SJames Wright PetscFunctionBeginUser; 7391f639d2SJames Wright PetscCheck(label, comm, PETSC_ERR_ARG_BADPTR, "Label should be non-NULL"); 7491f639d2SJames Wright 7591f639d2SJames Wright { 7691f639d2SJames Wright size_t num_elements; 7791f639d2SJames Wright const PetscScalar *label_values; 78b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorGetContextDoubleRead(op, label, &num_elements, &label_values)); 7991f639d2SJames Wright PetscCheck(num_elements == 1, comm, PETSC_ERR_SUP, "%s does not support labels with more than 1 value. Label has %zu values", __func__, 8091f639d2SJames Wright num_elements); 8191f639d2SJames Wright label_value = *label_values; 82b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorRestoreContextDoubleRead(op, label, &label_values)); 83c2cb7fc8SJames Wright } 8491f639d2SJames Wright 8591f639d2SJames Wright if (label_value != update_value) { 86b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetContextDouble(op, label, &update_value)); 87c2cb7fc8SJames Wright } 88d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 89c2cb7fc8SJames Wright } 90c2cb7fc8SJames Wright 91a515125bSLeila Ghaffari // RHS (Explicit time-stepper) function setup 92a515125bSLeila Ghaffari // This is the RHS of the ODE, given as u_t = G(t,u) 93a515125bSLeila Ghaffari // This function takes in a state vector Q and writes into G 94a515125bSLeila Ghaffari PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) { 95a515125bSLeila Ghaffari User user = *(User *)user_data; 9691f639d2SJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)ts); 97fd969b44SJames Wright PetscScalar dt; 98da5fe0e4SJames Wright Vec Q_loc = user->Q_loc; 99a515125bSLeila Ghaffari 100*06f41313SJames Wright PetscFunctionBeginUser; 101e2f84137SJeremy L Thompson // Update time dependent data 102c996854bSJames Wright PetscCall(UpdateBoundaryValues(user, Q_loc, t)); 103b4c37c5cSJames Wright if (user->phys->solution_time_label) PetscCall(UpdateContextLabel(user->ceed, comm, t, user->op_rhs_ctx->op, user->phys->solution_time_label)); 1042b916ea7SJeremy L Thompson PetscCall(TSGetTimeStep(ts, &dt)); 105b4c37c5cSJames Wright if (user->phys->timestep_size_label) PetscCall(UpdateContextLabel(user->ceed, comm, dt, user->op_rhs_ctx->op, user->phys->timestep_size_label)); 106a515125bSLeila Ghaffari 107da5fe0e4SJames Wright PetscCall(ApplyCeedOperatorGlobalToGlobal(Q, G, user->op_rhs_ctx)); 108a515125bSLeila Ghaffari 109f3fcf8f4SJames Wright // Inverse of the lumped mass matrix 110cc9aa765SJames Wright PetscCall(VecPointwiseMult(G, G, user->M_inv)); 111d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 112a515125bSLeila Ghaffari } 113a515125bSLeila Ghaffari 114c5e9980aSAdeleke O. Bankole // Surface forces function setup 115c5e9980aSAdeleke O. Bankole static PetscErrorCode Surface_Forces_NS(DM dm, Vec G_loc, PetscInt num_walls, const PetscInt walls[], PetscScalar *reaction_force) { 116c5e9980aSAdeleke O. Bankole DMLabel face_label; 117c5e9980aSAdeleke O. Bankole const PetscScalar *g; 1182004e3acSAdeleke O. Bankole PetscInt dof, dim = 3; 119c5e9980aSAdeleke O. Bankole MPI_Comm comm; 1202004e3acSAdeleke O. Bankole PetscSection s; 121c5e9980aSAdeleke O. Bankole 122c5e9980aSAdeleke O. Bankole PetscFunctionBeginUser; 123c5e9980aSAdeleke O. Bankole PetscCall(PetscArrayzero(reaction_force, num_walls * dim)); 124c5e9980aSAdeleke O. Bankole PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 125c5e9980aSAdeleke O. Bankole PetscCall(DMGetLabel(dm, "Face Sets", &face_label)); 126c5e9980aSAdeleke O. Bankole PetscCall(VecGetArrayRead(G_loc, &g)); 127c5e9980aSAdeleke O. Bankole for (PetscInt w = 0; w < num_walls; w++) { 128c5e9980aSAdeleke O. Bankole const PetscInt wall = walls[w]; 129c5e9980aSAdeleke O. Bankole IS wall_is; 1302004e3acSAdeleke O. Bankole PetscCall(DMGetLocalSection(dm, &s)); 131c5e9980aSAdeleke O. Bankole PetscCall(DMLabelGetStratumIS(face_label, wall, &wall_is)); 132c5e9980aSAdeleke O. Bankole if (wall_is) { // There exist such points on this process 133c5e9980aSAdeleke O. Bankole PetscInt num_points; 1342004e3acSAdeleke O. Bankole PetscInt num_comp = 0; 135c5e9980aSAdeleke O. Bankole const PetscInt *points; 1362004e3acSAdeleke O. Bankole PetscCall(PetscSectionGetFieldComponents(s, 0, &num_comp)); 137c5e9980aSAdeleke O. Bankole PetscCall(ISGetSize(wall_is, &num_points)); 138c5e9980aSAdeleke O. Bankole PetscCall(ISGetIndices(wall_is, &points)); 139c5e9980aSAdeleke O. Bankole for (PetscInt i = 0; i < num_points; i++) { 140c5e9980aSAdeleke O. Bankole const PetscInt p = points[i]; 141c5e9980aSAdeleke O. Bankole const StateConservative *r; 142c5e9980aSAdeleke O. Bankole PetscCall(DMPlexPointLocalRead(dm, p, g, &r)); 1432004e3acSAdeleke O. Bankole PetscCall(PetscSectionGetDof(s, p, &dof)); 1442004e3acSAdeleke O. Bankole for (PetscInt node = 0; node < dof / num_comp; node++) { 145c5e9980aSAdeleke O. Bankole for (PetscInt j = 0; j < 3; j++) { 1462004e3acSAdeleke O. Bankole reaction_force[w * dim + j] -= r[node].momentum[j]; 1472004e3acSAdeleke O. Bankole } 148c5e9980aSAdeleke O. Bankole } 149c5e9980aSAdeleke O. Bankole } 150c5e9980aSAdeleke O. Bankole PetscCall(ISRestoreIndices(wall_is, &points)); 151c5e9980aSAdeleke O. Bankole } 152c5e9980aSAdeleke O. Bankole PetscCall(ISDestroy(&wall_is)); 153c5e9980aSAdeleke O. Bankole } 154c5e9980aSAdeleke O. Bankole PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, reaction_force, dim * num_walls, MPIU_SCALAR, MPI_SUM, comm)); 155c5e9980aSAdeleke O. Bankole // Restore Vectors 156c5e9980aSAdeleke O. Bankole PetscCall(VecRestoreArrayRead(G_loc, &g)); 157d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 158c5e9980aSAdeleke O. Bankole } 159c5e9980aSAdeleke O. Bankole 160a515125bSLeila Ghaffari // Implicit time-stepper function setup 1612b916ea7SJeremy L Thompson PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) { 162a515125bSLeila Ghaffari User user = *(User *)user_data; 16391f639d2SJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)ts); 164fd969b44SJames Wright PetscScalar dt; 165e2f84137SJeremy L Thompson Vec Q_loc = user->Q_loc, Q_dot_loc = user->Q_dot_loc, G_loc; 166a515125bSLeila Ghaffari PetscMemType q_mem_type, q_dot_mem_type, g_mem_type; 167a515125bSLeila Ghaffari 168*06f41313SJames Wright PetscFunctionBeginUser; 169e2f84137SJeremy L Thompson // Get local vectors 170c5e9980aSAdeleke O. Bankole PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); 171e2f84137SJeremy L Thompson 172e2f84137SJeremy L Thompson // Update time dependent data 173c996854bSJames Wright PetscCall(UpdateBoundaryValues(user, Q_loc, t)); 174b4c37c5cSJames Wright if (user->phys->solution_time_label) PetscCall(UpdateContextLabel(user->ceed, comm, t, user->op_ifunction, user->phys->solution_time_label)); 1752b916ea7SJeremy L Thompson PetscCall(TSGetTimeStep(ts, &dt)); 176b4c37c5cSJames Wright if (user->phys->timestep_size_label) PetscCall(UpdateContextLabel(user->ceed, comm, dt, user->op_ifunction, user->phys->timestep_size_label)); 177a515125bSLeila Ghaffari 178a515125bSLeila Ghaffari // Global-to-local 17906108310SJames Wright PetscCall(DMGlobalToLocalBegin(user->dm, Q, INSERT_VALUES, Q_loc)); 18006108310SJames Wright PetscCall(DMGlobalToLocalBegin(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc)); 18106108310SJames Wright PetscCall(DMGlobalToLocalEnd(user->dm, Q, INSERT_VALUES, Q_loc)); 18206108310SJames Wright PetscCall(DMGlobalToLocalEnd(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc)); 183a515125bSLeila Ghaffari 184a515125bSLeila Ghaffari // Place PETSc vectors in CEED vectors 185fd969b44SJames Wright PetscCall(VecReadP2C(Q_loc, &q_mem_type, user->q_ceed)); 186fd969b44SJames Wright PetscCall(VecReadP2C(Q_dot_loc, &q_dot_mem_type, user->q_dot_ceed)); 187fd969b44SJames Wright PetscCall(VecP2C(G_loc, &g_mem_type, user->g_ceed)); 188a515125bSLeila Ghaffari 189a515125bSLeila Ghaffari // Apply CEED operator 1907eedc94cSJames Wright PetscCall(PetscLogEventBegin(FLUIDS_CeedOperatorApply, Q, G, 0, 0)); 1917eedc94cSJames Wright PetscCall(PetscLogGpuTimeBegin()); 192b4c37c5cSJames Wright PetscCallCeed(user->ceed, CeedOperatorApply(user->op_ifunction, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE)); 1937eedc94cSJames Wright PetscCall(PetscLogGpuTimeEnd()); 1947eedc94cSJames Wright PetscCall(PetscLogEventEnd(FLUIDS_CeedOperatorApply, Q, G, 0, 0)); 195a515125bSLeila Ghaffari 196a515125bSLeila Ghaffari // Restore vectors 197fd969b44SJames Wright PetscCall(VecReadC2P(user->q_ceed, q_mem_type, Q_loc)); 198fd969b44SJames Wright PetscCall(VecReadC2P(user->q_dot_ceed, q_dot_mem_type, Q_dot_loc)); 199fd969b44SJames Wright PetscCall(VecC2P(user->g_ceed, g_mem_type, G_loc)); 200a515125bSLeila Ghaffari 20101ab89c1SJames Wright if (user->app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) { 20242454adaSJames Wright PetscCall(SgsDDModelApplyIFunction(user, Q_loc, G_loc)); 20301ab89c1SJames Wright } 2049c678832SJames Wright 205a515125bSLeila Ghaffari // Local-to-Global 2062b916ea7SJeremy L Thompson PetscCall(VecZeroEntries(G)); 2072b916ea7SJeremy L Thompson PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G)); 208a515125bSLeila Ghaffari 209a515125bSLeila Ghaffari // Restore vectors 210c5e9980aSAdeleke O. Bankole PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); 211d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 212a515125bSLeila Ghaffari } 213a515125bSLeila Ghaffari 2142b916ea7SJeremy L Thompson static PetscErrorCode FormPreallocation(User user, PetscBool pbdiagonal, Mat J, CeedVector *coo_values) { 215b107fddaSJed Brown PetscCount ncoo; 216defe8520SJames Wright PetscInt *rows_petsc, *cols_petsc; 217b107fddaSJed Brown 218b107fddaSJed Brown PetscFunctionBeginUser; 219b107fddaSJed Brown if (pbdiagonal) { 220b107fddaSJed Brown CeedSize l_size; 221b4c37c5cSJames Wright PetscCallCeed(user->ceed, CeedOperatorGetActiveVectorLengths(user->op_ijacobian, &l_size, NULL)); 222b107fddaSJed Brown ncoo = l_size * 5; 223defe8520SJames Wright rows_petsc = malloc(ncoo * sizeof(rows_petsc)); 224defe8520SJames Wright cols_petsc = malloc(ncoo * sizeof(cols_petsc)); 225b107fddaSJed Brown for (PetscCount n = 0; n < l_size / 5; n++) { 226b107fddaSJed Brown for (PetscInt i = 0; i < 5; i++) { 227b107fddaSJed Brown for (PetscInt j = 0; j < 5; j++) { 228defe8520SJames Wright rows_petsc[(n * 5 + i) * 5 + j] = n * 5 + i; 229defe8520SJames Wright cols_petsc[(n * 5 + i) * 5 + j] = n * 5 + j; 230b107fddaSJed Brown } 231b107fddaSJed Brown } 232b107fddaSJed Brown } 233b107fddaSJed Brown } else { 234defe8520SJames Wright CeedInt *rows_ceed, *cols_ceed; 235b4c37c5cSJames Wright PetscCallCeed(user->ceed, CeedOperatorLinearAssembleSymbolic(user->op_ijacobian, &ncoo, &rows_ceed, &cols_ceed)); 236defe8520SJames Wright PetscCall(IntArrayC2P(ncoo, &rows_ceed, &rows_petsc)); 237defe8520SJames Wright PetscCall(IntArrayC2P(ncoo, &cols_ceed, &cols_petsc)); 238b107fddaSJed Brown } 239defe8520SJames Wright PetscCall(MatSetPreallocationCOOLocal(J, ncoo, rows_petsc, cols_petsc)); 240defe8520SJames Wright free(rows_petsc); 241defe8520SJames Wright free(cols_petsc); 242b4c37c5cSJames Wright PetscCallCeed(user->ceed, CeedVectorCreate(user->ceed, ncoo, coo_values)); 243d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 244b107fddaSJed Brown } 245b107fddaSJed Brown 2462b916ea7SJeremy L Thompson static PetscErrorCode FormSetValues(User user, PetscBool pbdiagonal, Mat J, CeedVector coo_values) { 247b107fddaSJed Brown CeedMemType mem_type = CEED_MEM_HOST; 248b107fddaSJed Brown const PetscScalar *values; 249b107fddaSJed Brown MatType mat_type; 250b107fddaSJed Brown 251b107fddaSJed Brown PetscFunctionBeginUser; 252b107fddaSJed Brown PetscCall(MatGetType(J, &mat_type)); 2532b916ea7SJeremy L Thompson if (strstr(mat_type, "kokkos") || strstr(mat_type, "cusparse")) mem_type = CEED_MEM_DEVICE; 254cb315d14SJames Wright if (pbdiagonal) { 2557eedc94cSJames Wright PetscCall(PetscLogEventBegin(FLUIDS_CeedOperatorAssemblePointBlockDiagonal, J, 0, 0, 0)); 2567eedc94cSJames Wright PetscCall(PetscLogGpuTimeBegin()); 257b4c37c5cSJames Wright PetscCallCeed(user->ceed, CeedOperatorLinearAssemblePointBlockDiagonal(user->op_ijacobian, coo_values, CEED_REQUEST_IMMEDIATE)); 2587eedc94cSJames Wright PetscCall(PetscLogGpuTimeEnd()); 2597eedc94cSJames Wright PetscCall(PetscLogEventEnd(FLUIDS_CeedOperatorAssemblePointBlockDiagonal, J, 0, 0, 0)); 260b107fddaSJed Brown } else { 2617eedc94cSJames Wright PetscCall(PetscLogEventBegin(FLUIDS_CeedOperatorAssemble, J, 0, 0, 0)); 2627eedc94cSJames Wright PetscCall(PetscLogGpuTimeBegin()); 263b4c37c5cSJames Wright PetscCallCeed(user->ceed, CeedOperatorLinearAssemble(user->op_ijacobian, coo_values)); 2647eedc94cSJames Wright PetscCall(PetscLogGpuTimeEnd()); 2657eedc94cSJames Wright PetscCall(PetscLogEventEnd(FLUIDS_CeedOperatorAssemble, J, 0, 0, 0)); 266b107fddaSJed Brown } 267b4c37c5cSJames Wright PetscCallCeed(user->ceed, CeedVectorGetArrayRead(coo_values, mem_type, &values)); 268b107fddaSJed Brown PetscCall(MatSetValuesCOO(J, values, INSERT_VALUES)); 269b4c37c5cSJames Wright PetscCallCeed(user->ceed, CeedVectorRestoreArrayRead(coo_values, &values)); 270d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 271b107fddaSJed Brown } 272b107fddaSJed Brown 2732b916ea7SJeremy L Thompson PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, PetscReal shift, Mat J, Mat J_pre, void *user_data) { 274f0b65372SJed Brown User user = *(User *)user_data; 275b4c37c5cSJames Wright Ceed ceed = user->ceed; 27604855949SJed Brown PetscBool J_is_shell, J_is_mffd, J_pre_is_shell; 277*06f41313SJames Wright 278f0b65372SJed Brown PetscFunctionBeginUser; 279b4c37c5cSJames Wright if (user->phys->ijacobian_time_shift_label) 280b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetContextDouble(user->op_ijacobian, user->phys->ijacobian_time_shift_label, &shift)); 28104855949SJed Brown PetscCall(PetscObjectTypeCompare((PetscObject)J, MATMFFD, &J_is_mffd)); 282f0b65372SJed Brown PetscCall(PetscObjectTypeCompare((PetscObject)J, MATSHELL, &J_is_shell)); 2832b916ea7SJeremy L Thompson PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATSHELL, &J_pre_is_shell)); 284f0b65372SJed Brown if (!user->matrices_set_up) { 285f0b65372SJed Brown if (J_is_shell) { 286f9028c3cSJames Wright OperatorApplyContext op_ijacobian_ctx; 287f9028c3cSJames Wright OperatorApplyContextCreate(user->dm, user->dm, user->ceed, user->op_ijacobian, user->q_ceed, user->g_ceed, user->Q_dot_loc, NULL, 288f9028c3cSJames Wright &op_ijacobian_ctx); 289f9028c3cSJames Wright PetscCall(MatShellSetContext(J, op_ijacobian_ctx)); 290f9028c3cSJames Wright PetscCall(MatShellSetContextDestroy(J, (PetscErrorCode(*)(void *))OperatorApplyContextDestroy)); 291f9028c3cSJames Wright PetscCall(MatShellSetOperation(J, MATOP_MULT, (void (*)(void))MatMult_Ceed)); 292f9028c3cSJames Wright PetscCall(MatShellSetOperation(J, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiag_Ceed)); 293f0b65372SJed Brown PetscCall(MatSetUp(J)); 294f0b65372SJed Brown } 295f0b65372SJed Brown if (!J_pre_is_shell) { 2962b916ea7SJeremy L Thompson PetscCall(FormPreallocation(user, user->app_ctx->pmat_pbdiagonal, J_pre, &user->coo_values_pmat)); 297b107fddaSJed Brown } 29804855949SJed Brown if (J != J_pre && !J_is_shell && !J_is_mffd) { 299b107fddaSJed Brown PetscCall(FormPreallocation(user, PETSC_FALSE, J, &user->coo_values_amat)); 300b107fddaSJed Brown } 301f0b65372SJed Brown user->matrices_set_up = true; 302f0b65372SJed Brown } 303f0b65372SJed Brown if (!J_pre_is_shell) { 3042b916ea7SJeremy L Thompson PetscCall(FormSetValues(user, user->app_ctx->pmat_pbdiagonal, J_pre, user->coo_values_pmat)); 305f0b65372SJed Brown } 30604855949SJed Brown if (user->coo_values_amat) { 30704855949SJed Brown PetscCall(FormSetValues(user, PETSC_FALSE, J, user->coo_values_amat)); 30804855949SJed Brown } else if (J_is_mffd) { 30904855949SJed Brown PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 31004855949SJed Brown PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 31104855949SJed Brown } 312d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 313f0b65372SJed Brown } 314f0b65372SJed Brown 3152b916ea7SJeremy L Thompson PetscErrorCode WriteOutput(User user, Vec Q, PetscInt step_no, PetscScalar time) { 316a515125bSLeila Ghaffari Vec Q_loc; 317a515125bSLeila Ghaffari char file_path[PETSC_MAX_PATH_LEN]; 318a515125bSLeila Ghaffari PetscViewer viewer; 319a515125bSLeila Ghaffari 320*06f41313SJames Wright PetscFunctionBeginUser; 321852e5969SJed Brown if (user->app_ctx->checkpoint_vtk) { 322a515125bSLeila Ghaffari // Set up output 3237538d537SJames Wright PetscCall(DMGetLocalVector(user->dm, &Q_loc)); 3247538d537SJames Wright PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec")); 3257538d537SJames Wright PetscCall(VecZeroEntries(Q_loc)); 3267538d537SJames Wright PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc)); 327a515125bSLeila Ghaffari 328a515125bSLeila Ghaffari // Output 329852e5969SJed Brown PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no)); 3307538d537SJames Wright 3312b916ea7SJeremy L Thompson PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path, FILE_MODE_WRITE, &viewer)); 3327538d537SJames Wright PetscCall(VecView(Q_loc, viewer)); 3337538d537SJames Wright PetscCall(PetscViewerDestroy(&viewer)); 334a515125bSLeila Ghaffari if (user->dm_viz) { 335a515125bSLeila Ghaffari Vec Q_refined, Q_refined_loc; 336a515125bSLeila Ghaffari char file_path_refined[PETSC_MAX_PATH_LEN]; 337a515125bSLeila Ghaffari PetscViewer viewer_refined; 338a515125bSLeila Ghaffari 3397538d537SJames Wright PetscCall(DMGetGlobalVector(user->dm_viz, &Q_refined)); 3407538d537SJames Wright PetscCall(DMGetLocalVector(user->dm_viz, &Q_refined_loc)); 3417538d537SJames Wright PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined")); 3427538d537SJames Wright 3437538d537SJames Wright PetscCall(MatInterpolate(user->interp_viz, Q, Q_refined)); 3447538d537SJames Wright PetscCall(VecZeroEntries(Q_refined_loc)); 3452b916ea7SJeremy L Thompson PetscCall(DMGlobalToLocal(user->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc)); 3467538d537SJames Wright 347852e5969SJed Brown PetscCall( 348852e5969SJed Brown PetscSNPrintf(file_path_refined, sizeof file_path_refined, "%s/nsrefined-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no)); 3497538d537SJames Wright 3502b916ea7SJeremy L Thompson PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined), file_path_refined, FILE_MODE_WRITE, &viewer_refined)); 3517538d537SJames Wright PetscCall(VecView(Q_refined_loc, viewer_refined)); 3527538d537SJames Wright PetscCall(DMRestoreLocalVector(user->dm_viz, &Q_refined_loc)); 3537538d537SJames Wright PetscCall(DMRestoreGlobalVector(user->dm_viz, &Q_refined)); 3547538d537SJames Wright PetscCall(PetscViewerDestroy(&viewer_refined)); 355a515125bSLeila Ghaffari } 3567538d537SJames Wright PetscCall(DMRestoreLocalVector(user->dm, &Q_loc)); 357852e5969SJed Brown } 358a515125bSLeila Ghaffari 359a515125bSLeila Ghaffari // Save data in a binary file for continuation of simulations 36091a36801SJames Wright if (user->app_ctx->add_stepnum2bin) { 361852e5969SJed Brown PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", user->app_ctx->output_dir, step_no)); 36291a36801SJames Wright } else { 3632b916ea7SJeremy L Thompson PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", user->app_ctx->output_dir)); 36491a36801SJames Wright } 3652b916ea7SJeremy L Thompson PetscCall(PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, &viewer)); 3667538d537SJames Wright 367e1233009SJames Wright PetscInt32 token = PetscDefined(USE_64BIT_INDICES) ? FLUIDS_FILE_TOKEN_64 : FLUIDS_FILE_TOKEN_32; 368e1233009SJames Wright PetscCall(PetscViewerBinaryWrite(viewer, &token, 1, PETSC_INT32)); 3699293eaa1SJed Brown PetscCall(PetscViewerBinaryWrite(viewer, &step_no, 1, PETSC_INT)); 3709293eaa1SJed Brown time /= user->units->second; // Dimensionalize time back 3719293eaa1SJed Brown PetscCall(PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL)); 3727538d537SJames Wright PetscCall(VecView(Q, viewer)); 3737538d537SJames Wright PetscCall(PetscViewerDestroy(&viewer)); 374d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3757538d537SJames Wright } 3767538d537SJames Wright 377c5e9980aSAdeleke O. Bankole // CSV Monitor 378c5e9980aSAdeleke O. Bankole PetscErrorCode TSMonitor_WallForce(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) { 379c5e9980aSAdeleke O. Bankole User user = ctx; 380c5e9980aSAdeleke O. Bankole Vec G_loc; 381c5e9980aSAdeleke O. Bankole PetscInt num_wall = user->app_ctx->wall_forces.num_wall, dim = 3; 382c5e9980aSAdeleke O. Bankole const PetscInt *walls = user->app_ctx->wall_forces.walls; 383c5e9980aSAdeleke O. Bankole PetscViewer viewer = user->app_ctx->wall_forces.viewer; 384c5e9980aSAdeleke O. Bankole PetscViewerFormat format = user->app_ctx->wall_forces.viewer_format; 385c5e9980aSAdeleke O. Bankole PetscScalar *reaction_force; 386c5e9980aSAdeleke O. Bankole PetscBool iascii; 387c5e9980aSAdeleke O. Bankole 388c5e9980aSAdeleke O. Bankole PetscFunctionBeginUser; 389d949ddfcSJames Wright if (!viewer) PetscFunctionReturn(PETSC_SUCCESS); 390c5e9980aSAdeleke O. Bankole PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); 391c5e9980aSAdeleke O. Bankole PetscCall(PetscMalloc1(num_wall * dim, &reaction_force)); 392c5e9980aSAdeleke O. Bankole PetscCall(Surface_Forces_NS(user->dm, G_loc, num_wall, walls, reaction_force)); 393c5e9980aSAdeleke O. Bankole PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); 394c5e9980aSAdeleke O. Bankole 395c5e9980aSAdeleke O. Bankole PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 396c5e9980aSAdeleke O. Bankole 397c5e9980aSAdeleke O. Bankole if (iascii) { 398c5e9980aSAdeleke O. Bankole if (format == PETSC_VIEWER_ASCII_CSV && !user->app_ctx->wall_forces.header_written) { 399c5e9980aSAdeleke O. Bankole PetscCall(PetscViewerASCIIPrintf(viewer, "Step,Time,Wall,ForceX,ForceY,ForceZ\n")); 400c5e9980aSAdeleke O. Bankole user->app_ctx->wall_forces.header_written = PETSC_TRUE; 401c5e9980aSAdeleke O. Bankole } 402c5e9980aSAdeleke O. Bankole for (PetscInt w = 0; w < num_wall; w++) { 403c5e9980aSAdeleke O. Bankole PetscInt wall = walls[w]; 404c5e9980aSAdeleke O. Bankole if (format == PETSC_VIEWER_ASCII_CSV) { 405c5e9980aSAdeleke O. Bankole PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT ",%g,%" PetscInt_FMT ",%g,%g,%g\n", step_no, time, wall, 406c5e9980aSAdeleke O. Bankole reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2])); 407c5e9980aSAdeleke O. Bankole 408c5e9980aSAdeleke O. Bankole } else { 409c5e9980aSAdeleke O. Bankole PetscCall(PetscViewerASCIIPrintf(viewer, "Wall %" PetscInt_FMT " Forces: Force_x = %12g, Force_y = %12g, Force_z = %12g\n", wall, 410c5e9980aSAdeleke O. Bankole reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2])); 411c5e9980aSAdeleke O. Bankole } 412c5e9980aSAdeleke O. Bankole } 413c5e9980aSAdeleke O. Bankole } 414c5e9980aSAdeleke O. Bankole PetscCall(PetscFree(reaction_force)); 415d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 416c5e9980aSAdeleke O. Bankole } 417c5e9980aSAdeleke O. Bankole 4187538d537SJames Wright // User provided TS Monitor 4192b916ea7SJeremy L Thompson PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) { 4207538d537SJames Wright User user = ctx; 4217538d537SJames Wright 422*06f41313SJames Wright PetscFunctionBeginUser; 423852e5969SJed Brown // Print every 'checkpoint_interval' steps 424c539088bSJames Wright if (user->app_ctx->checkpoint_interval <= 0 || step_no % user->app_ctx->checkpoint_interval != 0 || 425e419654dSJeremy L Thompson (user->app_ctx->cont_steps == step_no && step_no != 0)) { 426d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 427e419654dSJeremy L Thompson } 4287538d537SJames Wright 4297538d537SJames Wright PetscCall(WriteOutput(user, Q, step_no, time)); 430d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 431a515125bSLeila Ghaffari } 432a515125bSLeila Ghaffari 433a515125bSLeila Ghaffari // TS: Create, setup, and solve 4342b916ea7SJeremy L Thompson PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, Vec *Q, PetscScalar *f_time, TS *ts) { 435a515125bSLeila Ghaffari MPI_Comm comm = user->comm; 436a515125bSLeila Ghaffari TSAdapt adapt; 437a515125bSLeila Ghaffari PetscScalar final_time; 438a515125bSLeila Ghaffari 439*06f41313SJames Wright PetscFunctionBeginUser; 4402b916ea7SJeremy L Thompson PetscCall(TSCreate(comm, ts)); 4412b916ea7SJeremy L Thompson PetscCall(TSSetDM(*ts, dm)); 442a515125bSLeila Ghaffari if (phys->implicit) { 4432b916ea7SJeremy L Thompson PetscCall(TSSetType(*ts, TSBDF)); 444a515125bSLeila Ghaffari if (user->op_ifunction) { 4452b916ea7SJeremy L Thompson PetscCall(TSSetIFunction(*ts, NULL, IFunction_NS, &user)); 446a515125bSLeila Ghaffari } else { // Implicit integrators can fall back to using an RHSFunction 4472b916ea7SJeremy L Thompson PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user)); 448a515125bSLeila Ghaffari } 449f0b65372SJed Brown if (user->op_ijacobian) { 4502b916ea7SJeremy L Thompson PetscCall(DMTSSetIJacobian(dm, FormIJacobian_NS, &user)); 451b107fddaSJed Brown if (app_ctx->amat_type) { 452b107fddaSJed Brown Mat Pmat, Amat; 4532b916ea7SJeremy L Thompson PetscCall(DMCreateMatrix(dm, &Pmat)); 4542b916ea7SJeremy L Thompson PetscCall(DMSetMatType(dm, app_ctx->amat_type)); 4552b916ea7SJeremy L Thompson PetscCall(DMCreateMatrix(dm, &Amat)); 4562b916ea7SJeremy L Thompson PetscCall(TSSetIJacobian(*ts, Amat, Pmat, NULL, NULL)); 4572b916ea7SJeremy L Thompson PetscCall(MatDestroy(&Amat)); 4582b916ea7SJeremy L Thompson PetscCall(MatDestroy(&Pmat)); 459b107fddaSJed Brown } 460f0b65372SJed Brown } 461a515125bSLeila Ghaffari } else { 462da5fe0e4SJames Wright PetscCheck(user->op_rhs_ctx, comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction"); 4632b916ea7SJeremy L Thompson PetscCall(TSSetType(*ts, TSRK)); 4642b916ea7SJeremy L Thompson PetscCall(TSRKSetType(*ts, TSRK5F)); 4652b916ea7SJeremy L Thompson PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user)); 466a515125bSLeila Ghaffari } 4672b916ea7SJeremy L Thompson PetscCall(TSSetMaxTime(*ts, 500. * user->units->second)); 4682b916ea7SJeremy L Thompson PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER)); 46922387d3aSJames Wright if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(TSSetErrorIfStepFails(*ts, PETSC_FALSE)); 4702b916ea7SJeremy L Thompson PetscCall(TSSetTimeStep(*ts, 1.e-2 * user->units->second)); 4710e1e9333SJames Wright if (app_ctx->test_type != TESTTYPE_NONE) { 4722b916ea7SJeremy L Thompson PetscCall(TSSetMaxSteps(*ts, 10)); 4732b916ea7SJeremy L Thompson } 4742b916ea7SJeremy L Thompson PetscCall(TSGetAdapt(*ts, &adapt)); 4752b916ea7SJeremy L Thompson PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * user->units->second, 1.e2 * user->units->second)); 4762b916ea7SJeremy L Thompson PetscCall(TSSetFromOptions(*ts)); 47791f639d2SJames Wright user->time_bc_set = -1.0; // require all BCs be updated 478c26b555cSJames Wright if (app_ctx->cont_steps) { // continue from previous timestep data 479a515125bSLeila Ghaffari PetscInt count; 480a515125bSLeila Ghaffari PetscViewer viewer; 4812b916ea7SJeremy L Thompson 4829293eaa1SJed Brown if (app_ctx->cont_time <= 0) { // Legacy files did not include step number and time 4832b916ea7SJeremy L Thompson PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_time_file, FILE_MODE_READ, &viewer)); 4849293eaa1SJed Brown PetscCall(PetscViewerBinaryRead(viewer, &app_ctx->cont_time, 1, &count, PETSC_REAL)); 4852b916ea7SJeremy L Thompson PetscCall(PetscViewerDestroy(&viewer)); 4869293eaa1SJed Brown PetscCheck(app_ctx->cont_steps != -1, comm, PETSC_ERR_ARG_INCOMP, 4879293eaa1SJed Brown "-continue step number not specified, but checkpoint file does not contain a step number (likely written by older code version)"); 4889293eaa1SJed Brown } 4899293eaa1SJed Brown PetscCall(TSSetTime(*ts, app_ctx->cont_time * user->units->second)); 49074a6f4ddSJed Brown PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps)); 491a515125bSLeila Ghaffari } 4920e1e9333SJames Wright if (app_ctx->test_type == TESTTYPE_NONE) { 4932b916ea7SJeremy L Thompson PetscCall(TSMonitorSet(*ts, TSMonitor_NS, user, NULL)); 4940e1e9333SJames Wright } 495c5e9980aSAdeleke O. Bankole if (app_ctx->wall_forces.viewer) { 496c5e9980aSAdeleke O. Bankole PetscCall(TSMonitorSet(*ts, TSMonitor_WallForce, user, NULL)); 497c5e9980aSAdeleke O. Bankole } 498c931fa59SJames Wright if (app_ctx->turb_spanstats_enable) { 49991933550SJames Wright PetscCall(TSMonitorSet(*ts, TSMonitor_TurbulenceStatistics, user, NULL)); 500b8daee98SJames Wright CeedScalar previous_time = app_ctx->cont_time * user->units->second; 501b4c37c5cSJames Wright PetscCallCeed(user->ceed, 502b4c37c5cSJames Wright CeedOperatorSetContextDouble(user->spanstats.op_stats_collect_ctx->op, user->spanstats.previous_time_label, &previous_time)); 503b0488d1fSJames Wright } 50488b07121SJames Wright if (app_ctx->diff_filter_monitor) PetscCall(TSMonitorSet(*ts, TSMonitor_DifferentialFilter, user, NULL)); 505a515125bSLeila Ghaffari 506a515125bSLeila Ghaffari // Solve 50774a6f4ddSJed Brown PetscReal start_time; 50874a6f4ddSJed Brown PetscInt start_step; 5092b916ea7SJeremy L Thompson PetscCall(TSGetTime(*ts, &start_time)); 51074a6f4ddSJed Brown PetscCall(TSGetStepNumber(*ts, &start_step)); 51191982731SJeremy L Thompson 512df4304b5SJed Brown PetscCall(PetscLogDefaultBegin()); // So we can use PetscLogStageGetPerfInfo without -log_view 51391982731SJeremy L Thompson PetscPreLoadBegin(PETSC_FALSE, "Fluids Solve"); 51491982731SJeremy L Thompson PetscCall(TSSetTime(*ts, start_time)); 51574a6f4ddSJed Brown PetscCall(TSSetStepNumber(*ts, start_step)); 51691982731SJeremy L Thompson if (PetscPreLoadingOn) { 51791982731SJeremy L Thompson // LCOV_EXCL_START 51891982731SJeremy L Thompson SNES snes; 51991982731SJeremy L Thompson Vec Q_preload; 52091982731SJeremy L Thompson PetscReal rtol; 52191982731SJeremy L Thompson PetscCall(VecDuplicate(*Q, &Q_preload)); 52291982731SJeremy L Thompson PetscCall(VecCopy(*Q, Q_preload)); 52391982731SJeremy L Thompson PetscCall(TSGetSNES(*ts, &snes)); 52491982731SJeremy L Thompson PetscCall(SNESGetTolerances(snes, NULL, &rtol, NULL, NULL, NULL)); 5252b916ea7SJeremy L Thompson PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, .99, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 52622c1b34eSJames Wright PetscCall(TSSetSolution(*ts, Q_preload)); 52791982731SJeremy L Thompson PetscCall(TSStep(*ts)); 5282b916ea7SJeremy L Thompson PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, rtol, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 52991982731SJeremy L Thompson PetscCall(VecDestroy(&Q_preload)); 53091982731SJeremy L Thompson // LCOV_EXCL_STOP 53191982731SJeremy L Thompson } else { 5322b916ea7SJeremy L Thompson PetscCall(PetscBarrier((PetscObject)*ts)); 5332b916ea7SJeremy L Thompson PetscCall(TSSolve(*ts, *Q)); 53491982731SJeremy L Thompson } 53591982731SJeremy L Thompson PetscPreLoadEnd(); 53691982731SJeremy L Thompson 53791982731SJeremy L Thompson PetscCall(TSGetSolveTime(*ts, &final_time)); 538a515125bSLeila Ghaffari *f_time = final_time; 53991982731SJeremy L Thompson 5400e1e9333SJames Wright if (app_ctx->test_type == TESTTYPE_NONE) { 5417538d537SJames Wright PetscInt step_no; 5427538d537SJames Wright PetscCall(TSGetStepNumber(*ts, &step_no)); 543b0488d1fSJames Wright if (user->app_ctx->checkpoint_interval > 0 || user->app_ctx->checkpoint_interval == -1) { 5447538d537SJames Wright PetscCall(WriteOutput(user, *Q, step_no, final_time)); 5457538d537SJames Wright } 5467538d537SJames Wright 5477eedc94cSJames Wright PetscLogStage stage_id; 548df4304b5SJed Brown PetscEventPerfInfo stage_perf; 54991982731SJeremy L Thompson 55091982731SJeremy L Thompson PetscCall(PetscLogStageGetId("Fluids Solve", &stage_id)); 551df4304b5SJed Brown PetscCall(PetscLogStageGetPerfInfo(stage_id, &stage_perf)); 552df4304b5SJed Brown PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time taken for solution (sec): %g\n", stage_perf.time)); 553a515125bSLeila Ghaffari } 554d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 555a515125bSLeila Ghaffari } 556