1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, Lawrence Livermore National Security, LLC and other CEED contributors. 2cfb075a4SJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3cfb075a4SJames Wright // 4cfb075a4SJames Wright // SPDX-License-Identifier: BSD-2-Clause 5cfb075a4SJames Wright // 6cfb075a4SJames Wright // This file is part of CEED: http://github.com/ceed 7cfb075a4SJames Wright 8cfb075a4SJames Wright #include "../navierstokes.h" 9cfb075a4SJames Wright 10cfb075a4SJames Wright #include <petscsection.h> 11cfb075a4SJames Wright #include "../qfunctions/setupgeo.h" 12cfb075a4SJames Wright #include "../qfunctions/setupgeo2d.h" 13cfb075a4SJames Wright 14cfb075a4SJames Wright /** 15cfb075a4SJames Wright * @brief Get number of components of quadrature data for domain 16cfb075a4SJames Wright * 17cfb075a4SJames Wright * @param[in] dm DM where quadrature data would be used 18cfb075a4SJames Wright * @param[out] q_data_size Number of components of quadrature data 19cfb075a4SJames Wright */ 20cfb075a4SJames Wright PetscErrorCode QDataGetNumComponents(DM dm, CeedInt *q_data_size) { 21cfb075a4SJames Wright PetscInt num_comp_x, dim; 22cfb075a4SJames Wright 23cfb075a4SJames Wright PetscFunctionBeginUser; 24cfb075a4SJames Wright PetscCall(DMGetDimension(dm, &dim)); 25cfb075a4SJames Wright { // Get number of coordinate components 26cfb075a4SJames Wright DM dm_coord; 27cfb075a4SJames Wright PetscSection section_coord; 28cfb075a4SJames Wright PetscInt field = 0; // Default field has the coordinates 29cfb075a4SJames Wright PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 30cfb075a4SJames Wright PetscCall(DMGetLocalSection(dm_coord, §ion_coord)); 31cfb075a4SJames Wright PetscCall(PetscSectionGetFieldComponents(section_coord, field, &num_comp_x)); 32cfb075a4SJames Wright } 33cfb075a4SJames Wright switch (dim) { 34cfb075a4SJames Wright case 2: 35cfb075a4SJames Wright switch (num_comp_x) { 36cfb075a4SJames Wright case 2: 37cfb075a4SJames Wright *q_data_size = 5; 38cfb075a4SJames Wright break; 39cfb075a4SJames Wright case 3: 40cfb075a4SJames Wright *q_data_size = 7; 41cfb075a4SJames Wright break; 42cfb075a4SJames Wright default: 43cfb075a4SJames Wright SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, 44cfb075a4SJames Wright "QData not valid for DM of dimension %" PetscInt_FMT " and coordinates with dimension %" PetscInt_FMT, dim, num_comp_x); 45cfb075a4SJames Wright break; 46cfb075a4SJames Wright } 47cfb075a4SJames Wright break; 48cfb075a4SJames Wright case 3: 49cfb075a4SJames Wright *q_data_size = 10; 50cfb075a4SJames Wright break; 51cfb075a4SJames Wright default: 52cfb075a4SJames Wright SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, 53cfb075a4SJames Wright "QData not valid for DM of dimension %" PetscInt_FMT " and coordinates with dimension %" PetscInt_FMT, dim, num_comp_x); 54cfb075a4SJames Wright break; 55cfb075a4SJames Wright } 56cfb075a4SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 57cfb075a4SJames Wright } 58cfb075a4SJames Wright 59cfb075a4SJames Wright /** 60cfb075a4SJames Wright * @brief Create quadrature data for domain 61cfb075a4SJames Wright * 62cfb075a4SJames Wright * @param[in] ceed Ceed object quadrature data will be used with 63cfb075a4SJames Wright * @param[in] dm DM where quadrature data would be used 64cfb075a4SJames Wright * @param[in] domain_label DMLabel that quadrature data would be used one 65cfb075a4SJames Wright * @param[in] label_value Value of label 66cfb075a4SJames Wright * @param[in] elem_restr_x CeedElemRestriction of the coordinates (must match `domain_label` and `label_value` selections) 67cfb075a4SJames Wright * @param[in] basis_x CeedBasis of the coordinates 68cfb075a4SJames Wright * @param[in] x_coord CeedVector of the coordinates 69cfb075a4SJames Wright * @param[out] elem_restr_qd CeedElemRestriction of the quadrature data 70cfb075a4SJames Wright * @param[out] q_data CeedVector of the quadrature data 71cfb075a4SJames Wright * @param[out] q_data_size number of components of quadrature data 72cfb075a4SJames Wright */ 73cfb075a4SJames Wright PetscErrorCode QDataGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedElemRestriction elem_restr_x, CeedBasis basis_x, 74cfb075a4SJames Wright CeedVector x_coord, CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size) { 75cfb075a4SJames Wright CeedQFunction qf_setup; 76cfb075a4SJames Wright CeedOperator op_setup; 77cfb075a4SJames Wright CeedInt num_comp_x; 78cfb075a4SJames Wright PetscInt dim, height = 0; 79cfb075a4SJames Wright 80cfb075a4SJames Wright PetscFunctionBeginUser; 81cfb075a4SJames Wright PetscCall(QDataGetNumComponents(dm, q_data_size)); 82cfb075a4SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_x, &num_comp_x)); 83cfb075a4SJames Wright PetscCall(DMGetDimension(dm, &dim)); 84cfb075a4SJames Wright switch (dim) { 85cfb075a4SJames Wright case 2: 86cfb075a4SJames Wright switch (num_comp_x) { 87cfb075a4SJames Wright case 2: 88cfb075a4SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Setup2d, Setup2d_loc, &qf_setup)); 89cfb075a4SJames Wright break; 90cfb075a4SJames Wright case 3: 91cfb075a4SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Setup2D_3Dcoords, Setup2D_3Dcoords_loc, &qf_setup)); 92cfb075a4SJames Wright break; 93cfb075a4SJames Wright } 94cfb075a4SJames Wright break; 95cfb075a4SJames Wright case 3: 96cfb075a4SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Setup, Setup_loc, &qf_setup)); 97cfb075a4SJames Wright break; 98cfb075a4SJames Wright } 99cfb075a4SJames Wright 100cfb075a4SJames Wright // -- Create QFunction for quadrature data 101cfb075a4SJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_setup, 0)); 102cfb075a4SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup, "dx", num_comp_x * (dim - height), CEED_EVAL_GRAD)); 103cfb075a4SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup, "weight", 1, CEED_EVAL_WEIGHT)); 104cfb075a4SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_setup, "surface qdata", *q_data_size, CEED_EVAL_NONE)); 105cfb075a4SJames Wright 106cfb075a4SJames Wright PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, *q_data_size, elem_restr_qd)); 107cfb075a4SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(*elem_restr_qd, q_data, NULL)); 108cfb075a4SJames Wright 109cfb075a4SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_setup, NULL, NULL, &op_setup)); 110cfb075a4SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "dx", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE)); 111cfb075a4SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "weight", CEED_ELEMRESTRICTION_NONE, basis_x, CEED_VECTOR_NONE)); 112cfb075a4SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "surface qdata", *elem_restr_qd, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 113cfb075a4SJames Wright 114cfb075a4SJames Wright PetscCallCeed(ceed, CeedOperatorApply(op_setup, x_coord, *q_data, CEED_REQUEST_IMMEDIATE)); 115cfb075a4SJames Wright 116cfb075a4SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup)); 117cfb075a4SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_setup)); 118cfb075a4SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 119cfb075a4SJames Wright } 120cfb075a4SJames Wright 121cfb075a4SJames Wright /** 122cfb075a4SJames Wright * @brief Get number of components of quadrature data for boundary of domain 123cfb075a4SJames Wright * 124cfb075a4SJames Wright * @param[in] dm DM where quadrature data would be used 125cfb075a4SJames Wright * @param[out] q_data_size Number of components of quadrature data 126cfb075a4SJames Wright */ 127cfb075a4SJames Wright PetscErrorCode QDataBoundaryGetNumComponents(DM dm, CeedInt *q_data_size) { 128cfb075a4SJames Wright PetscInt dim; 129cfb075a4SJames Wright 130cfb075a4SJames Wright PetscFunctionBeginUser; 131cfb075a4SJames Wright PetscCall(DMGetDimension(dm, &dim)); 132cfb075a4SJames Wright switch (dim) { 133cfb075a4SJames Wright case 2: 134cfb075a4SJames Wright *q_data_size = 3; 135cfb075a4SJames Wright break; 136cfb075a4SJames Wright case 3: 137cfb075a4SJames Wright *q_data_size = 10; 138cfb075a4SJames Wright break; 139cfb075a4SJames Wright default: 140cfb075a4SJames Wright SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "QDataBoundary not valid for DM of dimension %" PetscInt_FMT, dim); 141cfb075a4SJames Wright break; 142cfb075a4SJames Wright } 143cfb075a4SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 144cfb075a4SJames Wright } 145cfb075a4SJames Wright 146cfb075a4SJames Wright /** 147cfb075a4SJames Wright * @brief Create quadrature data for boundary of domain 148cfb075a4SJames Wright * 149cfb075a4SJames Wright * @param[in] ceed Ceed object quadrature data will be used with 150cfb075a4SJames Wright * @param[in] dm DM where quadrature data would be used 151cfb075a4SJames Wright * @param[in] domain_label DMLabel that quadrature data would be used one 152cfb075a4SJames Wright * @param[in] label_value Value of label 153cfb075a4SJames Wright * @param[in] elem_restr_x CeedElemRestriction of the coordinates (must match `domain_label` and `label_value` selections) 154cfb075a4SJames Wright * @param[in] basis_x CeedBasis of the coordinates 155cfb075a4SJames Wright * @param[in] x_coord CeedVector of the coordinates 156cfb075a4SJames Wright * @param[out] elem_restr_qd CeedElemRestriction of the quadrature data 157cfb075a4SJames Wright * @param[out] q_data CeedVector of the quadrature data 158cfb075a4SJames Wright * @param[out] q_data_size number of components of quadrature data 159cfb075a4SJames Wright */ 160cfb075a4SJames Wright PetscErrorCode QDataBoundaryGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedElemRestriction elem_restr_x, CeedBasis basis_x, 161cfb075a4SJames Wright CeedVector x_coord, CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size) { 162cfb075a4SJames Wright CeedQFunction qf_setup_sur; 163cfb075a4SJames Wright CeedOperator op_setup_sur; 164cfb075a4SJames Wright CeedInt num_comp_x; 165cfb075a4SJames Wright PetscInt dim, height = 1; 166cfb075a4SJames Wright 167cfb075a4SJames Wright PetscFunctionBeginUser; 168cfb075a4SJames Wright PetscCall(QDataBoundaryGetNumComponents(dm, q_data_size)); 169cfb075a4SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_x, &num_comp_x)); 170cfb075a4SJames Wright PetscCall(DMGetDimension(dm, &dim)); 171cfb075a4SJames Wright switch (dim) { 172cfb075a4SJames Wright case 2: 173cfb075a4SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, SetupBoundary2d, SetupBoundary2d_loc, &qf_setup_sur)); 174cfb075a4SJames Wright break; 175cfb075a4SJames Wright case 3: 176cfb075a4SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, SetupBoundary, SetupBoundary_loc, &qf_setup_sur)); 177cfb075a4SJames Wright break; 178cfb075a4SJames Wright } 179cfb075a4SJames Wright 180cfb075a4SJames Wright // -- Create QFunction for quadrature data 181cfb075a4SJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_setup_sur, 0)); 182cfb075a4SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup_sur, "dx", num_comp_x * (dim - height), CEED_EVAL_GRAD)); 183cfb075a4SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup_sur, "weight", 1, CEED_EVAL_WEIGHT)); 184cfb075a4SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_setup_sur, "surface qdata", *q_data_size, CEED_EVAL_NONE)); 185cfb075a4SJames Wright 186cfb075a4SJames Wright PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, *q_data_size, elem_restr_qd)); 187cfb075a4SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(*elem_restr_qd, q_data, NULL)); 188cfb075a4SJames Wright 189cfb075a4SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_setup_sur, NULL, NULL, &op_setup_sur)); 190cfb075a4SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "dx", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE)); 191cfb075a4SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, basis_x, CEED_VECTOR_NONE)); 192cfb075a4SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "surface qdata", *elem_restr_qd, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 193cfb075a4SJames Wright 194cfb075a4SJames Wright PetscCallCeed(ceed, CeedOperatorApply(op_setup_sur, x_coord, *q_data, CEED_REQUEST_IMMEDIATE)); 195cfb075a4SJames Wright 196cfb075a4SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup_sur)); 197cfb075a4SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_setup_sur)); 198cfb075a4SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 199cfb075a4SJames Wright } 200