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 Vec M_loc; 21a515125bSLeila Ghaffari CeedQFunction qf_mass; 22a515125bSLeila Ghaffari CeedOperator op_mass; 23a515125bSLeila Ghaffari CeedVector m_ceed, ones_vec; 24a515125bSLeila Ghaffari CeedInt num_comp_q, q_data_size; 25a515125bSLeila Ghaffari PetscFunctionBeginUser; 26a515125bSLeila Ghaffari 27a515125bSLeila Ghaffari // CEED Restriction 28a515125bSLeila Ghaffari CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q); 29a515125bSLeila Ghaffari CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &q_data_size); 30a515125bSLeila Ghaffari CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &m_ceed, NULL); 31a515125bSLeila Ghaffari CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &ones_vec, NULL); 32a515125bSLeila Ghaffari CeedVectorSetValue(ones_vec, 1.0); 33a515125bSLeila Ghaffari 34a515125bSLeila Ghaffari // CEED QFunction 359f59f36eSJames Wright PetscCall(CreateMassQFunction(ceed, num_comp_q, q_data_size, &qf_mass)); 36a515125bSLeila Ghaffari 37a515125bSLeila Ghaffari // CEED Operator 38a515125bSLeila Ghaffari CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass); 399f59f36eSJames Wright CeedOperatorSetField(op_mass, "u", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 402b916ea7SJeremy L Thompson CeedOperatorSetField(op_mass, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data); 412b916ea7SJeremy L Thompson CeedOperatorSetField(op_mass, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 42a515125bSLeila Ghaffari 43a515125bSLeila Ghaffari // Place PETSc vector in CEED vector 44a515125bSLeila Ghaffari PetscMemType m_mem_type; 452b916ea7SJeremy L Thompson PetscCall(DMGetLocalVector(dm, &M_loc)); 46fd969b44SJames Wright PetscCall(VecP2C(M_loc, &m_mem_type, m_ceed)); 47a515125bSLeila Ghaffari 48a515125bSLeila Ghaffari // Apply CEED Operator 49a515125bSLeila Ghaffari CeedOperatorApply(op_mass, ones_vec, m_ceed, CEED_REQUEST_IMMEDIATE); 50a515125bSLeila Ghaffari 51a515125bSLeila Ghaffari // Restore vectors 52fd969b44SJames Wright PetscCall(VecC2P(m_ceed, m_mem_type, M_loc)); 53a515125bSLeila Ghaffari 54a515125bSLeila Ghaffari // Local-to-Global 552b916ea7SJeremy L Thompson PetscCall(VecZeroEntries(M)); 562b916ea7SJeremy L Thompson PetscCall(DMLocalToGlobal(dm, M_loc, ADD_VALUES, M)); 572b916ea7SJeremy L Thompson PetscCall(DMRestoreLocalVector(dm, &M_loc)); 58a515125bSLeila Ghaffari 59a515125bSLeila Ghaffari // Invert diagonally lumped mass vector for RHS function 602b916ea7SJeremy L Thompson PetscCall(VecReciprocal(M)); 61a515125bSLeila Ghaffari 62a515125bSLeila Ghaffari // Cleanup 63a515125bSLeila Ghaffari CeedVectorDestroy(&ones_vec); 64a515125bSLeila Ghaffari CeedVectorDestroy(&m_ceed); 65a515125bSLeila Ghaffari CeedQFunctionDestroy(&qf_mass); 66a515125bSLeila Ghaffari CeedOperatorDestroy(&op_mass); 67a515125bSLeila Ghaffari 68a515125bSLeila Ghaffari PetscFunctionReturn(0); 69a515125bSLeila Ghaffari } 70a515125bSLeila Ghaffari 71c996854bSJames Wright // Insert Boundary values if it's a new time 72c996854bSJames Wright PetscErrorCode UpdateBoundaryValues(User user, Vec Q_loc, PetscReal t) { 73c996854bSJames Wright PetscFunctionBeginUser; 74c996854bSJames Wright if (user->time_bc_set != t) { 75c996854bSJames Wright PetscCall(DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Q_loc, t, NULL, NULL, NULL)); 76c996854bSJames Wright user->time_bc_set = t; 77c996854bSJames Wright } 78c996854bSJames Wright PetscFunctionReturn(0); 79c996854bSJames Wright } 80c996854bSJames Wright 8191f639d2SJames Wright // @brief Update the context label value to new value if necessary. 8291f639d2SJames Wright // @note This only supports labels with scalar label values (ie. not arrays) 8391f639d2SJames Wright PetscErrorCode UpdateContextLabel(MPI_Comm comm, PetscScalar update_value, CeedOperator op, CeedContextFieldLabel label) { 8491f639d2SJames Wright PetscScalar label_value; 85c2cb7fc8SJames Wright 8691f639d2SJames Wright PetscFunctionBeginUser; 8791f639d2SJames Wright PetscCheck(label, comm, PETSC_ERR_ARG_BADPTR, "Label should be non-NULL"); 8891f639d2SJames Wright 8991f639d2SJames Wright { 9091f639d2SJames Wright size_t num_elements; 9191f639d2SJames Wright const PetscScalar *label_values; 9291f639d2SJames Wright CeedOperatorGetContextDoubleRead(op, label, &num_elements, &label_values); 9391f639d2SJames Wright PetscCheck(num_elements == 1, comm, PETSC_ERR_SUP, "%s does not support labels with more than 1 value. Label has %zu values", __func__, 9491f639d2SJames Wright num_elements); 9591f639d2SJames Wright label_value = *label_values; 9691f639d2SJames Wright CeedOperatorRestoreContextDoubleRead(op, label, &label_values); 97c2cb7fc8SJames Wright } 9891f639d2SJames Wright 9991f639d2SJames Wright if (label_value != update_value) { 10091f639d2SJames Wright CeedOperatorSetContextDouble(op, label, &update_value); 101c2cb7fc8SJames Wright } 102c2cb7fc8SJames Wright PetscFunctionReturn(0); 103c2cb7fc8SJames Wright } 104c2cb7fc8SJames Wright 105a515125bSLeila Ghaffari // RHS (Explicit time-stepper) function setup 106a515125bSLeila Ghaffari // This is the RHS of the ODE, given as u_t = G(t,u) 107a515125bSLeila Ghaffari // This function takes in a state vector Q and writes into G 108a515125bSLeila Ghaffari PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data) { 109a515125bSLeila Ghaffari User user = *(User *)user_data; 11091f639d2SJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)ts); 111fd969b44SJames Wright PetscScalar dt; 112e2f84137SJeremy L Thompson Vec Q_loc = user->Q_loc, G_loc; 113a515125bSLeila Ghaffari PetscMemType q_mem_type, g_mem_type; 114a515125bSLeila Ghaffari PetscFunctionBeginUser; 115a515125bSLeila Ghaffari 116e2f84137SJeremy L Thompson // Get local vector 1172b916ea7SJeremy L Thompson PetscCall(DMGetLocalVector(user->dm, &G_loc)); 118e2f84137SJeremy L Thompson 119e2f84137SJeremy L Thompson // Update time dependent data 120c996854bSJames Wright PetscCall(UpdateBoundaryValues(user, Q_loc, t)); 12191f639d2SJames Wright if (user->phys->solution_time_label) PetscCall(UpdateContextLabel(comm, t, user->op_rhs, user->phys->solution_time_label)); 1222b916ea7SJeremy L Thompson PetscCall(TSGetTimeStep(ts, &dt)); 12391f639d2SJames Wright if (user->phys->timestep_size_label) PetscCall(UpdateContextLabel(comm, dt, user->op_rhs, user->phys->timestep_size_label)); 124a515125bSLeila Ghaffari 125a515125bSLeila Ghaffari // Global-to-local 1262b916ea7SJeremy L Thompson PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc)); 127a515125bSLeila Ghaffari 128a515125bSLeila Ghaffari // Place PETSc vectors in CEED vectors 129fd969b44SJames Wright PetscCall(VecReadP2C(Q_loc, &q_mem_type, user->q_ceed)); 130fd969b44SJames Wright PetscCall(VecP2C(G_loc, &g_mem_type, user->g_ceed)); 131a515125bSLeila Ghaffari 132a515125bSLeila Ghaffari // Apply CEED operator 1332b916ea7SJeremy L Thompson CeedOperatorApply(user->op_rhs, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE); 134a515125bSLeila Ghaffari 135a515125bSLeila Ghaffari // Restore vectors 136fd969b44SJames Wright PetscCall(VecReadC2P(user->q_ceed, q_mem_type, Q_loc)); 137fd969b44SJames Wright PetscCall(VecC2P(user->g_ceed, g_mem_type, G_loc)); 138a515125bSLeila Ghaffari 139a515125bSLeila Ghaffari // Local-to-Global 1402b916ea7SJeremy L Thompson PetscCall(VecZeroEntries(G)); 1412b916ea7SJeremy L Thompson PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G)); 142a515125bSLeila Ghaffari 143a515125bSLeila Ghaffari // Inverse of the lumped mass matrix (M is Minv) 1442b916ea7SJeremy L Thompson PetscCall(VecPointwiseMult(G, G, user->M)); 145a515125bSLeila Ghaffari 146a515125bSLeila Ghaffari // Restore vectors 1472b916ea7SJeremy L Thompson PetscCall(DMRestoreLocalVector(user->dm, &G_loc)); 148a515125bSLeila Ghaffari 149a515125bSLeila Ghaffari PetscFunctionReturn(0); 150a515125bSLeila Ghaffari } 151a515125bSLeila Ghaffari 152c5e9980aSAdeleke O. Bankole // Surface forces function setup 153c5e9980aSAdeleke O. Bankole static PetscErrorCode Surface_Forces_NS(DM dm, Vec G_loc, PetscInt num_walls, const PetscInt walls[], PetscScalar *reaction_force) { 154c5e9980aSAdeleke O. Bankole DMLabel face_label; 155c5e9980aSAdeleke O. Bankole const PetscScalar *g; 1562004e3acSAdeleke O. Bankole PetscInt dof, dim = 3; 157c5e9980aSAdeleke O. Bankole MPI_Comm comm; 1582004e3acSAdeleke O. Bankole PetscSection s; 159c5e9980aSAdeleke O. Bankole 160c5e9980aSAdeleke O. Bankole PetscFunctionBeginUser; 161c5e9980aSAdeleke O. Bankole PetscCall(PetscArrayzero(reaction_force, num_walls * dim)); 162c5e9980aSAdeleke O. Bankole PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 163c5e9980aSAdeleke O. Bankole PetscCall(DMGetLabel(dm, "Face Sets", &face_label)); 164c5e9980aSAdeleke O. Bankole PetscCall(VecGetArrayRead(G_loc, &g)); 165c5e9980aSAdeleke O. Bankole for (PetscInt w = 0; w < num_walls; w++) { 166c5e9980aSAdeleke O. Bankole const PetscInt wall = walls[w]; 167c5e9980aSAdeleke O. Bankole IS wall_is; 1682004e3acSAdeleke O. Bankole PetscCall(DMGetLocalSection(dm, &s)); 169c5e9980aSAdeleke O. Bankole PetscCall(DMLabelGetStratumIS(face_label, wall, &wall_is)); 170c5e9980aSAdeleke O. Bankole if (wall_is) { // There exist such points on this process 171c5e9980aSAdeleke O. Bankole PetscInt num_points; 1722004e3acSAdeleke O. Bankole PetscInt num_comp = 0; 173c5e9980aSAdeleke O. Bankole const PetscInt *points; 1742004e3acSAdeleke O. Bankole PetscCall(PetscSectionGetFieldComponents(s, 0, &num_comp)); 175c5e9980aSAdeleke O. Bankole PetscCall(ISGetSize(wall_is, &num_points)); 176c5e9980aSAdeleke O. Bankole PetscCall(ISGetIndices(wall_is, &points)); 177c5e9980aSAdeleke O. Bankole for (PetscInt i = 0; i < num_points; i++) { 178c5e9980aSAdeleke O. Bankole const PetscInt p = points[i]; 179c5e9980aSAdeleke O. Bankole const StateConservative *r; 180c5e9980aSAdeleke O. Bankole PetscCall(DMPlexPointLocalRead(dm, p, g, &r)); 1812004e3acSAdeleke O. Bankole PetscCall(PetscSectionGetDof(s, p, &dof)); 1822004e3acSAdeleke O. Bankole for (PetscInt node = 0; node < dof / num_comp; node++) { 183c5e9980aSAdeleke O. Bankole for (PetscInt j = 0; j < 3; j++) { 1842004e3acSAdeleke O. Bankole reaction_force[w * dim + j] -= r[node].momentum[j]; 1852004e3acSAdeleke O. Bankole } 186c5e9980aSAdeleke O. Bankole } 187c5e9980aSAdeleke O. Bankole } 188c5e9980aSAdeleke O. Bankole PetscCall(ISRestoreIndices(wall_is, &points)); 189c5e9980aSAdeleke O. Bankole } 190c5e9980aSAdeleke O. Bankole PetscCall(ISDestroy(&wall_is)); 191c5e9980aSAdeleke O. Bankole } 192c5e9980aSAdeleke O. Bankole PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, reaction_force, dim * num_walls, MPIU_SCALAR, MPI_SUM, comm)); 193c5e9980aSAdeleke O. Bankole // Restore Vectors 194c5e9980aSAdeleke O. Bankole PetscCall(VecRestoreArrayRead(G_loc, &g)); 195c5e9980aSAdeleke O. Bankole 196c5e9980aSAdeleke O. Bankole PetscFunctionReturn(0); 197c5e9980aSAdeleke O. Bankole } 198c5e9980aSAdeleke O. Bankole 199a515125bSLeila Ghaffari // Implicit time-stepper function setup 2002b916ea7SJeremy L Thompson PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data) { 201a515125bSLeila Ghaffari User user = *(User *)user_data; 20291f639d2SJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)ts); 203fd969b44SJames Wright PetscScalar dt; 204e2f84137SJeremy L Thompson Vec Q_loc = user->Q_loc, Q_dot_loc = user->Q_dot_loc, G_loc; 205a515125bSLeila Ghaffari PetscMemType q_mem_type, q_dot_mem_type, g_mem_type; 206a515125bSLeila Ghaffari PetscFunctionBeginUser; 207a515125bSLeila Ghaffari 208e2f84137SJeremy L Thompson // Get local vectors 209c5e9980aSAdeleke O. Bankole PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); 210e2f84137SJeremy L Thompson 211e2f84137SJeremy L Thompson // Update time dependent data 212c996854bSJames Wright PetscCall(UpdateBoundaryValues(user, Q_loc, t)); 21391f639d2SJames Wright if (user->phys->solution_time_label) PetscCall(UpdateContextLabel(comm, t, user->op_ifunction, user->phys->solution_time_label)); 2142b916ea7SJeremy L Thompson PetscCall(TSGetTimeStep(ts, &dt)); 21591f639d2SJames Wright if (user->phys->timestep_size_label) PetscCall(UpdateContextLabel(comm, dt, user->op_ifunction, user->phys->timestep_size_label)); 216a515125bSLeila Ghaffari 217a515125bSLeila Ghaffari // Global-to-local 21806108310SJames Wright PetscCall(DMGlobalToLocalBegin(user->dm, Q, INSERT_VALUES, Q_loc)); 21906108310SJames Wright PetscCall(DMGlobalToLocalBegin(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc)); 22006108310SJames Wright PetscCall(DMGlobalToLocalEnd(user->dm, Q, INSERT_VALUES, Q_loc)); 22106108310SJames Wright PetscCall(DMGlobalToLocalEnd(user->dm, Q_dot, INSERT_VALUES, Q_dot_loc)); 222a515125bSLeila Ghaffari 223a515125bSLeila Ghaffari // Place PETSc vectors in CEED vectors 224fd969b44SJames Wright PetscCall(VecReadP2C(Q_loc, &q_mem_type, user->q_ceed)); 225fd969b44SJames Wright PetscCall(VecReadP2C(Q_dot_loc, &q_dot_mem_type, user->q_dot_ceed)); 226fd969b44SJames Wright PetscCall(VecP2C(G_loc, &g_mem_type, user->g_ceed)); 227a515125bSLeila Ghaffari 228a515125bSLeila Ghaffari // Apply CEED operator 2292b916ea7SJeremy L Thompson CeedOperatorApply(user->op_ifunction, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE); 230a515125bSLeila Ghaffari 231a515125bSLeila Ghaffari // Restore vectors 232fd969b44SJames Wright PetscCall(VecReadC2P(user->q_ceed, q_mem_type, Q_loc)); 233fd969b44SJames Wright PetscCall(VecReadC2P(user->q_dot_ceed, q_dot_mem_type, Q_dot_loc)); 234fd969b44SJames Wright PetscCall(VecC2P(user->g_ceed, g_mem_type, G_loc)); 235a515125bSLeila Ghaffari 23601ab89c1SJames Wright if (user->app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) { 23701ab89c1SJames Wright PetscCall(SGS_DD_ModelApplyIFunction(user, Q_loc, G_loc)); 23801ab89c1SJames Wright } 2399c678832SJames Wright 240a515125bSLeila Ghaffari // Local-to-Global 2412b916ea7SJeremy L Thompson PetscCall(VecZeroEntries(G)); 2422b916ea7SJeremy L Thompson PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G)); 243a515125bSLeila Ghaffari 244a515125bSLeila Ghaffari // Restore vectors 245c5e9980aSAdeleke O. Bankole PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); 246a515125bSLeila Ghaffari 247a515125bSLeila Ghaffari PetscFunctionReturn(0); 248a515125bSLeila Ghaffari } 249a515125bSLeila Ghaffari 250f0b65372SJed Brown static PetscErrorCode MatMult_NS_IJacobian(Mat J, Vec Q, Vec G) { 251f0b65372SJed Brown User user; 252f0b65372SJed Brown const PetscScalar *q; 253f0b65372SJed Brown PetscScalar *g; 254f0b65372SJed Brown PetscMemType q_mem_type, g_mem_type; 255f0b65372SJed Brown PetscFunctionBeginUser; 2562b916ea7SJeremy L Thompson PetscCall(MatShellGetContext(J, &user)); 257e2f84137SJeremy L Thompson Vec Q_loc = user->Q_dot_loc, // Note - Q_dot_loc has zero BCs 258e2f84137SJeremy L Thompson G_loc; 259e2f84137SJeremy L Thompson 260f0b65372SJed Brown // Get local vectors 2612b916ea7SJeremy L Thompson PetscCall(DMGetLocalVector(user->dm, &G_loc)); 262f0b65372SJed Brown 263f0b65372SJed Brown // Global-to-local 2642b916ea7SJeremy L Thompson PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc)); 265f0b65372SJed Brown 266f0b65372SJed Brown // Place PETSc vectors in CEED vectors 2672b916ea7SJeremy L Thompson PetscCall(VecGetArrayReadAndMemType(Q_loc, &q, &q_mem_type)); 2682b916ea7SJeremy L Thompson PetscCall(VecGetArrayAndMemType(G_loc, &g, &g_mem_type)); 2692b916ea7SJeremy L Thompson CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER, (PetscScalar *)q); 270f0b65372SJed Brown CeedVectorSetArray(user->g_ceed, MemTypeP2C(g_mem_type), CEED_USE_POINTER, g); 271f0b65372SJed Brown 272f0b65372SJed Brown // Apply CEED operator 2732b916ea7SJeremy L Thompson CeedOperatorApply(user->op_ijacobian, user->q_ceed, user->g_ceed, CEED_REQUEST_IMMEDIATE); 274f0b65372SJed Brown 275f0b65372SJed Brown // Restore vectors 276f0b65372SJed Brown CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL); 277f0b65372SJed Brown CeedVectorTakeArray(user->g_ceed, MemTypeP2C(g_mem_type), NULL); 2782b916ea7SJeremy L Thompson PetscCall(VecRestoreArrayReadAndMemType(Q_loc, &q)); 2792b916ea7SJeremy L Thompson PetscCall(VecRestoreArrayAndMemType(G_loc, &g)); 280f0b65372SJed Brown 281f0b65372SJed Brown // Local-to-Global 2822b916ea7SJeremy L Thompson PetscCall(VecZeroEntries(G)); 2832b916ea7SJeremy L Thompson PetscCall(DMLocalToGlobal(user->dm, G_loc, ADD_VALUES, G)); 284f0b65372SJed Brown 285f0b65372SJed Brown // Restore vectors 2862b916ea7SJeremy L Thompson PetscCall(DMRestoreLocalVector(user->dm, &G_loc)); 287f0b65372SJed Brown PetscFunctionReturn(0); 288f0b65372SJed Brown } 289f0b65372SJed Brown 290f0b65372SJed Brown PetscErrorCode MatGetDiagonal_NS_IJacobian(Mat A, Vec D) { 291f0b65372SJed Brown User user; 292f0b65372SJed Brown Vec D_loc; 293f0b65372SJed Brown PetscScalar *d; 294f0b65372SJed Brown PetscMemType mem_type; 295f0b65372SJed Brown 296f0b65372SJed Brown PetscFunctionBeginUser; 297f0b65372SJed Brown MatShellGetContext(A, &user); 298f0b65372SJed Brown PetscCall(DMGetLocalVector(user->dm, &D_loc)); 299f0b65372SJed Brown PetscCall(VecGetArrayAndMemType(D_loc, &d, &mem_type)); 300f0b65372SJed Brown CeedVectorSetArray(user->g_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, d); 3012b916ea7SJeremy L Thompson CeedOperatorLinearAssembleDiagonal(user->op_ijacobian, user->g_ceed, CEED_REQUEST_IMMEDIATE); 302f0b65372SJed Brown CeedVectorTakeArray(user->g_ceed, MemTypeP2C(mem_type), NULL); 303f0b65372SJed Brown PetscCall(VecRestoreArrayAndMemType(D_loc, &d)); 304f0b65372SJed Brown PetscCall(VecZeroEntries(D)); 305f0b65372SJed Brown PetscCall(DMLocalToGlobal(user->dm, D_loc, ADD_VALUES, D)); 306f0b65372SJed Brown PetscCall(DMRestoreLocalVector(user->dm, &D_loc)); 307f0b65372SJed Brown VecViewFromOptions(D, NULL, "-diag_vec_view"); 308f0b65372SJed Brown PetscFunctionReturn(0); 309f0b65372SJed Brown } 310f0b65372SJed Brown 3112b916ea7SJeremy L Thompson static PetscErrorCode FormPreallocation(User user, PetscBool pbdiagonal, Mat J, CeedVector *coo_values) { 312b107fddaSJed Brown PetscCount ncoo; 313b107fddaSJed Brown PetscInt *rows, *cols; 314b107fddaSJed Brown 315b107fddaSJed Brown PetscFunctionBeginUser; 316b107fddaSJed Brown if (pbdiagonal) { 317b107fddaSJed Brown CeedSize l_size; 318b107fddaSJed Brown CeedOperatorGetActiveVectorLengths(user->op_ijacobian, &l_size, NULL); 319b107fddaSJed Brown ncoo = l_size * 5; 320b107fddaSJed Brown rows = malloc(ncoo * sizeof(rows[0])); 321b107fddaSJed Brown cols = malloc(ncoo * sizeof(cols[0])); 322b107fddaSJed Brown for (PetscCount n = 0; n < l_size / 5; n++) { 323b107fddaSJed Brown for (PetscInt i = 0; i < 5; i++) { 324b107fddaSJed Brown for (PetscInt j = 0; j < 5; j++) { 325b107fddaSJed Brown rows[(n * 5 + i) * 5 + j] = n * 5 + i; 326b107fddaSJed Brown cols[(n * 5 + i) * 5 + j] = n * 5 + j; 327b107fddaSJed Brown } 328b107fddaSJed Brown } 329b107fddaSJed Brown } 330b107fddaSJed Brown } else { 3312b916ea7SJeremy L Thompson PetscCall(CeedOperatorLinearAssembleSymbolic(user->op_ijacobian, &ncoo, &rows, &cols)); 332b107fddaSJed Brown } 333b107fddaSJed Brown PetscCall(MatSetPreallocationCOOLocal(J, ncoo, rows, cols)); 334b107fddaSJed Brown free(rows); 335b107fddaSJed Brown free(cols); 336b107fddaSJed Brown CeedVectorCreate(user->ceed, ncoo, coo_values); 337b107fddaSJed Brown PetscFunctionReturn(0); 338b107fddaSJed Brown } 339b107fddaSJed Brown 3402b916ea7SJeremy L Thompson static PetscErrorCode FormSetValues(User user, PetscBool pbdiagonal, Mat J, CeedVector coo_values) { 341b107fddaSJed Brown CeedMemType mem_type = CEED_MEM_HOST; 342b107fddaSJed Brown const PetscScalar *values; 343b107fddaSJed Brown MatType mat_type; 344b107fddaSJed Brown 345b107fddaSJed Brown PetscFunctionBeginUser; 346b107fddaSJed Brown PetscCall(MatGetType(J, &mat_type)); 3472b916ea7SJeremy L Thompson if (strstr(mat_type, "kokkos") || strstr(mat_type, "cusparse")) mem_type = CEED_MEM_DEVICE; 348b107fddaSJed Brown if (user->app_ctx->pmat_pbdiagonal) { 3492b916ea7SJeremy L Thompson CeedOperatorLinearAssemblePointBlockDiagonal(user->op_ijacobian, coo_values, CEED_REQUEST_IMMEDIATE); 350b107fddaSJed Brown } else { 351b107fddaSJed Brown CeedOperatorLinearAssemble(user->op_ijacobian, coo_values); 352b107fddaSJed Brown } 353b107fddaSJed Brown CeedVectorGetArrayRead(coo_values, mem_type, &values); 354b107fddaSJed Brown PetscCall(MatSetValuesCOO(J, values, INSERT_VALUES)); 355b107fddaSJed Brown CeedVectorRestoreArrayRead(coo_values, &values); 356b107fddaSJed Brown PetscFunctionReturn(0); 357b107fddaSJed Brown } 358b107fddaSJed Brown 3592b916ea7SJeremy L Thompson PetscErrorCode FormIJacobian_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, PetscReal shift, Mat J, Mat J_pre, void *user_data) { 360f0b65372SJed Brown User user = *(User *)user_data; 36104855949SJed Brown PetscBool J_is_shell, J_is_mffd, J_pre_is_shell; 362f0b65372SJed Brown PetscFunctionBeginUser; 363a5f46a7bSJeremy L Thompson if (user->phys->ijacobian_time_shift_label) CeedOperatorSetContextDouble(user->op_ijacobian, user->phys->ijacobian_time_shift_label, &shift); 36404855949SJed Brown PetscCall(PetscObjectTypeCompare((PetscObject)J, MATMFFD, &J_is_mffd)); 365f0b65372SJed Brown PetscCall(PetscObjectTypeCompare((PetscObject)J, MATSHELL, &J_is_shell)); 3662b916ea7SJeremy L Thompson PetscCall(PetscObjectTypeCompare((PetscObject)J_pre, MATSHELL, &J_pre_is_shell)); 367f0b65372SJed Brown if (!user->matrices_set_up) { 368f0b65372SJed Brown if (J_is_shell) { 369f0b65372SJed Brown PetscCall(MatShellSetContext(J, user)); 3702b916ea7SJeremy L Thompson PetscCall(MatShellSetOperation(J, MATOP_MULT, (void (*)(void))MatMult_NS_IJacobian)); 3712b916ea7SJeremy L Thompson PetscCall(MatShellSetOperation(J, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_NS_IJacobian)); 372f0b65372SJed Brown PetscCall(MatSetUp(J)); 373f0b65372SJed Brown } 374f0b65372SJed Brown if (!J_pre_is_shell) { 3752b916ea7SJeremy L Thompson PetscCall(FormPreallocation(user, user->app_ctx->pmat_pbdiagonal, J_pre, &user->coo_values_pmat)); 376b107fddaSJed Brown } 37704855949SJed Brown if (J != J_pre && !J_is_shell && !J_is_mffd) { 378b107fddaSJed Brown PetscCall(FormPreallocation(user, PETSC_FALSE, J, &user->coo_values_amat)); 379b107fddaSJed Brown } 380f0b65372SJed Brown user->matrices_set_up = true; 381f0b65372SJed Brown } 382f0b65372SJed Brown if (!J_pre_is_shell) { 3832b916ea7SJeremy L Thompson PetscCall(FormSetValues(user, user->app_ctx->pmat_pbdiagonal, J_pre, user->coo_values_pmat)); 384f0b65372SJed Brown } 38504855949SJed Brown if (user->coo_values_amat) { 38604855949SJed Brown PetscCall(FormSetValues(user, PETSC_FALSE, J, user->coo_values_amat)); 38704855949SJed Brown } else if (J_is_mffd) { 38804855949SJed Brown PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 38904855949SJed Brown PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 39004855949SJed Brown } 391f0b65372SJed Brown PetscFunctionReturn(0); 392f0b65372SJed Brown } 393f0b65372SJed Brown 3942b916ea7SJeremy L Thompson PetscErrorCode WriteOutput(User user, Vec Q, PetscInt step_no, PetscScalar time) { 395a515125bSLeila Ghaffari Vec Q_loc; 396a515125bSLeila Ghaffari char file_path[PETSC_MAX_PATH_LEN]; 397a515125bSLeila Ghaffari PetscViewer viewer; 398a515125bSLeila Ghaffari PetscFunctionBeginUser; 399a515125bSLeila Ghaffari 400852e5969SJed Brown if (user->app_ctx->checkpoint_vtk) { 401a515125bSLeila Ghaffari // Set up output 4027538d537SJames Wright PetscCall(DMGetLocalVector(user->dm, &Q_loc)); 4037538d537SJames Wright PetscCall(PetscObjectSetName((PetscObject)Q_loc, "StateVec")); 4047538d537SJames Wright PetscCall(VecZeroEntries(Q_loc)); 4057538d537SJames Wright PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc)); 406a515125bSLeila Ghaffari 407a515125bSLeila Ghaffari // Output 408852e5969SJed Brown PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no)); 4097538d537SJames Wright 4102b916ea7SJeremy L Thompson PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), file_path, FILE_MODE_WRITE, &viewer)); 4117538d537SJames Wright PetscCall(VecView(Q_loc, viewer)); 4127538d537SJames Wright PetscCall(PetscViewerDestroy(&viewer)); 413a515125bSLeila Ghaffari if (user->dm_viz) { 414a515125bSLeila Ghaffari Vec Q_refined, Q_refined_loc; 415a515125bSLeila Ghaffari char file_path_refined[PETSC_MAX_PATH_LEN]; 416a515125bSLeila Ghaffari PetscViewer viewer_refined; 417a515125bSLeila Ghaffari 4187538d537SJames Wright PetscCall(DMGetGlobalVector(user->dm_viz, &Q_refined)); 4197538d537SJames Wright PetscCall(DMGetLocalVector(user->dm_viz, &Q_refined_loc)); 4207538d537SJames Wright PetscCall(PetscObjectSetName((PetscObject)Q_refined_loc, "Refined")); 4217538d537SJames Wright 4227538d537SJames Wright PetscCall(MatInterpolate(user->interp_viz, Q, Q_refined)); 4237538d537SJames Wright PetscCall(VecZeroEntries(Q_refined_loc)); 4242b916ea7SJeremy L Thompson PetscCall(DMGlobalToLocal(user->dm_viz, Q_refined, INSERT_VALUES, Q_refined_loc)); 4257538d537SJames Wright 426852e5969SJed Brown PetscCall( 427852e5969SJed Brown PetscSNPrintf(file_path_refined, sizeof file_path_refined, "%s/nsrefined-%03" PetscInt_FMT ".vtu", user->app_ctx->output_dir, step_no)); 4287538d537SJames Wright 4292b916ea7SJeremy L Thompson PetscCall(PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q_refined), file_path_refined, FILE_MODE_WRITE, &viewer_refined)); 4307538d537SJames Wright PetscCall(VecView(Q_refined_loc, viewer_refined)); 4317538d537SJames Wright PetscCall(DMRestoreLocalVector(user->dm_viz, &Q_refined_loc)); 4327538d537SJames Wright PetscCall(DMRestoreGlobalVector(user->dm_viz, &Q_refined)); 4337538d537SJames Wright PetscCall(PetscViewerDestroy(&viewer_refined)); 434a515125bSLeila Ghaffari } 4357538d537SJames Wright PetscCall(DMRestoreLocalVector(user->dm, &Q_loc)); 436852e5969SJed Brown } 437a515125bSLeila Ghaffari 438a515125bSLeila Ghaffari // Save data in a binary file for continuation of simulations 43991a36801SJames Wright if (user->app_ctx->add_stepnum2bin) { 440852e5969SJed Brown PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution-%" PetscInt_FMT ".bin", user->app_ctx->output_dir, step_no)); 44191a36801SJames Wright } else { 4422b916ea7SJeremy L Thompson PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/ns-solution.bin", user->app_ctx->output_dir)); 44391a36801SJames Wright } 4442b916ea7SJeremy L Thompson PetscCall(PetscViewerBinaryOpen(user->comm, file_path, FILE_MODE_WRITE, &viewer)); 4457538d537SJames Wright 4469293eaa1SJed Brown PetscInt token = FLUIDS_FILE_TOKEN; 4479293eaa1SJed Brown PetscCall(PetscViewerBinaryWrite(viewer, &token, 1, PETSC_INT)); 4489293eaa1SJed Brown PetscCall(PetscViewerBinaryWrite(viewer, &step_no, 1, PETSC_INT)); 4499293eaa1SJed Brown time /= user->units->second; // Dimensionalize time back 4509293eaa1SJed Brown PetscCall(PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL)); 4517538d537SJames Wright PetscCall(VecView(Q, viewer)); 4527538d537SJames Wright PetscCall(PetscViewerDestroy(&viewer)); 4537538d537SJames Wright PetscFunctionReturn(0); 4547538d537SJames Wright } 4557538d537SJames Wright 456c5e9980aSAdeleke O. Bankole // CSV Monitor 457c5e9980aSAdeleke O. Bankole PetscErrorCode TSMonitor_WallForce(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) { 458c5e9980aSAdeleke O. Bankole User user = ctx; 459c5e9980aSAdeleke O. Bankole Vec G_loc; 460c5e9980aSAdeleke O. Bankole PetscInt num_wall = user->app_ctx->wall_forces.num_wall, dim = 3; 461c5e9980aSAdeleke O. Bankole const PetscInt *walls = user->app_ctx->wall_forces.walls; 462c5e9980aSAdeleke O. Bankole PetscViewer viewer = user->app_ctx->wall_forces.viewer; 463c5e9980aSAdeleke O. Bankole PetscViewerFormat format = user->app_ctx->wall_forces.viewer_format; 464c5e9980aSAdeleke O. Bankole PetscScalar *reaction_force; 465c5e9980aSAdeleke O. Bankole PetscBool iascii; 466c5e9980aSAdeleke O. Bankole 467c5e9980aSAdeleke O. Bankole PetscFunctionBeginUser; 468c5e9980aSAdeleke O. Bankole if (!viewer) PetscFunctionReturn(0); 469c5e9980aSAdeleke O. Bankole PetscCall(DMGetNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); 470c5e9980aSAdeleke O. Bankole PetscCall(PetscMalloc1(num_wall * dim, &reaction_force)); 471c5e9980aSAdeleke O. Bankole PetscCall(Surface_Forces_NS(user->dm, G_loc, num_wall, walls, reaction_force)); 472c5e9980aSAdeleke O. Bankole PetscCall(DMRestoreNamedLocalVector(user->dm, "ResidualLocal", &G_loc)); 473c5e9980aSAdeleke O. Bankole 474c5e9980aSAdeleke O. Bankole PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 475c5e9980aSAdeleke O. Bankole 476c5e9980aSAdeleke O. Bankole if (iascii) { 477c5e9980aSAdeleke O. Bankole if (format == PETSC_VIEWER_ASCII_CSV && !user->app_ctx->wall_forces.header_written) { 478c5e9980aSAdeleke O. Bankole PetscCall(PetscViewerASCIIPrintf(viewer, "Step,Time,Wall,ForceX,ForceY,ForceZ\n")); 479c5e9980aSAdeleke O. Bankole user->app_ctx->wall_forces.header_written = PETSC_TRUE; 480c5e9980aSAdeleke O. Bankole } 481c5e9980aSAdeleke O. Bankole for (PetscInt w = 0; w < num_wall; w++) { 482c5e9980aSAdeleke O. Bankole PetscInt wall = walls[w]; 483c5e9980aSAdeleke O. Bankole if (format == PETSC_VIEWER_ASCII_CSV) { 484c5e9980aSAdeleke O. Bankole PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT ",%g,%" PetscInt_FMT ",%g,%g,%g\n", step_no, time, wall, 485c5e9980aSAdeleke O. Bankole reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2])); 486c5e9980aSAdeleke O. Bankole 487c5e9980aSAdeleke O. Bankole } else { 488c5e9980aSAdeleke O. Bankole PetscCall(PetscViewerASCIIPrintf(viewer, "Wall %" PetscInt_FMT " Forces: Force_x = %12g, Force_y = %12g, Force_z = %12g\n", wall, 489c5e9980aSAdeleke O. Bankole reaction_force[w * dim + 0], reaction_force[w * dim + 1], reaction_force[w * dim + 2])); 490c5e9980aSAdeleke O. Bankole } 491c5e9980aSAdeleke O. Bankole } 492c5e9980aSAdeleke O. Bankole } 493c5e9980aSAdeleke O. Bankole PetscCall(PetscFree(reaction_force)); 494c5e9980aSAdeleke O. Bankole PetscFunctionReturn(0); 495c5e9980aSAdeleke O. Bankole } 496c5e9980aSAdeleke O. Bankole 4977538d537SJames Wright // User provided TS Monitor 4982b916ea7SJeremy L Thompson PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx) { 4997538d537SJames Wright User user = ctx; 5007538d537SJames Wright PetscFunctionBeginUser; 5017538d537SJames Wright 502852e5969SJed Brown // Print every 'checkpoint_interval' steps 503c539088bSJames Wright if (user->app_ctx->checkpoint_interval <= 0 || step_no % user->app_ctx->checkpoint_interval != 0 || 504e419654dSJeremy L Thompson (user->app_ctx->cont_steps == step_no && step_no != 0)) { 505c539088bSJames Wright PetscFunctionReturn(0); 506e419654dSJeremy L Thompson } 5077538d537SJames Wright 5087538d537SJames Wright PetscCall(WriteOutput(user, Q, step_no, time)); 509a515125bSLeila Ghaffari 510a515125bSLeila Ghaffari PetscFunctionReturn(0); 511a515125bSLeila Ghaffari } 512a515125bSLeila Ghaffari 513a515125bSLeila Ghaffari // TS: Create, setup, and solve 5142b916ea7SJeremy L Thompson PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, Vec *Q, PetscScalar *f_time, TS *ts) { 515a515125bSLeila Ghaffari MPI_Comm comm = user->comm; 516a515125bSLeila Ghaffari TSAdapt adapt; 517a515125bSLeila Ghaffari PetscScalar final_time; 518a515125bSLeila Ghaffari PetscFunctionBeginUser; 519a515125bSLeila Ghaffari 5202b916ea7SJeremy L Thompson PetscCall(TSCreate(comm, ts)); 5212b916ea7SJeremy L Thompson PetscCall(TSSetDM(*ts, dm)); 522a515125bSLeila Ghaffari if (phys->implicit) { 5232b916ea7SJeremy L Thompson PetscCall(TSSetType(*ts, TSBDF)); 524a515125bSLeila Ghaffari if (user->op_ifunction) { 5252b916ea7SJeremy L Thompson PetscCall(TSSetIFunction(*ts, NULL, IFunction_NS, &user)); 526a515125bSLeila Ghaffari } else { // Implicit integrators can fall back to using an RHSFunction 5272b916ea7SJeremy L Thompson PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user)); 528a515125bSLeila Ghaffari } 529f0b65372SJed Brown if (user->op_ijacobian) { 5302b916ea7SJeremy L Thompson PetscCall(DMTSSetIJacobian(dm, FormIJacobian_NS, &user)); 531b107fddaSJed Brown if (app_ctx->amat_type) { 532b107fddaSJed Brown Mat Pmat, Amat; 5332b916ea7SJeremy L Thompson PetscCall(DMCreateMatrix(dm, &Pmat)); 5342b916ea7SJeremy L Thompson PetscCall(DMSetMatType(dm, app_ctx->amat_type)); 5352b916ea7SJeremy L Thompson PetscCall(DMCreateMatrix(dm, &Amat)); 5362b916ea7SJeremy L Thompson PetscCall(TSSetIJacobian(*ts, Amat, Pmat, NULL, NULL)); 5372b916ea7SJeremy L Thompson PetscCall(MatDestroy(&Amat)); 5382b916ea7SJeremy L Thompson PetscCall(MatDestroy(&Pmat)); 539b107fddaSJed Brown } 540f0b65372SJed Brown } 541a515125bSLeila Ghaffari } else { 542*5d28dccaSJames Wright PetscCheck(user->op_rhs, comm, PETSC_ERR_ARG_NULL, "Problem does not provide RHSFunction"); 5432b916ea7SJeremy L Thompson PetscCall(TSSetType(*ts, TSRK)); 5442b916ea7SJeremy L Thompson PetscCall(TSRKSetType(*ts, TSRK5F)); 5452b916ea7SJeremy L Thompson PetscCall(TSSetRHSFunction(*ts, NULL, RHS_NS, &user)); 546a515125bSLeila Ghaffari } 5472b916ea7SJeremy L Thompson PetscCall(TSSetMaxTime(*ts, 500. * user->units->second)); 5482b916ea7SJeremy L Thompson PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER)); 54922387d3aSJames Wright if (app_ctx->test_type == TESTTYPE_NONE) PetscCall(TSSetErrorIfStepFails(*ts, PETSC_FALSE)); 5502b916ea7SJeremy L Thompson PetscCall(TSSetTimeStep(*ts, 1.e-2 * user->units->second)); 5510e1e9333SJames Wright if (app_ctx->test_type != TESTTYPE_NONE) { 5522b916ea7SJeremy L Thompson PetscCall(TSSetMaxSteps(*ts, 10)); 5532b916ea7SJeremy L Thompson } 5542b916ea7SJeremy L Thompson PetscCall(TSGetAdapt(*ts, &adapt)); 5552b916ea7SJeremy L Thompson PetscCall(TSAdaptSetStepLimits(adapt, 1.e-12 * user->units->second, 1.e2 * user->units->second)); 5562b916ea7SJeremy L Thompson PetscCall(TSSetFromOptions(*ts)); 55791f639d2SJames Wright user->time_bc_set = -1.0; // require all BCs be updated 558a515125bSLeila Ghaffari if (!app_ctx->cont_steps) { // print initial condition 5590e1e9333SJames Wright if (app_ctx->test_type == TESTTYPE_NONE) { 5602b916ea7SJeremy L Thompson PetscCall(TSMonitor_NS(*ts, 0, 0., *Q, user)); 561a515125bSLeila Ghaffari } 562a515125bSLeila Ghaffari } else { // continue from time of last output 563a515125bSLeila Ghaffari PetscInt count; 564a515125bSLeila Ghaffari PetscViewer viewer; 5652b916ea7SJeremy L Thompson 5669293eaa1SJed Brown if (app_ctx->cont_time <= 0) { // Legacy files did not include step number and time 5672b916ea7SJeremy L Thompson PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_time_file, FILE_MODE_READ, &viewer)); 5689293eaa1SJed Brown PetscCall(PetscViewerBinaryRead(viewer, &app_ctx->cont_time, 1, &count, PETSC_REAL)); 5692b916ea7SJeremy L Thompson PetscCall(PetscViewerDestroy(&viewer)); 5709293eaa1SJed Brown PetscCheck(app_ctx->cont_steps != -1, comm, PETSC_ERR_ARG_INCOMP, 5719293eaa1SJed Brown "-continue step number not specified, but checkpoint file does not contain a step number (likely written by older code version)"); 5729293eaa1SJed Brown } 5739293eaa1SJed Brown PetscCall(TSSetTime(*ts, app_ctx->cont_time * user->units->second)); 57474a6f4ddSJed Brown PetscCall(TSSetStepNumber(*ts, app_ctx->cont_steps)); 575a515125bSLeila Ghaffari } 5760e1e9333SJames Wright if (app_ctx->test_type == TESTTYPE_NONE) { 5772b916ea7SJeremy L Thompson PetscCall(TSMonitorSet(*ts, TSMonitor_NS, user, NULL)); 5780e1e9333SJames Wright } 579c5e9980aSAdeleke O. Bankole if (app_ctx->wall_forces.viewer) { 580c5e9980aSAdeleke O. Bankole PetscCall(TSMonitorSet(*ts, TSMonitor_WallForce, user, NULL)); 581c5e9980aSAdeleke O. Bankole } 582c931fa59SJames Wright if (app_ctx->turb_spanstats_enable) { 583b0488d1fSJames Wright PetscCall(TSMonitorSet(*ts, TSMonitor_Statistics, user, NULL)); 584b8daee98SJames Wright CeedScalar previous_time = app_ctx->cont_time * user->units->second; 585a5f46a7bSJeremy L Thompson CeedOperatorSetContextDouble(user->spanstats.op_stats_collect, user->spanstats.previous_time_label, &previous_time); 586b0488d1fSJames Wright } 587a515125bSLeila Ghaffari 588a515125bSLeila Ghaffari // Solve 58974a6f4ddSJed Brown PetscReal start_time; 59074a6f4ddSJed Brown PetscInt start_step; 5912b916ea7SJeremy L Thompson PetscCall(TSGetTime(*ts, &start_time)); 59274a6f4ddSJed Brown PetscCall(TSGetStepNumber(*ts, &start_step)); 59391982731SJeremy L Thompson 59491982731SJeremy L Thompson PetscPreLoadBegin(PETSC_FALSE, "Fluids Solve"); 59591982731SJeremy L Thompson PetscCall(TSSetTime(*ts, start_time)); 59674a6f4ddSJed Brown PetscCall(TSSetStepNumber(*ts, start_step)); 59791982731SJeremy L Thompson if (PetscPreLoadingOn) { 59891982731SJeremy L Thompson // LCOV_EXCL_START 59991982731SJeremy L Thompson SNES snes; 60091982731SJeremy L Thompson Vec Q_preload; 60191982731SJeremy L Thompson PetscReal rtol; 60291982731SJeremy L Thompson PetscCall(VecDuplicate(*Q, &Q_preload)); 60391982731SJeremy L Thompson PetscCall(VecCopy(*Q, Q_preload)); 60491982731SJeremy L Thompson PetscCall(TSGetSNES(*ts, &snes)); 60591982731SJeremy L Thompson PetscCall(SNESGetTolerances(snes, NULL, &rtol, NULL, NULL, NULL)); 6062b916ea7SJeremy L Thompson PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, .99, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 60791982731SJeremy L Thompson PetscCall(TSSetSolution(*ts, *Q)); 60891982731SJeremy L Thompson PetscCall(TSStep(*ts)); 6092b916ea7SJeremy L Thompson PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, rtol, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 61091982731SJeremy L Thompson PetscCall(VecDestroy(&Q_preload)); 61191982731SJeremy L Thompson // LCOV_EXCL_STOP 61291982731SJeremy L Thompson } else { 6132b916ea7SJeremy L Thompson PetscCall(PetscBarrier((PetscObject)*ts)); 6142b916ea7SJeremy L Thompson PetscCall(TSSolve(*ts, *Q)); 61591982731SJeremy L Thompson } 61691982731SJeremy L Thompson PetscPreLoadEnd(); 61791982731SJeremy L Thompson 61891982731SJeremy L Thompson PetscCall(TSGetSolveTime(*ts, &final_time)); 619a515125bSLeila Ghaffari *f_time = final_time; 62091982731SJeremy L Thompson 6210e1e9333SJames Wright if (app_ctx->test_type == TESTTYPE_NONE) { 6227538d537SJames Wright PetscInt step_no; 6237538d537SJames Wright PetscCall(TSGetStepNumber(*ts, &step_no)); 624b0488d1fSJames Wright if (user->app_ctx->checkpoint_interval > 0 || user->app_ctx->checkpoint_interval == -1) { 6257538d537SJames Wright PetscCall(WriteOutput(user, *Q, step_no, final_time)); 6267538d537SJames Wright } 6277538d537SJames Wright 62891982731SJeremy L Thompson PetscLogEvent stage_id; 62991982731SJeremy L Thompson PetscStageLog stage_log; 63091982731SJeremy L Thompson 63191982731SJeremy L Thompson PetscCall(PetscLogStageGetId("Fluids Solve", &stage_id)); 63291982731SJeremy L Thompson PetscCall(PetscLogGetStageLog(&stage_log)); 6332b916ea7SJeremy L Thompson PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time taken for solution (sec): %g\n", stage_log->stageInfo[stage_id].perfInfo.time)); 634a515125bSLeila Ghaffari } 635a515125bSLeila Ghaffari PetscFunctionReturn(0); 636a515125bSLeila Ghaffari } 637