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