xref: /honee/src/qdata.c (revision d033c862e3e2cd31b8487d554f45c0b0685a669f)
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));
44*d033c862SJames Wright     //TODO Need to reduce across ranks to ensure all ranks are consistent (not a race condition though since the data is purely local)
45be29160dSJames Wright     PetscCallCeed(ceed, CeedVectorDestroy(&difference_cvec));
46e816a7e4SJames Wright     if (max_difference < 100 * CEED_EPSILON) {
47e816a7e4SJames Wright       q_data_stored_index = i;
48e816a7e4SJames Wright       break;
49e816a7e4SJames Wright     }
50e816a7e4SJames Wright   }
51e816a7e4SJames Wright 
52e816a7e4SJames Wright   if (q_data_stored_index == -1) {
53e816a7e4SJames Wright     q_data_stored_index = num_q_data_stored++;
54e816a7e4SJames Wright 
55e816a7e4SJames Wright     PetscCall(PetscRealloc(num_q_data_stored * sizeof(CeedVector), &q_data_vecs));
56e816a7e4SJames Wright     PetscCall(PetscRealloc(num_q_data_stored * sizeof(CeedElemRestriction), &q_data_restrictions));
57e816a7e4SJames Wright     q_data_vecs[q_data_stored_index]         = NULL;  // Must set to NULL for ReferenceCopy
58e816a7e4SJames Wright     q_data_restrictions[q_data_stored_index] = NULL;  // Must set to NULL for ReferenceCopy
59e816a7e4SJames Wright     PetscCallCeed(ceed, CeedVectorReferenceCopy(q_data_created, &q_data_vecs[q_data_stored_index]));
60e816a7e4SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionReferenceCopy(elem_restr_qd_created, &q_data_restrictions[q_data_stored_index]));
61e816a7e4SJames Wright   }
62e816a7e4SJames Wright   *q_data_stored        = NULL;  // Must set to NULL for ReferenceCopy
63e816a7e4SJames Wright   *elem_restr_qd_stored = NULL;  // Must set to NULL for ReferenceCopy
64e816a7e4SJames Wright   PetscCallCeed(ceed, CeedVectorReferenceCopy(q_data_vecs[q_data_stored_index], q_data_stored));
65e816a7e4SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionReferenceCopy(q_data_restrictions[q_data_stored_index], elem_restr_qd_stored));
66e816a7e4SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
67e816a7e4SJames Wright }
68e816a7e4SJames Wright 
69e816a7e4SJames Wright /**
70e816a7e4SJames Wright   @brief Clear stored QData objects
71e816a7e4SJames Wright **/
72e816a7e4SJames Wright PetscErrorCode QDataClearStoredData() {
73e816a7e4SJames Wright   PetscFunctionBeginUser;
74e816a7e4SJames Wright   for (PetscInt i = 0; i < num_q_data_stored; i++) {
75e816a7e4SJames Wright     Ceed ceed = CeedVectorReturnCeed(q_data_vecs[i]);
76e816a7e4SJames Wright 
77e816a7e4SJames Wright     PetscCallCeed(ceed, CeedVectorDestroy(&q_data_vecs[i]));
78e816a7e4SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&q_data_restrictions[i]));
79e816a7e4SJames Wright   }
80e816a7e4SJames Wright   PetscCall(PetscFree(q_data_vecs));
81e816a7e4SJames Wright   PetscCall(PetscFree(q_data_restrictions));
82e816a7e4SJames Wright   q_data_vecs         = NULL;
83e816a7e4SJames Wright   q_data_restrictions = NULL;
84e816a7e4SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
85e816a7e4SJames Wright }
86e816a7e4SJames Wright 
87c864c5abSJames Wright /**
88c864c5abSJames Wright  * @brief Get number of components of quadrature data for domain
89c864c5abSJames Wright  *
90c864c5abSJames Wright  * @param[in]  dm          DM where quadrature data would be used
91c864c5abSJames Wright  * @param[out] q_data_size Number of components of quadrature data
92c864c5abSJames Wright  */
93c864c5abSJames Wright PetscErrorCode QDataGetNumComponents(DM dm, CeedInt *q_data_size) {
94c864c5abSJames Wright   PetscInt num_comp_x, dim;
95c864c5abSJames Wright 
96c864c5abSJames Wright   PetscFunctionBeginUser;
97c864c5abSJames Wright   PetscCall(DMGetDimension(dm, &dim));
98bc87537dSJames Wright   PetscCall(DMGetCoordinateNumComps(dm, &num_comp_x));
99da8b59d6SJames Wright   *q_data_size = 0;
100c864c5abSJames Wright   switch (dim) {
101c864c5abSJames Wright     case 2:
102c864c5abSJames Wright       switch (num_comp_x) {
103c864c5abSJames Wright         case 2:
104c864c5abSJames Wright           *q_data_size = 5;
105c864c5abSJames Wright           break;
106c864c5abSJames Wright         case 3:
107c864c5abSJames Wright           *q_data_size = 7;
108c864c5abSJames Wright           break;
109c864c5abSJames Wright       }
110c864c5abSJames Wright       break;
111c864c5abSJames Wright     case 3:
112c864c5abSJames Wright       *q_data_size = 10;
113c864c5abSJames Wright       break;
114c864c5abSJames Wright   }
115da8b59d6SJames Wright   PetscCheck(*q_data_size, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP,
116da8b59d6SJames Wright              "QData not valid for DM of dimension %" PetscInt_FMT " and coordinates with dimension %" PetscInt_FMT, dim, num_comp_x);
117c864c5abSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
118c864c5abSJames Wright }
119c864c5abSJames Wright 
120c864c5abSJames Wright /**
121c864c5abSJames Wright  * @brief Create quadrature data for domain
122c864c5abSJames Wright  *
123c864c5abSJames Wright  * @param[in]  ceed          Ceed object quadrature data will be used with
124c864c5abSJames Wright  * @param[in]  dm            DM where quadrature data would be used
125c864c5abSJames Wright  * @param[in]  domain_label  DMLabel that quadrature data would be used one
126c864c5abSJames Wright  * @param[in]  label_value   Value of label
127c864c5abSJames Wright  * @param[in]  elem_restr_x  CeedElemRestriction of the coordinates (must match `domain_label` and `label_value` selections)
128c864c5abSJames Wright  * @param[in]  basis_x       CeedBasis of the coordinates
129c864c5abSJames Wright  * @param[in]  x_coord       CeedVector of the coordinates
130c864c5abSJames Wright  * @param[out] elem_restr_qd CeedElemRestriction of the quadrature data
131c864c5abSJames Wright  * @param[out] q_data        CeedVector of the quadrature data
1324aea4664SJames Wright  * @param[out] q_data_size   Number of components of quadrature data
133c864c5abSJames Wright  */
134c864c5abSJames Wright PetscErrorCode QDataGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedElemRestriction elem_restr_x, CeedBasis basis_x,
135c864c5abSJames Wright                         CeedVector x_coord, CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size) {
13687d3884fSJeremy L Thompson   CeedQFunction       qf_setup = NULL;
137c864c5abSJames Wright   CeedOperator        op_setup;
138e816a7e4SJames Wright   CeedVector          q_data_created;
139e816a7e4SJames Wright   CeedElemRestriction elem_restr_qd_created;
140c864c5abSJames Wright   CeedInt             num_comp_x;
141c864c5abSJames Wright   PetscInt            dim, height = 0;
142c864c5abSJames Wright 
143c864c5abSJames Wright   PetscFunctionBeginUser;
144c864c5abSJames Wright   PetscCall(QDataGetNumComponents(dm, q_data_size));
145c864c5abSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_x, &num_comp_x));
146c864c5abSJames Wright   PetscCall(DMGetDimension(dm, &dim));
147c864c5abSJames Wright   switch (dim) {
148c864c5abSJames Wright     case 2:
149c864c5abSJames Wright       switch (num_comp_x) {
150c864c5abSJames Wright         case 2:
151c864c5abSJames Wright           PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Setup2d, Setup2d_loc, &qf_setup));
152c864c5abSJames Wright           break;
153c864c5abSJames Wright         case 3:
154c864c5abSJames Wright           PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Setup2D_3Dcoords, Setup2D_3Dcoords_loc, &qf_setup));
155c864c5abSJames Wright           break;
156c864c5abSJames Wright       }
157c864c5abSJames Wright       break;
158c864c5abSJames Wright     case 3:
159c864c5abSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Setup, Setup_loc, &qf_setup));
160c864c5abSJames Wright       break;
161c864c5abSJames Wright   }
162da8b59d6SJames Wright   PetscCheck(qf_setup, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "%s not valid for DM of dimension %" PetscInt_FMT, __func__, dim);
163c864c5abSJames Wright 
164c864c5abSJames Wright   // -- Create QFunction for quadrature data
165c864c5abSJames Wright   PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_setup, 0));
166c864c5abSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup, "dx", num_comp_x * (dim - height), CEED_EVAL_GRAD));
167c864c5abSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup, "weight", 1, CEED_EVAL_WEIGHT));
168c864c5abSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_setup, "surface qdata", *q_data_size, CEED_EVAL_NONE));
169c864c5abSJames Wright 
170e816a7e4SJames Wright   PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, *q_data_size, &elem_restr_qd_created));
171e816a7e4SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_qd_created, &q_data_created, NULL));
172c864c5abSJames Wright 
173c864c5abSJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_setup, NULL, NULL, &op_setup));
174c864c5abSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "dx", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE));
175c864c5abSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "weight", CEED_ELEMRESTRICTION_NONE, basis_x, CEED_VECTOR_NONE));
176e816a7e4SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "surface qdata", elem_restr_qd_created, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
177c864c5abSJames Wright 
178e816a7e4SJames Wright   PetscCallCeed(ceed, CeedOperatorApply(op_setup, x_coord, q_data_created, CEED_REQUEST_IMMEDIATE));
179c864c5abSJames Wright 
180e816a7e4SJames Wright   PetscCall(QDataGetStored(q_data_created, elem_restr_qd_created, q_data, elem_restr_qd));
181e816a7e4SJames Wright 
182e816a7e4SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&q_data_created));
183e816a7e4SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd_created));
184c864c5abSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup));
185c864c5abSJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_setup));
186c864c5abSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
187c864c5abSJames Wright }
188c864c5abSJames Wright 
189c864c5abSJames Wright /**
190c864c5abSJames Wright  * @brief Get number of components of quadrature data for boundary of domain
191c864c5abSJames Wright  *
192c864c5abSJames Wright  * @param[in]  dm          DM where quadrature data would be used
193c864c5abSJames Wright  * @param[out] q_data_size Number of components of quadrature data
194c864c5abSJames Wright  */
195c864c5abSJames Wright PetscErrorCode QDataBoundaryGetNumComponents(DM dm, CeedInt *q_data_size) {
196c864c5abSJames Wright   PetscInt dim;
197c864c5abSJames Wright 
198c864c5abSJames Wright   PetscFunctionBeginUser;
199c864c5abSJames Wright   PetscCall(DMGetDimension(dm, &dim));
200da8b59d6SJames Wright   *q_data_size = 0;
201c864c5abSJames Wright   switch (dim) {
202c864c5abSJames Wright     case 2:
203c864c5abSJames Wright       *q_data_size = 3;
204c864c5abSJames Wright       break;
205c864c5abSJames Wright     case 3:
206c864c5abSJames Wright       *q_data_size = 10;
207c864c5abSJames Wright       break;
208c864c5abSJames Wright   }
209da8b59d6SJames Wright   PetscCheck(*q_data_size, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "%s not valid for DM of dimension %" PetscInt_FMT, __func__, dim);
210c864c5abSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
211c864c5abSJames Wright }
212c864c5abSJames Wright 
213c864c5abSJames Wright /**
214c864c5abSJames Wright  * @brief Create quadrature data for boundary of domain
215c864c5abSJames Wright  *
216c864c5abSJames Wright  * @param[in]  ceed          Ceed object quadrature data will be used with
217c864c5abSJames Wright  * @param[in]  dm            DM where quadrature data would be used
218c864c5abSJames Wright  * @param[in]  domain_label  DMLabel that quadrature data would be used one
219c864c5abSJames Wright  * @param[in]  label_value   Value of label
220c864c5abSJames Wright  * @param[in]  elem_restr_x  CeedElemRestriction of the coordinates (must match `domain_label` and `label_value` selections)
221c864c5abSJames Wright  * @param[in]  basis_x       CeedBasis of the coordinates
222c864c5abSJames Wright  * @param[in]  x_coord       CeedVector of the coordinates
223c864c5abSJames Wright  * @param[out] elem_restr_qd CeedElemRestriction of the quadrature data
224c864c5abSJames Wright  * @param[out] q_data        CeedVector of the quadrature data
2254aea4664SJames Wright  * @param[out] q_data_size   Number of components of quadrature data
226c864c5abSJames Wright  */
227c864c5abSJames Wright PetscErrorCode QDataBoundaryGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedElemRestriction elem_restr_x, CeedBasis basis_x,
228c864c5abSJames Wright                                 CeedVector x_coord, CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size) {
22987d3884fSJeremy L Thompson   CeedQFunction       qf_setup_sur = NULL;
230c864c5abSJames Wright   CeedOperator        op_setup_sur;
231e816a7e4SJames Wright   CeedVector          q_data_created;
232e816a7e4SJames Wright   CeedElemRestriction elem_restr_qd_created;
233c864c5abSJames Wright   CeedInt             num_comp_x;
234c864c5abSJames Wright   PetscInt            dim, height = 1;
235c864c5abSJames Wright 
236c864c5abSJames Wright   PetscFunctionBeginUser;
237c864c5abSJames Wright   PetscCall(QDataBoundaryGetNumComponents(dm, q_data_size));
238c864c5abSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_x, &num_comp_x));
239c864c5abSJames Wright   PetscCall(DMGetDimension(dm, &dim));
240c864c5abSJames Wright   switch (dim) {
241c864c5abSJames Wright     case 2:
242c864c5abSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, SetupBoundary2d, SetupBoundary2d_loc, &qf_setup_sur));
243c864c5abSJames Wright       break;
244c864c5abSJames Wright     case 3:
245c864c5abSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, SetupBoundary, SetupBoundary_loc, &qf_setup_sur));
246c864c5abSJames Wright       break;
247c864c5abSJames Wright   }
248da8b59d6SJames Wright   PetscCheck(qf_setup_sur, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "%s not valid for DM of dimension %" PetscInt_FMT, __func__, dim);
249c864c5abSJames Wright 
250c864c5abSJames Wright   // -- Create QFunction for quadrature data
251c864c5abSJames Wright   PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_setup_sur, 0));
252c864c5abSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup_sur, "dx", num_comp_x * (dim - height), CEED_EVAL_GRAD));
253c864c5abSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup_sur, "weight", 1, CEED_EVAL_WEIGHT));
254c864c5abSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_setup_sur, "surface qdata", *q_data_size, CEED_EVAL_NONE));
255c864c5abSJames Wright 
256e816a7e4SJames Wright   PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, *q_data_size, &elem_restr_qd_created));
257e816a7e4SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_qd_created, &q_data_created, NULL));
258c864c5abSJames Wright 
259c864c5abSJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_setup_sur, NULL, NULL, &op_setup_sur));
260c864c5abSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "dx", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE));
261c864c5abSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, basis_x, CEED_VECTOR_NONE));
262e816a7e4SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "surface qdata", elem_restr_qd_created, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
263c864c5abSJames Wright 
264e816a7e4SJames Wright   PetscCallCeed(ceed, CeedOperatorApply(op_setup_sur, x_coord, q_data_created, CEED_REQUEST_IMMEDIATE));
265c864c5abSJames Wright 
266e816a7e4SJames Wright   PetscCall(QDataGetStored(q_data_created, elem_restr_qd_created, q_data, elem_restr_qd));
267e816a7e4SJames Wright 
268e816a7e4SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&q_data_created));
269e816a7e4SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd_created));
270c864c5abSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup_sur));
271c864c5abSJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_setup_sur));
272c864c5abSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
273c864c5abSJames Wright }
2748c85b835SJames Wright 
2758c85b835SJames Wright /**
2768c85b835SJames Wright  * @brief Get number of components of quadrature data for boundary of domain
2778c85b835SJames Wright  *
2788c85b835SJames Wright  * @param[in]  dm          DM where quadrature data would be used
2798c85b835SJames Wright  * @param[out] q_data_size Number of components of quadrature data
2808c85b835SJames Wright  */
2818c85b835SJames Wright PetscErrorCode QDataBoundaryGradientGetNumComponents(DM dm, CeedInt *q_data_size) {
2828c85b835SJames Wright   PetscInt dim;
2838c85b835SJames Wright 
2848c85b835SJames Wright   PetscFunctionBeginUser;
2858c85b835SJames Wright   PetscCall(DMGetDimension(dm, &dim));
286da8b59d6SJames Wright   *q_data_size = 0;
2878c85b835SJames Wright   switch (dim) {
288da8b59d6SJames Wright     case 2:
289da8b59d6SJames Wright       *q_data_size = 7;
290da8b59d6SJames Wright       break;
2918c85b835SJames Wright     case 3:
2928c85b835SJames Wright       *q_data_size = 13;
2938c85b835SJames Wright       break;
2948c85b835SJames Wright   }
295da8b59d6SJames Wright   PetscCheck(*q_data_size, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "%s not valid for DM of dimension %" PetscInt_FMT, __func__, dim);
2968c85b835SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2978c85b835SJames Wright }
2988c85b835SJames Wright 
2998c85b835SJames Wright /**
3008c85b835SJames Wright   @brief Compute `CeedOperator` surface gradient QData
3018c85b835SJames Wright 
3028c85b835SJames Wright   Collective across MPI processes.
3038c85b835SJames Wright 
3048c85b835SJames Wright   @param[in]  ceed          `Ceed` object
3058c85b835SJames Wright   @param[in]  dm            `DMPlex` grid
3068c85b835SJames Wright   @param[in]  domain_label  `DMLabel` for surface
3078c85b835SJames Wright   @param[in]  label_value   `DMPlex` label value for surface
3088c85b835SJames Wright   @param[in]  x_coord       `CeedVector` for coordinates
30900dbc7b1SJames Wright   @param[out] elem_restr_qd `CeedElemRestriction` for QData
3108c85b835SJames Wright   @param[out] q_data        `CeedVector` holding QData
3118c85b835SJames Wright   @param[out] q_data_size   Number of QData components per quadrature point
3128c85b835SJames Wright **/
3138c85b835SJames Wright PetscErrorCode QDataBoundaryGradientGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedVector x_coord,
31400dbc7b1SJames Wright                                         CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size) {
3158c85b835SJames Wright   PetscInt            dim;
3168c85b835SJames Wright   const PetscInt      height_cell = 0, height_face = 1;
3178c85b835SJames Wright   CeedInt             num_comp_x;
318f17df9b6SJames Wright   CeedElemRestriction elem_restr_x_cell, elem_restr_x_face, elem_restr_qd_created;
319f17df9b6SJames Wright   CeedVector          q_data_created;
3208c85b835SJames Wright   CeedBasis           basis_x_cell_to_face, basis_x_face;
321da8b59d6SJames Wright   CeedQFunction       qf_setup_sur = NULL;
3228c85b835SJames Wright   CeedOperator        op_setup_sur;
3238c85b835SJames Wright 
3248c85b835SJames Wright   PetscFunctionBeginUser;
3258c85b835SJames Wright   PetscCall(CreateCoordinateBasisFromPlex(ceed, dm, domain_label, label_value, height_face, &basis_x_face));
3268c85b835SJames Wright   PetscCall(DMPlexCeedBasisCellToFaceCoordinateCreate(ceed, dm, domain_label, label_value, label_value, &basis_x_cell_to_face));
32700dbc7b1SJames Wright   PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height_cell, &elem_restr_x_cell));
3288c85b835SJames Wright   PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height_face, &elem_restr_x_face));
3298c85b835SJames Wright 
3308c85b835SJames Wright   PetscCall(QDataBoundaryGradientGetNumComponents(dm, q_data_size));
331f17df9b6SJames Wright   PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height_face, *q_data_size, &elem_restr_qd_created));
332f17df9b6SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_qd_created, &q_data_created, NULL));
3338c85b835SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_x_face, &num_comp_x));
3348c85b835SJames Wright 
335da8b59d6SJames Wright   PetscCall(DMGetDimension(dm, &dim));
336da8b59d6SJames Wright   switch (dim) {
337da8b59d6SJames Wright     case 2:
338da8b59d6SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Setup2DBoundaryGradient, Setup2DBoundaryGradient_loc, &qf_setup_sur));
339da8b59d6SJames Wright       break;
340da8b59d6SJames Wright     case 3:
3418c85b835SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, SetupBoundaryGradient, SetupBoundaryGradient_loc, &qf_setup_sur));
342da8b59d6SJames Wright       break;
343da8b59d6SJames Wright   }
344da8b59d6SJames Wright   PetscCheck(qf_setup_sur, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "%s not valid for DM of dimension %" PetscInt_FMT, __func__, dim);
345da8b59d6SJames Wright 
3468c85b835SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup_sur, "dx/dX cell", num_comp_x * (dim - height_cell), CEED_EVAL_GRAD));
3478c85b835SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup_sur, "dx/dX face", num_comp_x * (dim - height_face), CEED_EVAL_GRAD));
3488c85b835SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup_sur, "weight", 1, CEED_EVAL_WEIGHT));
3498c85b835SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_setup_sur, "q data", *q_data_size, CEED_EVAL_NONE));
3508c85b835SJames Wright 
3518c85b835SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_setup_sur, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_setup_sur));
35200dbc7b1SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "dx/dX cell", elem_restr_x_cell, basis_x_cell_to_face, CEED_VECTOR_ACTIVE));
3538c85b835SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "dx/dX face", elem_restr_x_face, basis_x_face, CEED_VECTOR_ACTIVE));
3548c85b835SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, basis_x_face, CEED_VECTOR_NONE));
355f17df9b6SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "q data", elem_restr_qd_created, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
3568c85b835SJames Wright 
357f17df9b6SJames Wright   PetscCallCeed(ceed, CeedOperatorApply(op_setup_sur, x_coord, q_data_created, CEED_REQUEST_IMMEDIATE));
358f17df9b6SJames Wright 
359f17df9b6SJames Wright   PetscCall(QDataGetStored(q_data_created, elem_restr_qd_created, q_data, elem_restr_qd));
3608c85b835SJames Wright 
3618c85b835SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup_sur));
3628c85b835SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_setup_sur));
363f17df9b6SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&q_data_created));
364f17df9b6SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd_created));
3658c85b835SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x_face));
36600dbc7b1SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x_cell));
3678c85b835SJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_face));
3688c85b835SJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_cell_to_face));
3698c85b835SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3708c85b835SJames Wright }
371