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