1dc936754SJeremy L Thompson // Copyright (c) 2017-2024, 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 /// Setup libCEED for Navier-Stokes example using PETSc 10a515125bSLeila Ghaffari 11e419654dSJeremy L Thompson #include <ceed.h> 12e419654dSJeremy L Thompson #include <petscdmplex.h> 13e419654dSJeremy L Thompson 14a515125bSLeila Ghaffari #include "../navierstokes.h" 15a515125bSLeila Ghaffari 16a78efa86SJames Wright // @brief Create CeedOperator for unstabilized mass KSP for explicit timestepping 17a78efa86SJames Wright static PetscErrorCode CreateKSPMassOperator_Unstabilized(User user, CeedOperator *op_mass) { 18a78efa86SJames Wright Ceed ceed = user->ceed; 19a78efa86SJames Wright CeedInt num_comp_q, q_data_size; 20a78efa86SJames Wright CeedQFunction qf_mass; 21a78efa86SJames Wright CeedElemRestriction elem_restr_q, elem_restr_qd_i; 22a78efa86SJames Wright CeedBasis basis_q; 23a78efa86SJames Wright CeedVector q_data; 24a78efa86SJames Wright 25a78efa86SJames Wright PetscFunctionBeginUser; 26a78efa86SJames Wright { // Get restriction and basis from the RHS function 27a78efa86SJames Wright CeedOperator *sub_ops; 28a78efa86SJames Wright CeedOperatorField field; 29a78efa86SJames Wright PetscInt sub_op_index = 0; // will be 0 for the volume op 30a78efa86SJames Wright 31a78efa86SJames Wright PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(user->op_rhs_ctx->op, &sub_ops)); 32a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &field)); 33a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_q)); 34a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetBasis(field, &basis_q)); 35a78efa86SJames Wright 36a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &field)); 37a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_qd_i)); 38a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetVector(field, &q_data)); 39a78efa86SJames Wright } 40a78efa86SJames Wright 41a78efa86SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q)); 42a78efa86SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd_i, &q_data_size)); 43a78efa86SJames Wright 44a78efa86SJames Wright PetscCall(CreateMassQFunction(ceed, num_comp_q, q_data_size, &qf_mass)); 45a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass)); 46a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "u", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 47a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data)); 48a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 49a78efa86SJames Wright 50a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 51a78efa86SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 52a78efa86SJames Wright } 53a78efa86SJames Wright 54a78efa86SJames Wright // @brief Create KSP to solve the inverse mass operator for explicit time stepping schemes 55991aef52SJames Wright static PetscErrorCode CreateKSPMass(User user, ProblemData problem) { 56a78efa86SJames Wright Ceed ceed = user->ceed; 57a78efa86SJames Wright DM dm = user->dm; 58a78efa86SJames Wright CeedOperator op_mass; 59a78efa86SJames Wright 60a78efa86SJames Wright PetscFunctionBeginUser; 61a78efa86SJames Wright if (problem->create_mass_operator) PetscCall(problem->create_mass_operator(user, &op_mass)); 62a78efa86SJames Wright else PetscCall(CreateKSPMassOperator_Unstabilized(user, &op_mass)); 63a78efa86SJames Wright 64a78efa86SJames Wright { // -- Setup KSP for mass operator 65a78efa86SJames Wright Mat mat_mass; 66a78efa86SJames Wright Vec Zeros_loc; 67a78efa86SJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)dm); 68a78efa86SJames Wright 69a78efa86SJames Wright PetscCall(DMCreateLocalVector(dm, &Zeros_loc)); 70a78efa86SJames Wright PetscCall(VecZeroEntries(Zeros_loc)); 71a78efa86SJames Wright PetscCall(MatCeedCreate(dm, dm, op_mass, NULL, &mat_mass)); 72a78efa86SJames Wright PetscCall(MatCeedSetLocalVectors(mat_mass, Zeros_loc, NULL)); 73a78efa86SJames Wright 74a78efa86SJames Wright PetscCall(KSPCreate(comm, &user->mass_ksp)); 75a78efa86SJames Wright PetscCall(KSPSetOptionsPrefix(user->mass_ksp, "mass_")); 76a78efa86SJames Wright { // lumped by default 77a78efa86SJames Wright PC pc; 78a78efa86SJames Wright PetscCall(KSPGetPC(user->mass_ksp, &pc)); 79a78efa86SJames Wright PetscCall(PCSetType(pc, PCJACOBI)); 80a78efa86SJames Wright PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM)); 81a78efa86SJames Wright PetscCall(KSPSetType(user->mass_ksp, KSPPREONLY)); 82a78efa86SJames Wright } 83a78efa86SJames Wright PetscCall(KSPSetFromOptions_WithMatCeed(user->mass_ksp, mat_mass)); 84a78efa86SJames Wright PetscCall(KSPSetFromOptions(user->mass_ksp)); 85a78efa86SJames Wright PetscCall(VecDestroy(&Zeros_loc)); 86a78efa86SJames Wright PetscCall(MatDestroy(&mat_mass)); 87a78efa86SJames Wright } 88a78efa86SJames Wright 89a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass)); 90a78efa86SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 91a78efa86SJames Wright } 92a78efa86SJames Wright 9336c6cbc8SJames Wright PetscErrorCode AddBCSubOperator(Ceed ceed, DM dm, CeedData ceed_data, DMLabel domain_label, PetscInt label_value, CeedInt height, CeedInt Q_sur, 942b916ea7SJeremy L Thompson CeedInt q_data_size_sur, CeedInt jac_data_size_sur, CeedQFunction qf_apply_bc, CeedQFunction qf_apply_bc_jacobian, 95368c645fSJames Wright CeedOperator *op_apply, CeedOperator *op_apply_ijacobian) { 9615c18037SJames Wright CeedVector q_data_sur, jac_data_sur = NULL; 972b916ea7SJeremy L Thompson CeedOperator op_setup_sur, op_apply_bc, op_apply_bc_jacobian = NULL; 9815c18037SJames Wright CeedElemRestriction elem_restr_x_sur, elem_restr_q_sur, elem_restr_qd_i_sur, elem_restr_jd_i_sur = NULL; 9915c18037SJames Wright CeedInt num_qpts_sur, dm_field = 0; 100368c645fSJames Wright 10106f41313SJames Wright PetscFunctionBeginUser; 102368c645fSJames Wright // --- Get number of quadrature points for the boundaries 103b4c37c5cSJames Wright PetscCallCeed(ceed, CeedBasisGetNumQuadraturePoints(ceed_data->basis_q_sur, &num_qpts_sur)); 104368c645fSJames Wright 105368c645fSJames Wright // ---- CEED Restriction 10615c18037SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, domain_label, label_value, height, dm_field, &elem_restr_q_sur)); 10715c18037SJames Wright PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height, &elem_restr_x_sur)); 10815c18037SJames Wright PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, q_data_size_sur, &elem_restr_qd_i_sur)); 109368c645fSJames Wright if (jac_data_size_sur > 0) { 110368c645fSJames Wright // State-dependent data will be passed from residual to Jacobian. This will be collocated. 11115c18037SJames Wright PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, jac_data_size_sur, &elem_restr_jd_i_sur)); 112b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_jd_i_sur, &jac_data_sur, NULL)); 113368c645fSJames Wright } 114368c645fSJames Wright 115368c645fSJames Wright // ---- CEED Vector 116defe8520SJames Wright CeedInt loc_num_elem_sur; 117b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumElements(elem_restr_q_sur, &loc_num_elem_sur)); 118b4c37c5cSJames Wright PetscCallCeed(ceed, CeedVectorCreate(ceed, q_data_size_sur * loc_num_elem_sur * num_qpts_sur, &q_data_sur)); 119368c645fSJames Wright 120368c645fSJames Wright // ---- CEED Operator 121368c645fSJames Wright // ----- CEED Operator for Setup (geometric factors) 122b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, ceed_data->qf_setup_sur, NULL, NULL, &op_setup_sur)); 123b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "dx", elem_restr_x_sur, ceed_data->basis_x_sur, CEED_VECTOR_ACTIVE)); 124b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, ceed_data->basis_x_sur, CEED_VECTOR_NONE)); 12558e1cbfdSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 126368c645fSJames Wright 127368c645fSJames Wright // ----- CEED Operator for Physics 128b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_apply_bc, NULL, NULL, &op_apply_bc)); 129b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "q", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE)); 130b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "Grad_q", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE)); 13158e1cbfdSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_NONE, q_data_sur)); 132b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "x", elem_restr_x_sur, ceed_data->basis_x_sur, ceed_data->x_coord)); 133b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "v", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE)); 134b4c37c5cSJames Wright if (elem_restr_jd_i_sur) 13558e1cbfdSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_NONE, jac_data_sur)); 136368c645fSJames Wright 1374b96a86bSJames Wright if (qf_apply_bc_jacobian && elem_restr_jd_i_sur) { 138b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_apply_bc_jacobian, NULL, NULL, &op_apply_bc_jacobian)); 139b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "dq", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE)); 140b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "Grad_dq", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE)); 14158e1cbfdSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_NONE, q_data_sur)); 142b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "x", elem_restr_x_sur, ceed_data->basis_x_sur, ceed_data->x_coord)); 14358e1cbfdSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_NONE, jac_data_sur)); 144b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "v", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE)); 145368c645fSJames Wright } 146368c645fSJames Wright 147368c645fSJames Wright // ----- Apply CEED operator for Setup 148b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorApply(op_setup_sur, ceed_data->x_coord, q_data_sur, CEED_REQUEST_IMMEDIATE)); 149368c645fSJames Wright 150368c645fSJames Wright // ----- Apply Sub-Operator for Physics 151b4c37c5cSJames Wright PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_apply, op_apply_bc)); 152b4c37c5cSJames Wright if (op_apply_bc_jacobian) PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_apply_ijacobian, op_apply_bc_jacobian)); 153368c645fSJames Wright 154368c645fSJames Wright // ----- Cleanup 155b4c37c5cSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data_sur)); 156b4c37c5cSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&jac_data_sur)); 157b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q_sur)); 158b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x_sur)); 159b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd_i_sur)); 160b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_jd_i_sur)); 161b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup_sur)); 162b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_apply_bc)); 163b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_apply_bc_jacobian)); 164d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 165368c645fSJames Wright } 166368c645fSJames Wright 167a515125bSLeila Ghaffari // Utility function to create CEED Composite Operator for the entire domain 1682b916ea7SJeremy L Thompson PetscErrorCode CreateOperatorForDomain(Ceed ceed, DM dm, SimpleBC bc, CeedData ceed_data, Physics phys, CeedOperator op_apply_vol, 1692b916ea7SJeremy L Thompson CeedOperator op_apply_ijacobian_vol, CeedInt height, CeedInt P_sur, CeedInt Q_sur, CeedInt q_data_size_sur, 1702b916ea7SJeremy L Thompson CeedInt jac_data_size_sur, CeedOperator *op_apply, CeedOperator *op_apply_ijacobian) { 171a515125bSLeila Ghaffari DMLabel domain_label; 172a515125bSLeila Ghaffari 17306f41313SJames Wright PetscFunctionBeginUser; 174a515125bSLeila Ghaffari // Create Composite Operaters 175b4c37c5cSJames Wright PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, op_apply)); 176b4c37c5cSJames Wright if (op_apply_ijacobian) PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, op_apply_ijacobian)); 177a515125bSLeila Ghaffari 178a515125bSLeila Ghaffari // --Apply Sub-Operator for the volume 179b4c37c5cSJames Wright PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_apply, op_apply_vol)); 180b4c37c5cSJames Wright if (op_apply_ijacobian) PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_apply_ijacobian, op_apply_ijacobian_vol)); 181a515125bSLeila Ghaffari 182a515125bSLeila Ghaffari // -- Create Sub-Operator for in/outflow BCs 1832b916ea7SJeremy L Thompson PetscCall(DMGetLabel(dm, "Face Sets", &domain_label)); 184a515125bSLeila Ghaffari 185002797a3SLeila Ghaffari // --- Create Sub-Operator for inflow boundaries 186002797a3SLeila Ghaffari for (CeedInt i = 0; i < bc->num_inflow; i++) { 1872b916ea7SJeremy L Thompson PetscCall(AddBCSubOperator(ceed, dm, ceed_data, domain_label, bc->inflows[i], height, Q_sur, q_data_size_sur, jac_data_size_sur, 1882b916ea7SJeremy L Thompson ceed_data->qf_apply_inflow, ceed_data->qf_apply_inflow_jacobian, op_apply, op_apply_ijacobian)); 189a515125bSLeila Ghaffari } 190002797a3SLeila Ghaffari // --- Create Sub-Operator for outflow boundaries 191002797a3SLeila Ghaffari for (CeedInt i = 0; i < bc->num_outflow; i++) { 1922b916ea7SJeremy L Thompson PetscCall(AddBCSubOperator(ceed, dm, ceed_data, domain_label, bc->outflows[i], height, Q_sur, q_data_size_sur, jac_data_size_sur, 1932b916ea7SJeremy L Thompson ceed_data->qf_apply_outflow, ceed_data->qf_apply_outflow_jacobian, op_apply, op_apply_ijacobian)); 1942556a851SJed Brown } 195df55ba5fSJames Wright // --- Create Sub-Operator for freestream boundaries 196df55ba5fSJames Wright for (CeedInt i = 0; i < bc->num_freestream; i++) { 1972b916ea7SJeremy L Thompson PetscCall(AddBCSubOperator(ceed, dm, ceed_data, domain_label, bc->freestreams[i], height, Q_sur, q_data_size_sur, jac_data_size_sur, 1982b916ea7SJeremy L Thompson ceed_data->qf_apply_freestream, ceed_data->qf_apply_freestream_jacobian, op_apply, op_apply_ijacobian)); 199002797a3SLeila Ghaffari } 2009ed3d70dSJames Wright // --- Create Sub-Operator for slip boundaries 2019ed3d70dSJames Wright for (CeedInt i = 0; i < bc->num_slip; i++) { 2029ed3d70dSJames Wright PetscCall(AddBCSubOperator(ceed, dm, ceed_data, domain_label, bc->slips[i], height, Q_sur, q_data_size_sur, jac_data_size_sur, 2039ed3d70dSJames Wright ceed_data->qf_apply_slip, ceed_data->qf_apply_slip_jacobian, op_apply, op_apply_ijacobian)); 2049ed3d70dSJames Wright } 20592ada588SJames Wright 20692ada588SJames Wright // ----- Get Context Labels for Operator 207b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(*op_apply, "solution time", &phys->solution_time_label)); 208b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(*op_apply, "timestep size", &phys->timestep_size_label)); 209d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 210a515125bSLeila Ghaffari } 211a515125bSLeila Ghaffari 2122b916ea7SJeremy L Thompson PetscErrorCode SetupBCQFunctions(Ceed ceed, PetscInt dim_sur, PetscInt num_comp_x, PetscInt num_comp_q, PetscInt q_data_size_sur, 2132b916ea7SJeremy L Thompson PetscInt jac_data_size_sur, ProblemQFunctionSpec apply_bc, ProblemQFunctionSpec apply_bc_jacobian, 21425988f00SJames Wright CeedQFunction *qf_apply_bc, CeedQFunction *qf_apply_bc_jacobian) { 21525988f00SJames Wright PetscFunctionBeginUser; 21625988f00SJames Wright if (apply_bc.qfunction) { 217b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, apply_bc.qfunction, apply_bc.qfunction_loc, qf_apply_bc)); 218b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_apply_bc, apply_bc.qfunction_context)); 219ebfabadfSJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf_apply_bc, 0)); 220b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "q", num_comp_q, CEED_EVAL_INTERP)); 221b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "Grad_q", num_comp_q * dim_sur, CEED_EVAL_GRAD)); 222b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE)); 223b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "x", num_comp_x, CEED_EVAL_INTERP)); 224b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc, "v", num_comp_q, CEED_EVAL_INTERP)); 225b4c37c5cSJames Wright if (jac_data_size_sur) PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE)); 22625988f00SJames Wright } 22725988f00SJames Wright if (apply_bc_jacobian.qfunction) { 228b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, apply_bc_jacobian.qfunction, apply_bc_jacobian.qfunction_loc, qf_apply_bc_jacobian)); 229b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_apply_bc_jacobian, apply_bc_jacobian.qfunction_context)); 230ebfabadfSJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf_apply_bc_jacobian, 0)); 231b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "dq", num_comp_q, CEED_EVAL_INTERP)); 232b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "Grad_dq", num_comp_q * dim_sur, CEED_EVAL_GRAD)); 233b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface qdata", q_data_size_sur, CEED_EVAL_NONE)); 234b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "x", num_comp_x, CEED_EVAL_INTERP)); 235b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE)); 236b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc_jacobian, "v", num_comp_q, CEED_EVAL_INTERP)); 23725988f00SJames Wright } 238d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 23925988f00SJames Wright } 24025988f00SJames Wright 241991aef52SJames Wright PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData problem, SimpleBC bc) { 242a515125bSLeila Ghaffari PetscFunctionBeginUser; 243a515125bSLeila Ghaffari const PetscInt num_comp_q = 5; 24494a7b3d2SKenneth E. Jansen const CeedInt dim = problem->dim, num_comp_x = problem->dim, q_data_size_vol = problem->q_data_size_vol; 24594a7b3d2SKenneth E. Jansen CeedInt jac_data_size_vol = num_comp_q + 6 + 3; 24694a7b3d2SKenneth E. Jansen 2477d8a615bSJames Wright if (problem->apply_vol_ifunction.qfunction && problem->uses_newtonian) { 24894a7b3d2SKenneth E. Jansen NewtonianIdealGasContext gas; 24994a7b3d2SKenneth E. Jansen PetscCallCeed(ceed, CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfunction_context, CEED_MEM_HOST, &gas)); 25094a7b3d2SKenneth E. Jansen jac_data_size_vol += (gas->idl_enable ? 1 : 0); 25194a7b3d2SKenneth E. Jansen PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfunction_context, &gas)); 25294a7b3d2SKenneth E. Jansen } 25394a7b3d2SKenneth E. Jansen 254752f40e3SJed Brown CeedElemRestriction elem_restr_jd_i; 255752f40e3SJed Brown CeedVector jac_data; 25667263decSKenneth E. Jansen CeedInt num_qpts; 25715c18037SJames Wright DMLabel domain_label = NULL; 25815c18037SJames Wright PetscInt label_value = 0, height = 0, dm_field = 0; 259a515125bSLeila Ghaffari 26067263decSKenneth E. Jansen DM dm_coord; 26167263decSKenneth E. Jansen PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 26267263decSKenneth E. Jansen 26315c18037SJames Wright PetscCall(CreateBasisFromPlex(ceed, dm, domain_label, label_value, height, dm_field, &ceed_data->basis_q)); 26415c18037SJames Wright PetscCall(CreateBasisFromPlex(ceed, dm_coord, domain_label, label_value, height, dm_field, &ceed_data->basis_x)); 265b4c37c5cSJames Wright PetscCallCeed(ceed, CeedBasisCreateProjection(ceed_data->basis_x, ceed_data->basis_q, &ceed_data->basis_xc)); 266b4c37c5cSJames Wright PetscCallCeed(ceed, CeedBasisGetNumQuadraturePoints(ceed_data->basis_q, &num_qpts)); 267a515125bSLeila Ghaffari 26815c18037SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, domain_label, label_value, height, 0, &ceed_data->elem_restr_q)); 26915c18037SJames Wright PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height, &ceed_data->elem_restr_x)); 27015c18037SJames Wright PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, q_data_size_vol, &ceed_data->elem_restr_qd_i)); 27115c18037SJames Wright PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, jac_data_size_vol, &elem_restr_jd_i)); 272*ce8bebb6SJames Wright 273b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_ceed, NULL)); 274b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_dot_ceed, NULL)); 275b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->g_ceed, NULL)); 276b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_x, &ceed_data->x_coord, NULL)); 277*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_qd_i, &ceed_data->q_data, NULL)); 278*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_jd_i, &jac_data, NULL)); 279a515125bSLeila Ghaffari 280*ce8bebb6SJames Wright { // -- Copy PETSc coordinate vector into CEED vector 281a515125bSLeila Ghaffari Vec X_loc; 282cac8aa24SJed Brown DM cdm; 283*ce8bebb6SJames Wright 2842b916ea7SJeremy L Thompson PetscCall(DMGetCellCoordinateDM(dm, &cdm)); 2852b916ea7SJeremy L Thompson if (cdm) { 2862b916ea7SJeremy L Thompson PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc)); 2872b916ea7SJeremy L Thompson } else { 2882b916ea7SJeremy L Thompson PetscCall(DMGetCoordinatesLocal(dm, &X_loc)); 289cac8aa24SJed Brown } 2902b916ea7SJeremy L Thompson PetscCall(VecScale(X_loc, problem->dm_scale)); 291a7dac1d5SJames Wright PetscCall(VecCopyPetscToCeed(X_loc, ceed_data->x_coord)); 292*ce8bebb6SJames Wright } 293a515125bSLeila Ghaffari 294*ce8bebb6SJames Wright { // -- Create quadrature data 295*ce8bebb6SJames Wright CeedQFunction qf_setup_vol; 296*ce8bebb6SJames Wright CeedOperator op_setup_vol; 2970b0430b7SJames Wright 298*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->setup_vol.qfunction, problem->setup_vol.qfunction_loc, &qf_setup_vol)); 299*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_setup_vol, problem->setup_vol.qfunction_context)); 300*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_setup_vol, 0)); 301*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup_vol, "dx", num_comp_x * dim, CEED_EVAL_GRAD)); 302*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup_vol, "weight", 1, CEED_EVAL_WEIGHT)); 303*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_setup_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE)); 304a515125bSLeila Ghaffari 305*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_setup_vol, NULL, NULL, &op_setup_vol)); 306*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_setup_vol, "dx", ceed_data->elem_restr_x, ceed_data->basis_x, CEED_VECTOR_ACTIVE)); 307*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_setup_vol, "weight", CEED_ELEMRESTRICTION_NONE, ceed_data->basis_x, CEED_VECTOR_NONE)); 308*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_setup_vol, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 309*ce8bebb6SJames Wright 310*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedOperatorApply(op_setup_vol, ceed_data->x_coord, ceed_data->q_data, CEED_REQUEST_IMMEDIATE)); 311*ce8bebb6SJames Wright 312*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_setup_vol)); 313*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup_vol)); 314*ce8bebb6SJames Wright } 315*ce8bebb6SJames Wright 316*ce8bebb6SJames Wright { // -- Create QFunction for ICs 317*ce8bebb6SJames Wright CeedQFunction qf_ics; 3188f18bb8bSJames Wright CeedOperator op_ics; 319*ce8bebb6SJames Wright 320*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->ics.qfunction, problem->ics.qfunction_loc, &qf_ics)); 321*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ics, problem->ics.qfunction_context)); 322*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ics, 0)); 323*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ics, "x", num_comp_x, CEED_EVAL_INTERP)); 324*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ics, "dx", num_comp_x * dim, CEED_EVAL_GRAD)); 325*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ics, "q0", num_comp_q, CEED_EVAL_NONE)); 326*ce8bebb6SJames Wright 327*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ics, NULL, NULL, &op_ics)); 328b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "x", ceed_data->elem_restr_x, ceed_data->basis_xc, CEED_VECTOR_ACTIVE)); 329b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "dx", ceed_data->elem_restr_x, ceed_data->basis_xc, CEED_VECTOR_ACTIVE)); 33058e1cbfdSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "q0", ceed_data->elem_restr_q, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 331b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_ics, "evaluation time", &user->phys->ics_time_label)); 3328f18bb8bSJames Wright PetscCall(OperatorApplyContextCreate(NULL, dm, user->ceed, op_ics, ceed_data->x_coord, NULL, NULL, user->Q_loc, &ceed_data->op_ics_ctx)); 333a515125bSLeila Ghaffari 334*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ics)); 335*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_ics)); 336*ce8bebb6SJames Wright } 337*ce8bebb6SJames Wright 338*ce8bebb6SJames Wright if (problem->apply_vol_rhs.qfunction) { 339*ce8bebb6SJames Wright CeedQFunction qf_rhs_vol; 340a515125bSLeila Ghaffari CeedOperator op; 341*ce8bebb6SJames Wright 342*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_rhs.qfunction, problem->apply_vol_rhs.qfunction_loc, &qf_rhs_vol)); 343*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_vol, problem->apply_vol_rhs.qfunction_context)); 344*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_rhs_vol, 0)); 345*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "q", num_comp_q, CEED_EVAL_INTERP)); 346*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 347*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE)); 348*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "x", num_comp_x, CEED_EVAL_INTERP)); 349*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_vol, "v", num_comp_q, CEED_EVAL_INTERP)); 350*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD)); 351*ce8bebb6SJames Wright 352*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_vol, NULL, NULL, &op)); 353b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 354b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 35558e1cbfdSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data)); 356b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord)); 357b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 358b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 359a515125bSLeila Ghaffari user->op_rhs_vol = op; 360*ce8bebb6SJames Wright 361*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_vol)); 362a515125bSLeila Ghaffari } 363a515125bSLeila Ghaffari 364*ce8bebb6SJames Wright if (problem->apply_vol_ifunction.qfunction) { 365*ce8bebb6SJames Wright CeedQFunction qf_ifunction_vol; 366a515125bSLeila Ghaffari CeedOperator op; 367*ce8bebb6SJames Wright 368*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ifunction.qfunction, problem->apply_vol_ifunction.qfunction_loc, 369*ce8bebb6SJames Wright &qf_ifunction_vol)); 370*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ifunction_vol, problem->apply_vol_ifunction.qfunction_context)); 371*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ifunction_vol, 0)); 372*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "q", num_comp_q, CEED_EVAL_INTERP)); 373*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 374*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "q dot", num_comp_q, CEED_EVAL_INTERP)); 375*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE)); 376*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "x", num_comp_x, CEED_EVAL_INTERP)); 377*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "v", num_comp_q, CEED_EVAL_INTERP)); 378*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD)); 379*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE)); 380*ce8bebb6SJames Wright 381*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ifunction_vol, NULL, NULL, &op)); 382b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 383b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 384b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op, "q dot", ceed_data->elem_restr_q, ceed_data->basis_q, user->q_dot_ceed)); 38558e1cbfdSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data)); 386b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord)); 387b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 388b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 38958e1cbfdSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op, "jac_data", elem_restr_jd_i, CEED_BASIS_NONE, jac_data)); 390752f40e3SJed Brown 391a515125bSLeila Ghaffari user->op_ifunction_vol = op; 392*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ifunction_vol)); 393a515125bSLeila Ghaffari } 394a515125bSLeila Ghaffari 395f0b65372SJed Brown CeedOperator op_ijacobian_vol = NULL; 396*ce8bebb6SJames Wright if (problem->apply_vol_ijacobian.qfunction) { 397*ce8bebb6SJames Wright CeedQFunction qf_ijacobian_vol; 398f0b65372SJed Brown CeedOperator op; 399*ce8bebb6SJames Wright 400*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ijacobian.qfunction, problem->apply_vol_ijacobian.qfunction_loc, 401*ce8bebb6SJames Wright &qf_ijacobian_vol)); 402*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ijacobian_vol, problem->apply_vol_ijacobian.qfunction_context)); 403*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ijacobian_vol, 0)); 404*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "dq", num_comp_q, CEED_EVAL_INTERP)); 405*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "Grad_dq", num_comp_q * dim, CEED_EVAL_GRAD)); 406*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE)); 407*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE)); 408*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ijacobian_vol, "v", num_comp_q, CEED_EVAL_INTERP)); 409*ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ijacobian_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD)); 410*ce8bebb6SJames Wright 411b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ijacobian_vol, NULL, NULL, &op)); 412b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op, "dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 413b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 41458e1cbfdSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data)); 41558e1cbfdSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op, "jac_data", elem_restr_jd_i, CEED_BASIS_NONE, jac_data)); 416b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 417b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 418f0b65372SJed Brown op_ijacobian_vol = op; 419*ce8bebb6SJames Wright 420b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ijacobian_vol)); 421f0b65372SJed Brown } 422f0b65372SJed Brown 423a515125bSLeila Ghaffari // ***************************************************************************** 424a515125bSLeila Ghaffari // Set up CEED objects for the exterior domain (surface) 425a515125bSLeila Ghaffari // ***************************************************************************** 42615c18037SJames Wright height = 1; 42715c18037SJames Wright CeedInt dim_sur = dim - height, P_sur = app_ctx->degree + 1, Q_sur = P_sur + app_ctx->q_extra; 4284b96a86bSJames Wright const CeedInt q_data_size_sur = problem->q_data_size_sur, jac_data_size_sur = user->phys->implicit ? problem->jac_data_size_sur : 0; 429a515125bSLeila Ghaffari 430a515125bSLeila Ghaffari // ----------------------------------------------------------------------------- 431a515125bSLeila Ghaffari // CEED Bases 432a515125bSLeila Ghaffari // ----------------------------------------------------------------------------- 43367263decSKenneth E. Jansen 434*ce8bebb6SJames Wright { 435*ce8bebb6SJames Wright DM dm_coord; 436*ce8bebb6SJames Wright 437*ce8bebb6SJames Wright PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 438*ce8bebb6SJames Wright DMLabel label = NULL; 43967263decSKenneth E. Jansen PetscInt face_id = 0; 44067263decSKenneth E. Jansen PetscInt field = 0; // Still want the normal, default field 44167263decSKenneth E. Jansen PetscCall(CreateBasisFromPlex(ceed, dm, label, face_id, height, field, &ceed_data->basis_q_sur)); 44267263decSKenneth E. Jansen PetscCall(CreateBasisFromPlex(ceed, dm_coord, label, face_id, height, field, &ceed_data->basis_x_sur)); 443*ce8bebb6SJames Wright } 444a515125bSLeila Ghaffari 445a515125bSLeila Ghaffari // ----------------------------------------------------------------------------- 446a515125bSLeila Ghaffari // CEED QFunctions 447a515125bSLeila Ghaffari // ----------------------------------------------------------------------------- 448a515125bSLeila Ghaffari // -- Create QFunction for quadrature data 449b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->setup_sur.qfunction, problem->setup_sur.qfunction_loc, &ceed_data->qf_setup_sur)); 45015a3537eSJed Brown if (problem->setup_sur.qfunction_context) { 451b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(ceed_data->qf_setup_sur, problem->setup_sur.qfunction_context)); 45215a3537eSJed Brown } 453ebfabadfSJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(ceed_data->qf_setup_sur, 0)); 454b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_setup_sur, "dx", num_comp_x * dim_sur, CEED_EVAL_GRAD)); 455b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_setup_sur, "weight", 1, CEED_EVAL_WEIGHT)); 456b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(ceed_data->qf_setup_sur, "surface qdata", q_data_size_sur, CEED_EVAL_NONE)); 457a515125bSLeila Ghaffari 4582b916ea7SJeremy L Thompson PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_inflow, 4592b916ea7SJeremy L Thompson problem->apply_inflow_jacobian, &ceed_data->qf_apply_inflow, &ceed_data->qf_apply_inflow_jacobian)); 4602b916ea7SJeremy L Thompson PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_outflow, 4612b916ea7SJeremy L Thompson problem->apply_outflow_jacobian, &ceed_data->qf_apply_outflow, &ceed_data->qf_apply_outflow_jacobian)); 4622b916ea7SJeremy L Thompson PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_freestream, 4632b916ea7SJeremy L Thompson problem->apply_freestream_jacobian, &ceed_data->qf_apply_freestream, &ceed_data->qf_apply_freestream_jacobian)); 4649ed3d70dSJames Wright PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_slip, 4659ed3d70dSJames Wright problem->apply_slip_jacobian, &ceed_data->qf_apply_slip, &ceed_data->qf_apply_slip_jacobian)); 466a515125bSLeila Ghaffari 467a515125bSLeila Ghaffari // ***************************************************************************** 468a515125bSLeila Ghaffari // CEED Operator Apply 469a515125bSLeila Ghaffari // ***************************************************************************** 470a515125bSLeila Ghaffari // -- Create and apply CEED Composite Operator for the entire domain 471a515125bSLeila Ghaffari if (!user->phys->implicit) { // RHS 472da5fe0e4SJames Wright CeedOperator op_rhs; 473da5fe0e4SJames Wright PetscCall(CreateOperatorForDomain(ceed, dm, bc, ceed_data, user->phys, user->op_rhs_vol, NULL, height, P_sur, Q_sur, q_data_size_sur, 0, &op_rhs, 474da5fe0e4SJames Wright NULL)); 475da5fe0e4SJames Wright PetscCall(OperatorApplyContextCreate(dm, dm, ceed, op_rhs, user->q_ceed, user->g_ceed, user->Q_loc, NULL, &user->op_rhs_ctx)); 476b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs)); 477a78efa86SJames Wright PetscCall(CreateKSPMass(user, problem)); 478c2376cc9SJames Wright PetscCheck(app_ctx->sgs_model_type == SGS_MODEL_NONE, user->comm, PETSC_ERR_SUP, "SGS modeling not implemented for explicit timestepping"); 479a515125bSLeila Ghaffari } else { // IFunction 480ebfabadfSJames Wright CeedOperator op_ijacobian = NULL; 481ebfabadfSJames Wright 4822b916ea7SJeremy L Thompson PetscCall(CreateOperatorForDomain(ceed, dm, bc, ceed_data, user->phys, user->op_ifunction_vol, op_ijacobian_vol, height, P_sur, Q_sur, 483ebfabadfSJames Wright q_data_size_sur, jac_data_size_sur, &user->op_ifunction, op_ijacobian_vol ? &op_ijacobian : NULL)); 484ebfabadfSJames Wright if (op_ijacobian) { 485ebfabadfSJames Wright PetscCall(MatCeedCreate(user->dm, user->dm, op_ijacobian, NULL, &user->mat_ijacobian)); 486ebfabadfSJames Wright PetscCall(MatCeedSetLocalVectors(user->mat_ijacobian, user->Q_dot_loc, NULL)); 487ebfabadfSJames Wright PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_ijacobian, "ijacobian time shift", &user->phys->ijacobian_time_shift_label)); 488ebfabadfSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_ijacobian)); 489f0b65372SJed Brown } 490ad494f68SJames Wright if (app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) PetscCall(SgsDDSetup(ceed, user, ceed_data, problem)); 4916d0190e2SJames Wright } 49291933550SJames Wright 493c2376cc9SJames Wright if (problem->use_strong_bc_ceed) PetscCall(SetupStrongBC_Ceed(ceed, ceed_data, dm, user, problem, bc)); 49491933550SJames Wright if (app_ctx->turb_spanstats_enable) PetscCall(TurbulenceStatisticsSetup(ceed, user, ceed_data, problem)); 495aa0b7f76SJames Wright if (app_ctx->diff_filter_monitor && !user->diff_filter) PetscCall(DifferentialFilterSetup(ceed, user, ceed_data, problem)); 4961c17f66aSJames Wright if (app_ctx->sgs_train_enable) PetscCall(SGS_DD_TrainingSetup(ceed, user, ceed_data, problem)); 49791933550SJames Wright 498b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_jd_i)); 499b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_ijacobian_vol)); 500b4c37c5cSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&jac_data)); 501d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 502a515125bSLeila Ghaffari } 503