1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 3c864c5abSJames Wright 4149fb536SJames Wright #include <navierstokes.h> 5c864c5abSJames Wright 6c864c5abSJames Wright #include <petscsection.h> 7c864c5abSJames Wright #include "../qfunctions/setupgeo.h" 8c864c5abSJames Wright #include "../qfunctions/setupgeo2d.h" 9c864c5abSJames Wright 10e816a7e4SJames Wright static CeedVector *q_data_vecs; 11e816a7e4SJames Wright static CeedElemRestriction *q_data_restrictions; 12e816a7e4SJames Wright static PetscInt num_q_data_stored; 13e816a7e4SJames Wright 14e816a7e4SJames Wright /** 15e816a7e4SJames Wright @brief Get stored QData objects that match created objects, if present 16e816a7e4SJames Wright 17e816a7e4SJames Wright If created objects are not present, they are added to the storage and returned in the output 18e816a7e4SJames Wright 19e816a7e4SJames Wright Not Collective 20e816a7e4SJames Wright 21e816a7e4SJames Wright @param[in] q_data_created Vector with created QData 22e816a7e4SJames Wright @param[in] elem_restr_qd_created Restriction for created QData 23e816a7e4SJames Wright @param[out] q_data_stored Vector from storage matching QData 24e816a7e4SJames Wright @param[out] elem_restr_qd_stored Restriction from storage matching QData 25e816a7e4SJames Wright **/ 26e816a7e4SJames Wright static PetscErrorCode QDataGetStored(CeedVector q_data_created, CeedElemRestriction elem_restr_qd_created, CeedVector *q_data_stored, 27e816a7e4SJames Wright CeedElemRestriction *elem_restr_qd_stored) { 28e816a7e4SJames Wright Ceed ceed = CeedVectorReturnCeed(q_data_created); 29e816a7e4SJames Wright CeedSize created_length, stored_length; 30e816a7e4SJames Wright PetscInt q_data_stored_index = -1; 31e816a7e4SJames Wright 32e816a7e4SJames Wright PetscFunctionBeginUser; 33e816a7e4SJames Wright PetscCallCeed(ceed, CeedVectorGetLength(q_data_created, &created_length)); 34e816a7e4SJames Wright for (PetscInt i = 0; i < num_q_data_stored; i++) { 35e816a7e4SJames Wright CeedVector difference_cvec; 36e816a7e4SJames Wright CeedScalar max_difference; 37e816a7e4SJames Wright 38e816a7e4SJames Wright PetscCallCeed(ceed, CeedVectorGetLength(q_data_vecs[0], &stored_length)); 39e816a7e4SJames Wright if (created_length != stored_length) continue; 40e816a7e4SJames Wright PetscCallCeed(ceed, CeedVectorCreate(ceed, stored_length, &difference_cvec)); 41e816a7e4SJames Wright PetscCallCeed(ceed, CeedVectorCopy(q_data_vecs[i], difference_cvec)); 42e816a7e4SJames Wright PetscCallCeed(ceed, CeedVectorAXPY(difference_cvec, -1, q_data_created)); 43e816a7e4SJames Wright PetscCallCeed(ceed, CeedVectorNorm(difference_cvec, CEED_NORM_MAX, &max_difference)); 44be29160dSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&difference_cvec)); 45e816a7e4SJames Wright if (max_difference < 100 * CEED_EPSILON) { 46e816a7e4SJames Wright q_data_stored_index = i; 47e816a7e4SJames Wright break; 48e816a7e4SJames Wright } 49e816a7e4SJames Wright } 50e816a7e4SJames Wright 51e816a7e4SJames Wright if (q_data_stored_index == -1) { 52e816a7e4SJames Wright q_data_stored_index = num_q_data_stored++; 53e816a7e4SJames Wright 54e816a7e4SJames Wright PetscCall(PetscRealloc(num_q_data_stored * sizeof(CeedVector), &q_data_vecs)); 55e816a7e4SJames Wright PetscCall(PetscRealloc(num_q_data_stored * sizeof(CeedElemRestriction), &q_data_restrictions)); 56e816a7e4SJames Wright q_data_vecs[q_data_stored_index] = NULL; // Must set to NULL for ReferenceCopy 57e816a7e4SJames Wright q_data_restrictions[q_data_stored_index] = NULL; // Must set to NULL for ReferenceCopy 58e816a7e4SJames Wright PetscCallCeed(ceed, CeedVectorReferenceCopy(q_data_created, &q_data_vecs[q_data_stored_index])); 59e816a7e4SJames Wright PetscCallCeed(ceed, CeedElemRestrictionReferenceCopy(elem_restr_qd_created, &q_data_restrictions[q_data_stored_index])); 60e816a7e4SJames Wright } 61e816a7e4SJames Wright *q_data_stored = NULL; // Must set to NULL for ReferenceCopy 62e816a7e4SJames Wright *elem_restr_qd_stored = NULL; // Must set to NULL for ReferenceCopy 63e816a7e4SJames Wright PetscCallCeed(ceed, CeedVectorReferenceCopy(q_data_vecs[q_data_stored_index], q_data_stored)); 64e816a7e4SJames Wright PetscCallCeed(ceed, CeedElemRestrictionReferenceCopy(q_data_restrictions[q_data_stored_index], elem_restr_qd_stored)); 65e816a7e4SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 66e816a7e4SJames Wright } 67e816a7e4SJames Wright 68e816a7e4SJames Wright /** 69e816a7e4SJames Wright @brief Clear stored QData objects 70e816a7e4SJames Wright **/ 71e816a7e4SJames Wright PetscErrorCode QDataClearStoredData() { 72e816a7e4SJames Wright PetscFunctionBeginUser; 73e816a7e4SJames Wright for (PetscInt i = 0; i < num_q_data_stored; i++) { 74e816a7e4SJames Wright Ceed ceed = CeedVectorReturnCeed(q_data_vecs[i]); 75e816a7e4SJames Wright 76e816a7e4SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data_vecs[i])); 77e816a7e4SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&q_data_restrictions[i])); 78e816a7e4SJames Wright } 79e816a7e4SJames Wright PetscCall(PetscFree(q_data_vecs)); 80e816a7e4SJames Wright PetscCall(PetscFree(q_data_restrictions)); 81e816a7e4SJames Wright q_data_vecs = NULL; 82e816a7e4SJames Wright q_data_restrictions = NULL; 83e816a7e4SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 84e816a7e4SJames Wright } 85e816a7e4SJames Wright 86c864c5abSJames Wright /** 87c864c5abSJames Wright * @brief Get number of components of quadrature data for domain 88c864c5abSJames Wright * 89c864c5abSJames Wright * @param[in] dm DM where quadrature data would be used 90c864c5abSJames Wright * @param[out] q_data_size Number of components of quadrature data 91c864c5abSJames Wright */ 92c864c5abSJames Wright PetscErrorCode QDataGetNumComponents(DM dm, CeedInt *q_data_size) { 93c864c5abSJames Wright PetscInt num_comp_x, dim; 94c864c5abSJames Wright 95c864c5abSJames Wright PetscFunctionBeginUser; 96c864c5abSJames Wright PetscCall(DMGetDimension(dm, &dim)); 97c864c5abSJames Wright { // Get number of coordinate components 98c864c5abSJames Wright DM dm_coord; 99c864c5abSJames Wright PetscSection section_coord; 100c864c5abSJames Wright PetscInt field = 0; // Default field has the coordinates 1014aea4664SJames Wright 102c864c5abSJames Wright PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 103c864c5abSJames Wright PetscCall(DMGetLocalSection(dm_coord, §ion_coord)); 104c864c5abSJames Wright PetscCall(PetscSectionGetFieldComponents(section_coord, field, &num_comp_x)); 105c864c5abSJames Wright } 106*da8b59d6SJames Wright *q_data_size = 0; 107c864c5abSJames Wright switch (dim) { 108c864c5abSJames Wright case 2: 109c864c5abSJames Wright switch (num_comp_x) { 110c864c5abSJames Wright case 2: 111c864c5abSJames Wright *q_data_size = 5; 112c864c5abSJames Wright break; 113c864c5abSJames Wright case 3: 114c864c5abSJames Wright *q_data_size = 7; 115c864c5abSJames Wright break; 116c864c5abSJames Wright } 117c864c5abSJames Wright break; 118c864c5abSJames Wright case 3: 119c864c5abSJames Wright *q_data_size = 10; 120c864c5abSJames Wright break; 121c864c5abSJames Wright } 122*da8b59d6SJames Wright PetscCheck(*q_data_size, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, 123*da8b59d6SJames Wright "QData not valid for DM of dimension %" PetscInt_FMT " and coordinates with dimension %" PetscInt_FMT, dim, num_comp_x); 124c864c5abSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 125c864c5abSJames Wright } 126c864c5abSJames Wright 127c864c5abSJames Wright /** 128c864c5abSJames Wright * @brief Create quadrature data for domain 129c864c5abSJames Wright * 130c864c5abSJames Wright * @param[in] ceed Ceed object quadrature data will be used with 131c864c5abSJames Wright * @param[in] dm DM where quadrature data would be used 132c864c5abSJames Wright * @param[in] domain_label DMLabel that quadrature data would be used one 133c864c5abSJames Wright * @param[in] label_value Value of label 134c864c5abSJames Wright * @param[in] elem_restr_x CeedElemRestriction of the coordinates (must match `domain_label` and `label_value` selections) 135c864c5abSJames Wright * @param[in] basis_x CeedBasis of the coordinates 136c864c5abSJames Wright * @param[in] x_coord CeedVector of the coordinates 137c864c5abSJames Wright * @param[out] elem_restr_qd CeedElemRestriction of the quadrature data 138c864c5abSJames Wright * @param[out] q_data CeedVector of the quadrature data 1394aea4664SJames Wright * @param[out] q_data_size Number of components of quadrature data 140c864c5abSJames Wright */ 141c864c5abSJames Wright PetscErrorCode QDataGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedElemRestriction elem_restr_x, CeedBasis basis_x, 142c864c5abSJames Wright CeedVector x_coord, CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size) { 14387d3884fSJeremy L Thompson CeedQFunction qf_setup = NULL; 144c864c5abSJames Wright CeedOperator op_setup; 145e816a7e4SJames Wright CeedVector q_data_created; 146e816a7e4SJames Wright CeedElemRestriction elem_restr_qd_created; 147c864c5abSJames Wright CeedInt num_comp_x; 148c864c5abSJames Wright PetscInt dim, height = 0; 149c864c5abSJames Wright 150c864c5abSJames Wright PetscFunctionBeginUser; 151c864c5abSJames Wright PetscCall(QDataGetNumComponents(dm, q_data_size)); 152c864c5abSJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_x, &num_comp_x)); 153c864c5abSJames Wright PetscCall(DMGetDimension(dm, &dim)); 154c864c5abSJames Wright switch (dim) { 155c864c5abSJames Wright case 2: 156c864c5abSJames Wright switch (num_comp_x) { 157c864c5abSJames Wright case 2: 158c864c5abSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Setup2d, Setup2d_loc, &qf_setup)); 159c864c5abSJames Wright break; 160c864c5abSJames Wright case 3: 161c864c5abSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Setup2D_3Dcoords, Setup2D_3Dcoords_loc, &qf_setup)); 162c864c5abSJames Wright break; 163c864c5abSJames Wright } 164c864c5abSJames Wright break; 165c864c5abSJames Wright case 3: 166c864c5abSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Setup, Setup_loc, &qf_setup)); 167c864c5abSJames Wright break; 168c864c5abSJames Wright } 169*da8b59d6SJames Wright PetscCheck(qf_setup, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "%s not valid for DM of dimension %" PetscInt_FMT, __func__, dim); 170c864c5abSJames Wright 171c864c5abSJames Wright // -- Create QFunction for quadrature data 172c864c5abSJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_setup, 0)); 173c864c5abSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup, "dx", num_comp_x * (dim - height), CEED_EVAL_GRAD)); 174c864c5abSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup, "weight", 1, CEED_EVAL_WEIGHT)); 175c864c5abSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_setup, "surface qdata", *q_data_size, CEED_EVAL_NONE)); 176c864c5abSJames Wright 177e816a7e4SJames Wright PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, *q_data_size, &elem_restr_qd_created)); 178e816a7e4SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_qd_created, &q_data_created, NULL)); 179c864c5abSJames Wright 180c864c5abSJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_setup, NULL, NULL, &op_setup)); 181c864c5abSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "dx", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE)); 182c864c5abSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "weight", CEED_ELEMRESTRICTION_NONE, basis_x, CEED_VECTOR_NONE)); 183e816a7e4SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "surface qdata", elem_restr_qd_created, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 184c864c5abSJames Wright 185e816a7e4SJames Wright PetscCallCeed(ceed, CeedOperatorApply(op_setup, x_coord, q_data_created, CEED_REQUEST_IMMEDIATE)); 186c864c5abSJames Wright 187e816a7e4SJames Wright PetscCall(QDataGetStored(q_data_created, elem_restr_qd_created, q_data, elem_restr_qd)); 188e816a7e4SJames Wright 189e816a7e4SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data_created)); 190e816a7e4SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd_created)); 191c864c5abSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup)); 192c864c5abSJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_setup)); 193c864c5abSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 194c864c5abSJames Wright } 195c864c5abSJames Wright 196c864c5abSJames Wright /** 197c864c5abSJames Wright * @brief Get number of components of quadrature data for boundary of domain 198c864c5abSJames Wright * 199c864c5abSJames Wright * @param[in] dm DM where quadrature data would be used 200c864c5abSJames Wright * @param[out] q_data_size Number of components of quadrature data 201c864c5abSJames Wright */ 202c864c5abSJames Wright PetscErrorCode QDataBoundaryGetNumComponents(DM dm, CeedInt *q_data_size) { 203c864c5abSJames Wright PetscInt dim; 204c864c5abSJames Wright 205c864c5abSJames Wright PetscFunctionBeginUser; 206c864c5abSJames Wright PetscCall(DMGetDimension(dm, &dim)); 207*da8b59d6SJames Wright *q_data_size = 0; 208c864c5abSJames Wright switch (dim) { 209c864c5abSJames Wright case 2: 210c864c5abSJames Wright *q_data_size = 3; 211c864c5abSJames Wright break; 212c864c5abSJames Wright case 3: 213c864c5abSJames Wright *q_data_size = 10; 214c864c5abSJames Wright break; 215c864c5abSJames Wright } 216*da8b59d6SJames Wright PetscCheck(*q_data_size, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "%s not valid for DM of dimension %" PetscInt_FMT, __func__, dim); 217c864c5abSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 218c864c5abSJames Wright } 219c864c5abSJames Wright 220c864c5abSJames Wright /** 221c864c5abSJames Wright * @brief Create quadrature data for boundary of domain 222c864c5abSJames Wright * 223c864c5abSJames Wright * @param[in] ceed Ceed object quadrature data will be used with 224c864c5abSJames Wright * @param[in] dm DM where quadrature data would be used 225c864c5abSJames Wright * @param[in] domain_label DMLabel that quadrature data would be used one 226c864c5abSJames Wright * @param[in] label_value Value of label 227c864c5abSJames Wright * @param[in] elem_restr_x CeedElemRestriction of the coordinates (must match `domain_label` and `label_value` selections) 228c864c5abSJames Wright * @param[in] basis_x CeedBasis of the coordinates 229c864c5abSJames Wright * @param[in] x_coord CeedVector of the coordinates 230c864c5abSJames Wright * @param[out] elem_restr_qd CeedElemRestriction of the quadrature data 231c864c5abSJames Wright * @param[out] q_data CeedVector of the quadrature data 2324aea4664SJames Wright * @param[out] q_data_size Number of components of quadrature data 233c864c5abSJames Wright */ 234c864c5abSJames Wright PetscErrorCode QDataBoundaryGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedElemRestriction elem_restr_x, CeedBasis basis_x, 235c864c5abSJames Wright CeedVector x_coord, CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size) { 23687d3884fSJeremy L Thompson CeedQFunction qf_setup_sur = NULL; 237c864c5abSJames Wright CeedOperator op_setup_sur; 238e816a7e4SJames Wright CeedVector q_data_created; 239e816a7e4SJames Wright CeedElemRestriction elem_restr_qd_created; 240c864c5abSJames Wright CeedInt num_comp_x; 241c864c5abSJames Wright PetscInt dim, height = 1; 242c864c5abSJames Wright 243c864c5abSJames Wright PetscFunctionBeginUser; 244c864c5abSJames Wright PetscCall(QDataBoundaryGetNumComponents(dm, q_data_size)); 245c864c5abSJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_x, &num_comp_x)); 246c864c5abSJames Wright PetscCall(DMGetDimension(dm, &dim)); 247c864c5abSJames Wright switch (dim) { 248c864c5abSJames Wright case 2: 249c864c5abSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, SetupBoundary2d, SetupBoundary2d_loc, &qf_setup_sur)); 250c864c5abSJames Wright break; 251c864c5abSJames Wright case 3: 252c864c5abSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, SetupBoundary, SetupBoundary_loc, &qf_setup_sur)); 253c864c5abSJames Wright break; 254c864c5abSJames Wright } 255*da8b59d6SJames Wright PetscCheck(qf_setup_sur, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "%s not valid for DM of dimension %" PetscInt_FMT, __func__, dim); 256c864c5abSJames Wright 257c864c5abSJames Wright // -- Create QFunction for quadrature data 258c864c5abSJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_setup_sur, 0)); 259c864c5abSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup_sur, "dx", num_comp_x * (dim - height), CEED_EVAL_GRAD)); 260c864c5abSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup_sur, "weight", 1, CEED_EVAL_WEIGHT)); 261c864c5abSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_setup_sur, "surface qdata", *q_data_size, CEED_EVAL_NONE)); 262c864c5abSJames Wright 263e816a7e4SJames Wright PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, *q_data_size, &elem_restr_qd_created)); 264e816a7e4SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_qd_created, &q_data_created, NULL)); 265c864c5abSJames Wright 266c864c5abSJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_setup_sur, NULL, NULL, &op_setup_sur)); 267c864c5abSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "dx", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE)); 268c864c5abSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, basis_x, CEED_VECTOR_NONE)); 269e816a7e4SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "surface qdata", elem_restr_qd_created, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 270c864c5abSJames Wright 271e816a7e4SJames Wright PetscCallCeed(ceed, CeedOperatorApply(op_setup_sur, x_coord, q_data_created, CEED_REQUEST_IMMEDIATE)); 272c864c5abSJames Wright 273e816a7e4SJames Wright PetscCall(QDataGetStored(q_data_created, elem_restr_qd_created, q_data, elem_restr_qd)); 274e816a7e4SJames Wright 275e816a7e4SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data_created)); 276e816a7e4SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd_created)); 277c864c5abSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup_sur)); 278c864c5abSJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_setup_sur)); 279c864c5abSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 280c864c5abSJames Wright } 2818c85b835SJames Wright 2828c85b835SJames Wright /** 2838c85b835SJames Wright * @brief Get number of components of quadrature data for boundary of domain 2848c85b835SJames Wright * 2858c85b835SJames Wright * @param[in] dm DM where quadrature data would be used 2868c85b835SJames Wright * @param[out] q_data_size Number of components of quadrature data 2878c85b835SJames Wright */ 2888c85b835SJames Wright PetscErrorCode QDataBoundaryGradientGetNumComponents(DM dm, CeedInt *q_data_size) { 2898c85b835SJames Wright PetscInt dim; 2908c85b835SJames Wright 2918c85b835SJames Wright PetscFunctionBeginUser; 2928c85b835SJames Wright PetscCall(DMGetDimension(dm, &dim)); 293*da8b59d6SJames Wright *q_data_size = 0; 2948c85b835SJames Wright switch (dim) { 295*da8b59d6SJames Wright case 2: 296*da8b59d6SJames Wright *q_data_size = 7; 297*da8b59d6SJames Wright break; 2988c85b835SJames Wright case 3: 2998c85b835SJames Wright *q_data_size = 13; 3008c85b835SJames Wright break; 3018c85b835SJames Wright } 302*da8b59d6SJames Wright PetscCheck(*q_data_size, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "%s not valid for DM of dimension %" PetscInt_FMT, __func__, dim); 3038c85b835SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3048c85b835SJames Wright } 3058c85b835SJames Wright 3068c85b835SJames Wright /** 3078c85b835SJames Wright @brief Compute `CeedOperator` surface gradient QData 3088c85b835SJames Wright 3098c85b835SJames Wright Collective across MPI processes. 3108c85b835SJames Wright 3118c85b835SJames Wright @param[in] ceed `Ceed` object 3128c85b835SJames Wright @param[in] dm `DMPlex` grid 3138c85b835SJames Wright @param[in] domain_label `DMLabel` for surface 3148c85b835SJames Wright @param[in] label_value `DMPlex` label value for surface 3158c85b835SJames Wright @param[in] x_coord `CeedVector` for coordinates 31600dbc7b1SJames Wright @param[out] elem_restr_qd `CeedElemRestriction` for QData 3178c85b835SJames Wright @param[out] q_data `CeedVector` holding QData 3188c85b835SJames Wright @param[out] q_data_size Number of QData components per quadrature point 3198c85b835SJames Wright **/ 3208c85b835SJames Wright PetscErrorCode QDataBoundaryGradientGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedVector x_coord, 32100dbc7b1SJames Wright CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size) { 3228c85b835SJames Wright PetscInt dim; 3238c85b835SJames Wright const PetscInt height_cell = 0, height_face = 1; 3248c85b835SJames Wright CeedInt num_comp_x; 325f17df9b6SJames Wright CeedElemRestriction elem_restr_x_cell, elem_restr_x_face, elem_restr_qd_created; 326f17df9b6SJames Wright CeedVector q_data_created; 3278c85b835SJames Wright CeedBasis basis_x_cell_to_face, basis_x_face; 328*da8b59d6SJames Wright CeedQFunction qf_setup_sur = NULL; 3298c85b835SJames Wright CeedOperator op_setup_sur; 3308c85b835SJames Wright 3318c85b835SJames Wright PetscFunctionBeginUser; 3328c85b835SJames Wright PetscCall(CreateCoordinateBasisFromPlex(ceed, dm, domain_label, label_value, height_face, &basis_x_face)); 3338c85b835SJames Wright PetscCall(DMPlexCeedBasisCellToFaceCoordinateCreate(ceed, dm, domain_label, label_value, label_value, &basis_x_cell_to_face)); 33400dbc7b1SJames Wright PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height_cell, &elem_restr_x_cell)); 3358c85b835SJames Wright PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height_face, &elem_restr_x_face)); 3368c85b835SJames Wright 3378c85b835SJames Wright PetscCall(QDataBoundaryGradientGetNumComponents(dm, q_data_size)); 338f17df9b6SJames Wright PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height_face, *q_data_size, &elem_restr_qd_created)); 339f17df9b6SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_qd_created, &q_data_created, NULL)); 3408c85b835SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_x_face, &num_comp_x)); 3418c85b835SJames Wright 342*da8b59d6SJames Wright PetscCall(DMGetDimension(dm, &dim)); 343*da8b59d6SJames Wright switch (dim) { 344*da8b59d6SJames Wright case 2: 345*da8b59d6SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Setup2DBoundaryGradient, Setup2DBoundaryGradient_loc, &qf_setup_sur)); 346*da8b59d6SJames Wright break; 347*da8b59d6SJames Wright case 3: 3488c85b835SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, SetupBoundaryGradient, SetupBoundaryGradient_loc, &qf_setup_sur)); 349*da8b59d6SJames Wright break; 350*da8b59d6SJames Wright } 351*da8b59d6SJames Wright PetscCheck(qf_setup_sur, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "%s not valid for DM of dimension %" PetscInt_FMT, __func__, dim); 352*da8b59d6SJames Wright 3538c85b835SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup_sur, "dx/dX cell", num_comp_x * (dim - height_cell), CEED_EVAL_GRAD)); 3548c85b835SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup_sur, "dx/dX face", num_comp_x * (dim - height_face), CEED_EVAL_GRAD)); 3558c85b835SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup_sur, "weight", 1, CEED_EVAL_WEIGHT)); 3568c85b835SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_setup_sur, "q data", *q_data_size, CEED_EVAL_NONE)); 3578c85b835SJames Wright 3588c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_setup_sur, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_setup_sur)); 35900dbc7b1SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "dx/dX cell", elem_restr_x_cell, basis_x_cell_to_face, CEED_VECTOR_ACTIVE)); 3608c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "dx/dX face", elem_restr_x_face, basis_x_face, CEED_VECTOR_ACTIVE)); 3618c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, basis_x_face, CEED_VECTOR_NONE)); 362f17df9b6SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "q data", elem_restr_qd_created, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 3638c85b835SJames Wright 364f17df9b6SJames Wright PetscCallCeed(ceed, CeedOperatorApply(op_setup_sur, x_coord, q_data_created, CEED_REQUEST_IMMEDIATE)); 365f17df9b6SJames Wright 366f17df9b6SJames Wright PetscCall(QDataGetStored(q_data_created, elem_restr_qd_created, q_data, elem_restr_qd)); 3678c85b835SJames Wright 3688c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup_sur)); 3698c85b835SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_setup_sur)); 370f17df9b6SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data_created)); 371f17df9b6SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd_created)); 3728c85b835SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x_face)); 37300dbc7b1SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x_cell)); 3748c85b835SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_face)); 3758c85b835SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_cell_to_face)); 3768c85b835SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3778c85b835SJames Wright } 378