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