15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, 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 /// Setup libCEED for Navier-Stokes example using PETSc 1077841947SLeila Ghaffari 1149aac155SJeremy L Thompson #include <ceed.h> 1249aac155SJeremy L Thompson #include <petscdmplex.h> 1349aac155SJeremy L Thompson 1477841947SLeila Ghaffari #include "../navierstokes.h" 1577841947SLeila Ghaffari 16b4e9a8f8SJames Wright // @brief Create CeedOperator for unstabilized mass KSP for explicit timestepping 17b4e9a8f8SJames Wright static PetscErrorCode CreateKSPMassOperator_Unstabilized(User user, CeedOperator *op_mass) { 18b4e9a8f8SJames Wright Ceed ceed = user->ceed; 19b4e9a8f8SJames Wright CeedInt num_comp_q, q_data_size; 20b4e9a8f8SJames Wright CeedQFunction qf_mass; 21b4e9a8f8SJames Wright CeedElemRestriction elem_restr_q, elem_restr_qd_i; 22b4e9a8f8SJames Wright CeedBasis basis_q; 23b4e9a8f8SJames Wright CeedVector q_data; 24b4e9a8f8SJames Wright 25b4e9a8f8SJames Wright PetscFunctionBeginUser; 26b4e9a8f8SJames Wright { // Get restriction and basis from the RHS function 27b4e9a8f8SJames Wright CeedOperator *sub_ops; 28b4e9a8f8SJames Wright CeedOperatorField field; 29b4e9a8f8SJames Wright PetscInt sub_op_index = 0; // will be 0 for the volume op 30b4e9a8f8SJames Wright 31b4e9a8f8SJames Wright PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(user->op_rhs_ctx->op, &sub_ops)); 32b4e9a8f8SJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &field)); 33b4e9a8f8SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_q)); 34b4e9a8f8SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetBasis(field, &basis_q)); 35b4e9a8f8SJames Wright 36b4e9a8f8SJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &field)); 37b4e9a8f8SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_qd_i)); 38b4e9a8f8SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetVector(field, &q_data)); 39b4e9a8f8SJames Wright } 40b4e9a8f8SJames Wright 41b4e9a8f8SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q)); 42b4e9a8f8SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd_i, &q_data_size)); 43b4e9a8f8SJames Wright 44b4e9a8f8SJames Wright PetscCall(CreateMassQFunction(ceed, num_comp_q, q_data_size, &qf_mass)); 45b4e9a8f8SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass)); 46b4e9a8f8SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "u", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 47b4e9a8f8SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data)); 48b4e9a8f8SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 49b4e9a8f8SJames Wright 50b4e9a8f8SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 51b4e9a8f8SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 52b4e9a8f8SJames Wright } 53b4e9a8f8SJames Wright 54b4e9a8f8SJames Wright // @brief Create KSP to solve the inverse mass operator for explicit time stepping schemes 55731c13d7SJames Wright static PetscErrorCode CreateKSPMass(User user, ProblemData problem) { 56b4e9a8f8SJames Wright Ceed ceed = user->ceed; 57b4e9a8f8SJames Wright DM dm = user->dm; 58b4e9a8f8SJames Wright CeedOperator op_mass; 59b4e9a8f8SJames Wright 60b4e9a8f8SJames Wright PetscFunctionBeginUser; 61b4e9a8f8SJames Wright if (problem->create_mass_operator) PetscCall(problem->create_mass_operator(user, &op_mass)); 62b4e9a8f8SJames Wright else PetscCall(CreateKSPMassOperator_Unstabilized(user, &op_mass)); 63b4e9a8f8SJames Wright 64b4e9a8f8SJames Wright { // -- Setup KSP for mass operator 65b4e9a8f8SJames Wright Mat mat_mass; 66b4e9a8f8SJames Wright Vec Zeros_loc; 67b4e9a8f8SJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)dm); 68b4e9a8f8SJames Wright 69b4e9a8f8SJames Wright PetscCall(DMCreateLocalVector(dm, &Zeros_loc)); 70b4e9a8f8SJames Wright PetscCall(VecZeroEntries(Zeros_loc)); 71b4e9a8f8SJames Wright PetscCall(MatCeedCreate(dm, dm, op_mass, NULL, &mat_mass)); 72b4e9a8f8SJames Wright PetscCall(MatCeedSetLocalVectors(mat_mass, Zeros_loc, NULL)); 73b4e9a8f8SJames Wright 74b4e9a8f8SJames Wright PetscCall(KSPCreate(comm, &user->mass_ksp)); 75b4e9a8f8SJames Wright PetscCall(KSPSetOptionsPrefix(user->mass_ksp, "mass_")); 76b4e9a8f8SJames Wright { // lumped by default 77b4e9a8f8SJames Wright PC pc; 78b4e9a8f8SJames Wright PetscCall(KSPGetPC(user->mass_ksp, &pc)); 79b4e9a8f8SJames Wright PetscCall(PCSetType(pc, PCJACOBI)); 80b4e9a8f8SJames Wright PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM)); 81b4e9a8f8SJames Wright PetscCall(KSPSetType(user->mass_ksp, KSPPREONLY)); 82b4e9a8f8SJames Wright } 83b4e9a8f8SJames Wright PetscCall(KSPSetFromOptions_WithMatCeed(user->mass_ksp, mat_mass)); 84b4e9a8f8SJames Wright PetscCall(KSPSetFromOptions(user->mass_ksp)); 85b4e9a8f8SJames Wright PetscCall(VecDestroy(&Zeros_loc)); 86b4e9a8f8SJames Wright PetscCall(MatDestroy(&mat_mass)); 87b4e9a8f8SJames Wright } 88b4e9a8f8SJames Wright 89b4e9a8f8SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass)); 90b4e9a8f8SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 91b4e9a8f8SJames Wright } 92b4e9a8f8SJames Wright 93e9c36be0SJames Wright static PetscErrorCode AddBCSubOperator(Ceed ceed, DM dm, CeedData ceed_data, DMLabel domain_label, PetscInt label_value, CeedInt height, 94e9c36be0SJames Wright CeedInt Q_sur, CeedInt q_data_size_sur, CeedInt jac_data_size_sur, CeedBasis basis_q_sur, 95e9c36be0SJames Wright CeedBasis basis_x_sur, CeedQFunction qf_apply_bc, CeedQFunction qf_apply_bc_jacobian, CeedOperator op_apply, 96e9c36be0SJames Wright CeedOperator op_apply_ijacobian) { 97bb85d312SJames Wright CeedVector q_data_sur, jac_data_sur = NULL; 98e9c36be0SJames Wright CeedOperator op_apply_bc, op_apply_bc_jacobian = NULL; 99bb85d312SJames Wright CeedElemRestriction elem_restr_x_sur, elem_restr_q_sur, elem_restr_qd_i_sur, elem_restr_jd_i_sur = NULL; 100e9c36be0SJames Wright PetscInt dm_field = 0; 1019eb9c072SJames Wright 102f17d818dSJames Wright PetscFunctionBeginUser; 103bb85d312SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, domain_label, label_value, height, dm_field, &elem_restr_q_sur)); 104bb85d312SJames Wright PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height, &elem_restr_x_sur)); 1059eb9c072SJames Wright if (jac_data_size_sur > 0) { 1069eb9c072SJames Wright // State-dependent data will be passed from residual to Jacobian. This will be collocated. 107bb85d312SJames Wright PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, jac_data_size_sur, &elem_restr_jd_i_sur)); 108a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_jd_i_sur, &jac_data_sur, NULL)); 1099eb9c072SJames Wright } 1109eb9c072SJames Wright 111*cfb075a4SJames Wright PetscCall(QDataBoundaryGet(ceed, dm, domain_label, label_value, elem_restr_x_sur, basis_x_sur, ceed_data->x_coord, &elem_restr_qd_i_sur, 112*cfb075a4SJames Wright &q_data_sur, &q_data_size_sur)); 113e9c36be0SJames Wright 114a5b0ec6fSJames Wright // CEED Operator for Physics 115a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_apply_bc, NULL, NULL, &op_apply_bc)); 116e9c36be0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "q", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE)); 117e9c36be0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "Grad_q", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE)); 118356036faSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_NONE, q_data_sur)); 119e9c36be0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "x", elem_restr_x_sur, basis_x_sur, ceed_data->x_coord)); 120e9c36be0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "v", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE)); 121a424bcd0SJames Wright if (elem_restr_jd_i_sur) 122356036faSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_NONE, jac_data_sur)); 1239eb9c072SJames Wright 124f21e6b1cSJames Wright if (qf_apply_bc_jacobian && elem_restr_jd_i_sur) { 125a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_apply_bc_jacobian, NULL, NULL, &op_apply_bc_jacobian)); 126e9c36be0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "dq", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE)); 127e9c36be0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "Grad_dq", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE)); 128356036faSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_NONE, q_data_sur)); 129e9c36be0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "x", elem_restr_x_sur, basis_x_sur, ceed_data->x_coord)); 130356036faSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_NONE, jac_data_sur)); 131e9c36be0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "v", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE)); 1329eb9c072SJames Wright } 1339eb9c072SJames Wright 134a5b0ec6fSJames Wright // Apply Sub-Operator for Physics 135e9c36be0SJames Wright PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_apply, op_apply_bc)); 136e9c36be0SJames Wright if (op_apply_bc_jacobian) PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_apply_ijacobian, op_apply_bc_jacobian)); 1379eb9c072SJames Wright 138a424bcd0SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data_sur)); 139a424bcd0SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&jac_data_sur)); 140a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q_sur)); 141a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x_sur)); 142a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd_i_sur)); 143a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_jd_i_sur)); 144a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_apply_bc)); 145a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_apply_bc_jacobian)); 146ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 1479eb9c072SJames Wright } 1489eb9c072SJames Wright 149e9c36be0SJames Wright static PetscErrorCode SetupBCQFunctions(Ceed ceed, PetscInt dim_sur, PetscInt num_comp_x, PetscInt num_comp_q, PetscInt q_data_size_sur, 1502b730f8bSJeremy L Thompson PetscInt jac_data_size_sur, ProblemQFunctionSpec apply_bc, ProblemQFunctionSpec apply_bc_jacobian, 1515b881d11SJames Wright CeedQFunction *qf_apply_bc, CeedQFunction *qf_apply_bc_jacobian) { 1525b881d11SJames Wright PetscFunctionBeginUser; 1535b881d11SJames Wright if (apply_bc.qfunction) { 154a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, apply_bc.qfunction, apply_bc.qfunction_loc, qf_apply_bc)); 155a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_apply_bc, apply_bc.qfunction_context)); 15691c97f41SJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf_apply_bc, 0)); 157a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "q", num_comp_q, CEED_EVAL_INTERP)); 158a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "Grad_q", num_comp_q * dim_sur, CEED_EVAL_GRAD)); 159a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE)); 160a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "x", num_comp_x, CEED_EVAL_INTERP)); 161a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc, "v", num_comp_q, CEED_EVAL_INTERP)); 162a424bcd0SJames Wright if (jac_data_size_sur) PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE)); 1635b881d11SJames Wright } 1645b881d11SJames Wright if (apply_bc_jacobian.qfunction) { 165a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, apply_bc_jacobian.qfunction, apply_bc_jacobian.qfunction_loc, qf_apply_bc_jacobian)); 166a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_apply_bc_jacobian, apply_bc_jacobian.qfunction_context)); 16791c97f41SJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf_apply_bc_jacobian, 0)); 168a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "dq", num_comp_q, CEED_EVAL_INTERP)); 169a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "Grad_dq", num_comp_q * dim_sur, CEED_EVAL_GRAD)); 170a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface qdata", q_data_size_sur, CEED_EVAL_NONE)); 171a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "x", num_comp_x, CEED_EVAL_INTERP)); 172a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE)); 173a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc_jacobian, "v", num_comp_q, CEED_EVAL_INTERP)); 1745b881d11SJames Wright } 175ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 1765b881d11SJames Wright } 1775b881d11SJames Wright 178a5b0ec6fSJames Wright // Utility function to add boundary operators to the composite operator 179e9c36be0SJames Wright static PetscErrorCode AddBCSubOperators(User user, Ceed ceed, DM dm, SimpleBC bc, ProblemData problem, CeedData ceed_data, CeedOperator op_apply, 180e9c36be0SJames Wright CeedOperator op_apply_ijacobian) { 181e9c36be0SJames Wright CeedInt height = 1, num_comp_q, num_comp_x; 182*cfb075a4SJames Wright CeedInt P_sur = user->app_ctx->degree + 1, Q_sur = P_sur + user->app_ctx->q_extra, dim_sur, q_data_size_sur; 183*cfb075a4SJames Wright const CeedInt jac_data_size_sur = user->phys->implicit ? problem->jac_data_size_sur : 0; 184e9c36be0SJames Wright PetscInt dim; 185e9c36be0SJames Wright DMLabel face_sets_label; 186e9c36be0SJames Wright CeedBasis basis_q_sur, basis_x_sur; 187e9c36be0SJames Wright 188e9c36be0SJames Wright PetscFunctionBeginUser; 189e9c36be0SJames Wright PetscCall(DMGetDimension(dm, &dim)); 190*cfb075a4SJames Wright PetscCall(QDataBoundaryGetNumComponents(dm, &q_data_size_sur)); 191e9c36be0SJames Wright dim_sur = dim - height; 192e9c36be0SJames Wright { // Get number of components and coordinate dimension from op_apply 193e9c36be0SJames Wright CeedOperator *sub_ops; 194e9c36be0SJames Wright CeedOperatorField field; 195e9c36be0SJames Wright PetscInt sub_op_index = 0; // will be 0 for the volume op 196e9c36be0SJames Wright CeedElemRestriction elem_restr_q, elem_restr_x; 197e9c36be0SJames Wright 198e9c36be0SJames Wright PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(op_apply, &sub_ops)); 199e9c36be0SJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &field)); 200e9c36be0SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_q)); 201e9c36be0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q)); 202e9c36be0SJames Wright 203e9c36be0SJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "x", &field)); 204e9c36be0SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_x)); 205e9c36be0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_x, &num_comp_x)); 206e9c36be0SJames Wright } 207e9c36be0SJames Wright 208e9c36be0SJames Wright { // Get bases 209e9c36be0SJames Wright DM dm_coord; 210e9c36be0SJames Wright 211e9c36be0SJames Wright PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 212e9c36be0SJames Wright DMLabel label = NULL; 213e9c36be0SJames Wright PetscInt label_value = 0; 214a5b0ec6fSJames Wright PetscInt field = 0; 215e9c36be0SJames Wright PetscCall(CreateBasisFromPlex(ceed, dm, label, label_value, height, field, &basis_q_sur)); 216e9c36be0SJames Wright PetscCall(CreateBasisFromPlex(ceed, dm_coord, label, label_value, height, field, &basis_x_sur)); 217e9c36be0SJames Wright } 218e9c36be0SJames Wright 219e9c36be0SJames Wright PetscCall(DMGetLabel(dm, "Face Sets", &face_sets_label)); 220e9c36be0SJames Wright 221e9c36be0SJames Wright { // --- Create Sub-Operator for inflow boundaries 222e9c36be0SJames Wright CeedQFunction qf_apply_inflow = NULL, qf_apply_inflow_jacobian = NULL; 223e9c36be0SJames Wright 224e9c36be0SJames Wright PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_inflow, 225e9c36be0SJames Wright problem->apply_inflow_jacobian, &qf_apply_inflow, &qf_apply_inflow_jacobian)); 226e9c36be0SJames Wright for (CeedInt i = 0; i < bc->num_inflow; i++) { 227e9c36be0SJames Wright PetscCall(AddBCSubOperator(ceed, dm, ceed_data, face_sets_label, bc->inflows[i], height, Q_sur, q_data_size_sur, jac_data_size_sur, basis_q_sur, 228e9c36be0SJames Wright basis_x_sur, qf_apply_inflow, qf_apply_inflow_jacobian, op_apply, op_apply_ijacobian)); 229e9c36be0SJames Wright } 230e9c36be0SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_inflow)); 231e9c36be0SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_inflow_jacobian)); 232e9c36be0SJames Wright } 233e9c36be0SJames Wright 234e9c36be0SJames Wright { // --- Create Sub-Operator for outflow boundaries 235e9c36be0SJames Wright CeedQFunction qf_apply_outflow = NULL, qf_apply_outflow_jacobian = NULL; 236e9c36be0SJames Wright 237e9c36be0SJames Wright PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_outflow, 238e9c36be0SJames Wright problem->apply_outflow_jacobian, &qf_apply_outflow, &qf_apply_outflow_jacobian)); 239e9c36be0SJames Wright for (CeedInt i = 0; i < bc->num_outflow; i++) { 240e9c36be0SJames Wright PetscCall(AddBCSubOperator(ceed, dm, ceed_data, face_sets_label, bc->outflows[i], height, Q_sur, q_data_size_sur, jac_data_size_sur, 241e9c36be0SJames Wright basis_q_sur, basis_x_sur, qf_apply_outflow, qf_apply_outflow_jacobian, op_apply, op_apply_ijacobian)); 242e9c36be0SJames Wright } 243e9c36be0SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_outflow)); 244e9c36be0SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_outflow_jacobian)); 245e9c36be0SJames Wright } 246e9c36be0SJames Wright 247e9c36be0SJames Wright { // --- Create Sub-Operator for freestream boundaries 248e9c36be0SJames Wright CeedQFunction qf_apply_freestream = NULL, qf_apply_freestream_jacobian = NULL; 249e9c36be0SJames Wright 250e9c36be0SJames Wright PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_freestream, 251e9c36be0SJames Wright problem->apply_freestream_jacobian, &qf_apply_freestream, &qf_apply_freestream_jacobian)); 252e9c36be0SJames Wright for (CeedInt i = 0; i < bc->num_freestream; i++) { 253e9c36be0SJames Wright PetscCall(AddBCSubOperator(ceed, dm, ceed_data, face_sets_label, bc->freestreams[i], height, Q_sur, q_data_size_sur, jac_data_size_sur, 254e9c36be0SJames Wright basis_q_sur, basis_x_sur, qf_apply_freestream, qf_apply_freestream_jacobian, op_apply, op_apply_ijacobian)); 255e9c36be0SJames Wright } 256e9c36be0SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_freestream)); 257e9c36be0SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_freestream_jacobian)); 258e9c36be0SJames Wright } 259e9c36be0SJames Wright 260e9c36be0SJames Wright { // --- Create Sub-Operator for slip boundaries 261e9c36be0SJames Wright CeedQFunction qf_apply_slip = NULL, qf_apply_slip_jacobian = NULL; 262e9c36be0SJames Wright 263e9c36be0SJames Wright PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_slip, 264e9c36be0SJames Wright problem->apply_slip_jacobian, &qf_apply_slip, &qf_apply_slip_jacobian)); 265e9c36be0SJames Wright for (CeedInt i = 0; i < bc->num_slip; i++) { 266e9c36be0SJames Wright PetscCall(AddBCSubOperator(ceed, dm, ceed_data, face_sets_label, bc->slips[i], height, Q_sur, q_data_size_sur, jac_data_size_sur, basis_q_sur, 267e9c36be0SJames Wright basis_x_sur, qf_apply_slip, qf_apply_slip_jacobian, op_apply, op_apply_ijacobian)); 268e9c36be0SJames Wright } 269e9c36be0SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_slip)); 270e9c36be0SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_slip_jacobian)); 271e9c36be0SJames Wright } 272e9c36be0SJames Wright 273e9c36be0SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_q_sur)); 274e9c36be0SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_sur)); 275e9c36be0SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 276e9c36be0SJames Wright } 277e9c36be0SJames Wright 278731c13d7SJames Wright PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData problem, SimpleBC bc) { 27977841947SLeila Ghaffari const PetscInt num_comp_q = 5; 280*cfb075a4SJames Wright const CeedInt dim = problem->dim, num_comp_x = problem->dim; 2811d2a9659SKenneth E. Jansen CeedInt jac_data_size_vol = num_comp_q + 6 + 3; 282a5b0ec6fSJames Wright CeedElemRestriction elem_restr_jd_i; 283a5b0ec6fSJames Wright CeedVector jac_data; 284a5b0ec6fSJames Wright CeedOperator op_ifunction_vol = NULL, op_rhs_vol = NULL, op_ijacobian_vol = NULL; 285a5b0ec6fSJames Wright 286a5b0ec6fSJames Wright PetscFunctionBeginUser; 2871d2a9659SKenneth E. Jansen 288752a08a7SJames Wright if (problem->apply_vol_ifunction.qfunction && problem->uses_newtonian) { 2891d2a9659SKenneth E. Jansen NewtonianIdealGasContext gas; 2901d2a9659SKenneth E. Jansen PetscCallCeed(ceed, CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfunction_context, CEED_MEM_HOST, &gas)); 2911d2a9659SKenneth E. Jansen jac_data_size_vol += (gas->idl_enable ? 1 : 0); 2921d2a9659SKenneth E. Jansen PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfunction_context, &gas)); 2931d2a9659SKenneth E. Jansen } 2941d2a9659SKenneth E. Jansen 295a5b0ec6fSJames Wright { // Create bases and element restrictions 296bb85d312SJames Wright DMLabel domain_label = NULL; 297bb85d312SJames Wright PetscInt label_value = 0, height = 0, dm_field = 0; 2980814d5a7SKenneth E. Jansen DM dm_coord; 2990814d5a7SKenneth E. Jansen 300a5b0ec6fSJames Wright PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 301bb85d312SJames Wright PetscCall(CreateBasisFromPlex(ceed, dm, domain_label, label_value, height, dm_field, &ceed_data->basis_q)); 302bb85d312SJames Wright PetscCall(CreateBasisFromPlex(ceed, dm_coord, domain_label, label_value, height, dm_field, &ceed_data->basis_x)); 30377841947SLeila Ghaffari 304bb85d312SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, domain_label, label_value, height, 0, &ceed_data->elem_restr_q)); 305bb85d312SJames Wright PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height, &ceed_data->elem_restr_x)); 306bb85d312SJames Wright PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, jac_data_size_vol, &elem_restr_jd_i)); 30764f98e98SJames Wright 308a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_ceed, NULL)); 309a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_dot_ceed, NULL)); 310a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->g_ceed, NULL)); 311a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_x, &ceed_data->x_coord, NULL)); 31264f98e98SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_jd_i, &jac_data, NULL)); 31377841947SLeila Ghaffari 31464f98e98SJames Wright { // -- Copy PETSc coordinate vector into CEED vector 31577841947SLeila Ghaffari Vec X_loc; 3163796c488SJed Brown DM cdm; 31764f98e98SJames Wright 3182b730f8bSJeremy L Thompson PetscCall(DMGetCellCoordinateDM(dm, &cdm)); 3192b730f8bSJeremy L Thompson if (cdm) { 3202b730f8bSJeremy L Thompson PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc)); 3212b730f8bSJeremy L Thompson } else { 3222b730f8bSJeremy L Thompson PetscCall(DMGetCoordinatesLocal(dm, &X_loc)); 3233796c488SJed Brown } 3242b730f8bSJeremy L Thompson PetscCall(VecScale(X_loc, problem->dm_scale)); 325d0593705SJames Wright PetscCall(VecCopyPetscToCeed(X_loc, ceed_data->x_coord)); 32664f98e98SJames Wright } 32777841947SLeila Ghaffari 328*cfb075a4SJames Wright PetscCall(QDataGet(ceed, dm, domain_label, label_value, ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord, 329*cfb075a4SJames Wright &ceed_data->elem_restr_qd_i, &ceed_data->q_data, &problem->q_data_size_vol)); 33064f98e98SJames Wright } 33164f98e98SJames Wright 33264f98e98SJames Wright { // -- Create QFunction for ICs 333e9c36be0SJames Wright CeedBasis basis_xc; 33464f98e98SJames Wright CeedQFunction qf_ics; 3355263e9c6SJames Wright CeedOperator op_ics; 33664f98e98SJames Wright 337e9c36be0SJames Wright PetscCallCeed(ceed, CeedBasisCreateProjection(ceed_data->basis_x, ceed_data->basis_q, &basis_xc)); 33864f98e98SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->ics.qfunction, problem->ics.qfunction_loc, &qf_ics)); 33964f98e98SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ics, problem->ics.qfunction_context)); 34064f98e98SJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ics, 0)); 34164f98e98SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ics, "x", num_comp_x, CEED_EVAL_INTERP)); 34264f98e98SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ics, "dx", num_comp_x * dim, CEED_EVAL_GRAD)); 34364f98e98SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ics, "q0", num_comp_q, CEED_EVAL_NONE)); 34464f98e98SJames Wright 34564f98e98SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ics, NULL, NULL, &op_ics)); 346e9c36be0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "x", ceed_data->elem_restr_x, basis_xc, CEED_VECTOR_ACTIVE)); 347e9c36be0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "dx", ceed_data->elem_restr_x, basis_xc, CEED_VECTOR_ACTIVE)); 348356036faSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "q0", ceed_data->elem_restr_q, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 349a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_ics, "evaluation time", &user->phys->ics_time_label)); 3505263e9c6SJames Wright PetscCall(OperatorApplyContextCreate(NULL, dm, user->ceed, op_ics, ceed_data->x_coord, NULL, NULL, user->Q_loc, &ceed_data->op_ics_ctx)); 35177841947SLeila Ghaffari 352e9c36be0SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_xc)); 35364f98e98SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ics)); 35464f98e98SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_ics)); 35564f98e98SJames Wright } 35664f98e98SJames Wright 35764f98e98SJames Wright if (problem->apply_vol_rhs.qfunction) { 35864f98e98SJames Wright CeedQFunction qf_rhs_vol; 35964f98e98SJames Wright 36064f98e98SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_rhs.qfunction, problem->apply_vol_rhs.qfunction_loc, &qf_rhs_vol)); 36164f98e98SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_vol, problem->apply_vol_rhs.qfunction_context)); 36264f98e98SJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_rhs_vol, 0)); 36364f98e98SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "q", num_comp_q, CEED_EVAL_INTERP)); 36464f98e98SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 365*cfb075a4SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "qdata", problem->q_data_size_vol, CEED_EVAL_NONE)); 36664f98e98SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "x", num_comp_x, CEED_EVAL_INTERP)); 36764f98e98SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_vol, "v", num_comp_q, CEED_EVAL_INTERP)); 36864f98e98SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD)); 36964f98e98SJames Wright 370a5b0ec6fSJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_vol, NULL, NULL, &op_rhs_vol)); 371a5b0ec6fSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 372a5b0ec6fSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 373a5b0ec6fSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data)); 374a5b0ec6fSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord)); 375a5b0ec6fSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 376a5b0ec6fSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 37764f98e98SJames Wright 37864f98e98SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_vol)); 37977841947SLeila Ghaffari } 38077841947SLeila Ghaffari 38164f98e98SJames Wright if (problem->apply_vol_ifunction.qfunction) { 38264f98e98SJames Wright CeedQFunction qf_ifunction_vol; 38364f98e98SJames Wright 38464f98e98SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ifunction.qfunction, problem->apply_vol_ifunction.qfunction_loc, 38564f98e98SJames Wright &qf_ifunction_vol)); 38664f98e98SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ifunction_vol, problem->apply_vol_ifunction.qfunction_context)); 38764f98e98SJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ifunction_vol, 0)); 38864f98e98SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "q", num_comp_q, CEED_EVAL_INTERP)); 38964f98e98SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 39064f98e98SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "q dot", num_comp_q, CEED_EVAL_INTERP)); 391*cfb075a4SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "qdata", problem->q_data_size_vol, CEED_EVAL_NONE)); 39264f98e98SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "x", num_comp_x, CEED_EVAL_INTERP)); 39364f98e98SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "v", num_comp_q, CEED_EVAL_INTERP)); 39464f98e98SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD)); 39564f98e98SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE)); 39664f98e98SJames Wright 397a5b0ec6fSJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ifunction_vol, NULL, NULL, &op_ifunction_vol)); 398a5b0ec6fSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 399a5b0ec6fSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 400a5b0ec6fSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "q dot", ceed_data->elem_restr_q, ceed_data->basis_q, user->q_dot_ceed)); 401a5b0ec6fSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data)); 402a5b0ec6fSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord)); 403a5b0ec6fSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 404a5b0ec6fSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 405a5b0ec6fSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "jac_data", elem_restr_jd_i, CEED_BASIS_NONE, jac_data)); 406a3ae0734SJed Brown 40764f98e98SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ifunction_vol)); 40877841947SLeila Ghaffari } 40977841947SLeila Ghaffari 41064f98e98SJames Wright if (problem->apply_vol_ijacobian.qfunction) { 41164f98e98SJames Wright CeedQFunction qf_ijacobian_vol; 41264f98e98SJames Wright 41364f98e98SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ijacobian.qfunction, problem->apply_vol_ijacobian.qfunction_loc, 41464f98e98SJames Wright &qf_ijacobian_vol)); 41564f98e98SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ijacobian_vol, problem->apply_vol_ijacobian.qfunction_context)); 41664f98e98SJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ijacobian_vol, 0)); 41764f98e98SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "dq", num_comp_q, CEED_EVAL_INTERP)); 41864f98e98SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "Grad_dq", num_comp_q * dim, CEED_EVAL_GRAD)); 419*cfb075a4SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "qdata", problem->q_data_size_vol, CEED_EVAL_NONE)); 42064f98e98SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE)); 42164f98e98SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ijacobian_vol, "v", num_comp_q, CEED_EVAL_INTERP)); 42264f98e98SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ijacobian_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD)); 42364f98e98SJames Wright 424a5b0ec6fSJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ijacobian_vol, NULL, NULL, &op_ijacobian_vol)); 425a5b0ec6fSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 426a5b0ec6fSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "Grad_dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 427a5b0ec6fSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data)); 428a5b0ec6fSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "jac_data", elem_restr_jd_i, CEED_BASIS_NONE, jac_data)); 429a5b0ec6fSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 430a5b0ec6fSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 43164f98e98SJames Wright 432a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ijacobian_vol)); 433e334ad8fSJed Brown } 434e334ad8fSJed Brown 43577841947SLeila Ghaffari // -- Create and apply CEED Composite Operator for the entire domain 43677841947SLeila Ghaffari if (!user->phys->implicit) { // RHS 4379ad5e8e4SJames Wright CeedOperator op_rhs; 438e9c36be0SJames Wright 439e9c36be0SJames Wright PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &op_rhs)); 440a5b0ec6fSJames Wright PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_rhs, op_rhs_vol)); 441e9c36be0SJames Wright PetscCall(AddBCSubOperators(user, ceed, dm, bc, problem, ceed_data, op_rhs, NULL)); 442e9c36be0SJames Wright 4439ad5e8e4SJames Wright PetscCall(OperatorApplyContextCreate(dm, dm, ceed, op_rhs, user->q_ceed, user->g_ceed, user->Q_loc, NULL, &user->op_rhs_ctx)); 444e9c36be0SJames Wright 445e9c36be0SJames Wright // ----- Get Context Labels for Operator 446e9c36be0SJames Wright PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_rhs, "solution time", &user->phys->solution_time_label)); 447e9c36be0SJames Wright PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_rhs, "timestep size", &user->phys->timestep_size_label)); 448e9c36be0SJames Wright 449a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs)); 450b4e9a8f8SJames Wright PetscCall(CreateKSPMass(user, problem)); 451577c3361SJames Wright PetscCheck(app_ctx->sgs_model_type == SGS_MODEL_NONE, user->comm, PETSC_ERR_SUP, "SGS modeling not implemented for explicit timestepping"); 45277841947SLeila Ghaffari } else { // IFunction 45391c97f41SJames Wright CeedOperator op_ijacobian = NULL; 45491c97f41SJames Wright 455e9c36be0SJames Wright // Create Composite Operaters 456e9c36be0SJames Wright PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &user->op_ifunction)); 457a5b0ec6fSJames Wright PetscCallCeed(ceed, CeedCompositeOperatorAddSub(user->op_ifunction, op_ifunction_vol)); 458e9c36be0SJames Wright if (op_ijacobian_vol) { 459e9c36be0SJames Wright PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &op_ijacobian)); 460e9c36be0SJames Wright PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_ijacobian, op_ijacobian_vol)); 461e9c36be0SJames Wright } 462e9c36be0SJames Wright PetscCall(AddBCSubOperators(user, ceed, dm, bc, problem, ceed_data, user->op_ifunction, op_ijacobian)); 463e9c36be0SJames Wright 464e9c36be0SJames Wright // ----- Get Context Labels for Operator 465e9c36be0SJames Wright PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(user->op_ifunction, "solution time", &user->phys->solution_time_label)); 466e9c36be0SJames Wright PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(user->op_ifunction, "timestep size", &user->phys->timestep_size_label)); 467e9c36be0SJames Wright 46891c97f41SJames Wright if (op_ijacobian) { 46991c97f41SJames Wright PetscCall(MatCeedCreate(user->dm, user->dm, op_ijacobian, NULL, &user->mat_ijacobian)); 47091c97f41SJames Wright PetscCall(MatCeedSetLocalVectors(user->mat_ijacobian, user->Q_dot_loc, NULL)); 47191c97f41SJames Wright PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_ijacobian, "ijacobian time shift", &user->phys->ijacobian_time_shift_label)); 47291c97f41SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_ijacobian)); 473e334ad8fSJed Brown } 4748f5ab23bSJames Wright if (app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) PetscCall(SgsDDSetup(ceed, user, ceed_data, problem)); 475dada6cc0SJames Wright } 476f5452247SJames Wright 477577c3361SJames Wright if (problem->use_strong_bc_ceed) PetscCall(SetupStrongBC_Ceed(ceed, ceed_data, dm, user, problem, bc)); 478f5452247SJames Wright if (app_ctx->turb_spanstats_enable) PetscCall(TurbulenceStatisticsSetup(ceed, user, ceed_data, problem)); 479ee3b213aSJames Wright if (app_ctx->diff_filter_monitor && !user->diff_filter) PetscCall(DifferentialFilterSetup(ceed, user, ceed_data, problem)); 48028134cfaSJames Wright if (app_ctx->sgs_train_enable) PetscCall(SGS_DD_TrainingSetup(ceed, user, ceed_data, problem)); 481f5452247SJames Wright 482a5b0ec6fSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&jac_data)); 483a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_jd_i)); 484a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_ijacobian_vol)); 485a5b0ec6fSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_ifunction_vol)); 486a5b0ec6fSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_vol)); 487ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 48877841947SLeila Ghaffari } 489