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