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
69018c49aSJames Wright #include <dm-utils.h>
7c864c5abSJames Wright #include <petscsection.h>
8c864c5abSJames Wright #include "../qfunctions/setupgeo.h"
9c864c5abSJames Wright #include "../qfunctions/setupgeo2d.h"
10c864c5abSJames Wright
11e816a7e4SJames Wright static CeedVector *q_data_vecs;
12e816a7e4SJames Wright static CeedElemRestriction *q_data_restrictions;
13e816a7e4SJames Wright static PetscInt num_q_data_stored;
14e816a7e4SJames Wright
15e816a7e4SJames Wright /**
16e816a7e4SJames Wright @brief Get stored QData objects that match created objects, if present
17e816a7e4SJames Wright
18e816a7e4SJames Wright If created objects are not present, they are added to the storage and returned in the output
19e816a7e4SJames Wright
20e816a7e4SJames Wright Not Collective
21e816a7e4SJames Wright
22e816a7e4SJames Wright @param[in] q_data_created Vector with created QData
23e816a7e4SJames Wright @param[in] elem_restr_qd_created Restriction for created QData
24e816a7e4SJames Wright @param[out] q_data_stored Vector from storage matching QData
25e816a7e4SJames Wright @param[out] elem_restr_qd_stored Restriction from storage matching QData
26e816a7e4SJames Wright **/
QDataGetStored(CeedVector q_data_created,CeedElemRestriction elem_restr_qd_created,CeedVector * q_data_stored,CeedElemRestriction * elem_restr_qd_stored)27e816a7e4SJames Wright static PetscErrorCode QDataGetStored(CeedVector q_data_created, CeedElemRestriction elem_restr_qd_created, CeedVector *q_data_stored,
28e816a7e4SJames Wright CeedElemRestriction *elem_restr_qd_stored) {
29e816a7e4SJames Wright Ceed ceed = CeedVectorReturnCeed(q_data_created);
30e816a7e4SJames Wright CeedSize created_length, stored_length;
31e816a7e4SJames Wright PetscInt q_data_stored_index = -1;
32e816a7e4SJames Wright
33e816a7e4SJames Wright PetscFunctionBeginUser;
34e816a7e4SJames Wright PetscCallCeed(ceed, CeedVectorGetLength(q_data_created, &created_length));
35e816a7e4SJames Wright for (PetscInt i = 0; i < num_q_data_stored; i++) {
36e816a7e4SJames Wright CeedVector difference_cvec;
37e816a7e4SJames Wright CeedScalar max_difference;
38e816a7e4SJames Wright
399018c49aSJames Wright PetscCallCeed(ceed, CeedVectorGetLength(q_data_vecs[i], &stored_length));
40e816a7e4SJames Wright if (created_length != stored_length) continue;
41e816a7e4SJames Wright PetscCallCeed(ceed, CeedVectorCreate(ceed, stored_length, &difference_cvec));
42e816a7e4SJames Wright PetscCallCeed(ceed, CeedVectorCopy(q_data_vecs[i], difference_cvec));
43e816a7e4SJames Wright PetscCallCeed(ceed, CeedVectorAXPY(difference_cvec, -1, q_data_created));
44e816a7e4SJames Wright PetscCallCeed(ceed, CeedVectorNorm(difference_cvec, CEED_NORM_MAX, &max_difference));
45d033c862SJames Wright //TODO Need to reduce across ranks to ensure all ranks are consistent (not a race condition though since the data is purely local)
46be29160dSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&difference_cvec));
47e816a7e4SJames Wright if (max_difference < 100 * CEED_EPSILON) {
48e816a7e4SJames Wright q_data_stored_index = i;
49e816a7e4SJames Wright break;
50e816a7e4SJames Wright }
51e816a7e4SJames Wright }
52e816a7e4SJames Wright
53e816a7e4SJames Wright if (q_data_stored_index == -1) {
54e816a7e4SJames Wright q_data_stored_index = num_q_data_stored++;
55e816a7e4SJames Wright
56e816a7e4SJames Wright PetscCall(PetscRealloc(num_q_data_stored * sizeof(CeedVector), &q_data_vecs));
57e816a7e4SJames Wright PetscCall(PetscRealloc(num_q_data_stored * sizeof(CeedElemRestriction), &q_data_restrictions));
58e816a7e4SJames Wright q_data_vecs[q_data_stored_index] = NULL; // Must set to NULL for ReferenceCopy
59e816a7e4SJames Wright q_data_restrictions[q_data_stored_index] = NULL; // Must set to NULL for ReferenceCopy
60e816a7e4SJames Wright PetscCallCeed(ceed, CeedVectorReferenceCopy(q_data_created, &q_data_vecs[q_data_stored_index]));
61e816a7e4SJames Wright PetscCallCeed(ceed, CeedElemRestrictionReferenceCopy(elem_restr_qd_created, &q_data_restrictions[q_data_stored_index]));
62e816a7e4SJames Wright }
63e816a7e4SJames Wright *q_data_stored = NULL; // Must set to NULL for ReferenceCopy
64e816a7e4SJames Wright *elem_restr_qd_stored = NULL; // Must set to NULL for ReferenceCopy
65e816a7e4SJames Wright PetscCallCeed(ceed, CeedVectorReferenceCopy(q_data_vecs[q_data_stored_index], q_data_stored));
66e816a7e4SJames Wright PetscCallCeed(ceed, CeedElemRestrictionReferenceCopy(q_data_restrictions[q_data_stored_index], elem_restr_qd_stored));
67e816a7e4SJames Wright PetscFunctionReturn(PETSC_SUCCESS);
68e816a7e4SJames Wright }
69e816a7e4SJames Wright
70e816a7e4SJames Wright /**
71e816a7e4SJames Wright @brief Clear stored QData objects
72e816a7e4SJames Wright **/
QDataClearStoredData()73e816a7e4SJames Wright PetscErrorCode QDataClearStoredData() {
74e816a7e4SJames Wright PetscFunctionBeginUser;
75e816a7e4SJames Wright for (PetscInt i = 0; i < num_q_data_stored; i++) {
76e816a7e4SJames Wright Ceed ceed = CeedVectorReturnCeed(q_data_vecs[i]);
77e816a7e4SJames Wright
78e816a7e4SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data_vecs[i]));
79e816a7e4SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&q_data_restrictions[i]));
80e816a7e4SJames Wright }
81e816a7e4SJames Wright PetscCall(PetscFree(q_data_vecs));
82e816a7e4SJames Wright PetscCall(PetscFree(q_data_restrictions));
83e816a7e4SJames Wright q_data_vecs = NULL;
84e816a7e4SJames Wright q_data_restrictions = NULL;
85e816a7e4SJames Wright PetscFunctionReturn(PETSC_SUCCESS);
86e816a7e4SJames Wright }
87e816a7e4SJames Wright
88c864c5abSJames Wright /**
89c864c5abSJames Wright * @brief Get number of components of quadrature data for domain
90c864c5abSJames Wright *
91c864c5abSJames Wright * @param[in] dm DM where quadrature data would be used
92c864c5abSJames Wright * @param[out] q_data_size Number of components of quadrature data
93c864c5abSJames Wright */
QDataGetNumComponents(DM dm,CeedInt * q_data_size)94c864c5abSJames Wright PetscErrorCode QDataGetNumComponents(DM dm, CeedInt *q_data_size) {
95c864c5abSJames Wright PetscInt num_comp_x, dim;
96c864c5abSJames Wright
97c864c5abSJames Wright PetscFunctionBeginUser;
98c864c5abSJames Wright PetscCall(DMGetDimension(dm, &dim));
99bc87537dSJames Wright PetscCall(DMGetCoordinateNumComps(dm, &num_comp_x));
100da8b59d6SJames Wright *q_data_size = 0;
101c864c5abSJames Wright switch (dim) {
102c864c5abSJames Wright case 2:
103c864c5abSJames Wright switch (num_comp_x) {
104c864c5abSJames Wright case 2:
105c864c5abSJames Wright *q_data_size = 5;
106c864c5abSJames Wright break;
107c864c5abSJames Wright case 3:
108c864c5abSJames Wright *q_data_size = 7;
109c864c5abSJames Wright break;
110c864c5abSJames Wright }
111c864c5abSJames Wright break;
112c864c5abSJames Wright case 3:
113c864c5abSJames Wright *q_data_size = 10;
114c864c5abSJames Wright break;
115c864c5abSJames Wright }
116da8b59d6SJames Wright PetscCheck(*q_data_size, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP,
117da8b59d6SJames Wright "QData not valid for DM of dimension %" PetscInt_FMT " and coordinates with dimension %" PetscInt_FMT, dim, num_comp_x);
118c864c5abSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
119c864c5abSJames Wright }
120c864c5abSJames Wright
121c864c5abSJames Wright /**
122c864c5abSJames Wright * @brief Create quadrature data for domain
123c864c5abSJames Wright *
124c864c5abSJames Wright * @param[in] ceed Ceed object quadrature data will be used with
125c864c5abSJames Wright * @param[in] dm DM where quadrature data would be used
126c864c5abSJames Wright * @param[in] domain_label DMLabel that quadrature data would be used one
127c864c5abSJames Wright * @param[in] label_value Value of label
128c864c5abSJames Wright * @param[out] elem_restr_qd CeedElemRestriction of the quadrature data
129c864c5abSJames Wright * @param[out] q_data CeedVector of the quadrature data
1304aea4664SJames Wright * @param[out] q_data_size Number of components of quadrature data
131c864c5abSJames Wright */
QDataGet(Ceed ceed,DM dm,DMLabel domain_label,PetscInt label_value,CeedElemRestriction * elem_restr_qd,CeedVector * q_data,CeedInt * q_data_size)1329018c49aSJames Wright PetscErrorCode QDataGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedElemRestriction *elem_restr_qd, CeedVector *q_data,
1339018c49aSJames Wright CeedInt *q_data_size) {
13487d3884fSJeremy L Thompson CeedQFunction qf_setup = NULL;
135c864c5abSJames Wright CeedOperator op_setup;
1369018c49aSJames Wright CeedVector q_data_created, x_coord;
1379018c49aSJames Wright CeedElemRestriction elem_restr_qd_created, elem_restr_x;
1389018c49aSJames Wright PetscInt dim, height = 0, num_comp_x;
1399018c49aSJames Wright CeedBasis basis_x;
140c864c5abSJames Wright
141c864c5abSJames Wright PetscFunctionBeginUser;
1429018c49aSJames Wright PetscCall(DMPlexCeedCoordinateCreateField(ceed, dm, domain_label, label_value, height, &elem_restr_x, &basis_x, &x_coord));
143c864c5abSJames Wright PetscCall(QDataGetNumComponents(dm, q_data_size));
1449018c49aSJames Wright PetscCall(DMGetCoordinateNumComps(dm, &num_comp_x));
145c864c5abSJames Wright PetscCall(DMGetDimension(dm, &dim));
146c864c5abSJames Wright switch (dim) {
147c864c5abSJames Wright case 2:
148c864c5abSJames Wright switch (num_comp_x) {
149c864c5abSJames Wright case 2:
150c864c5abSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Setup2d, Setup2d_loc, &qf_setup));
151c864c5abSJames Wright break;
152c864c5abSJames Wright case 3:
153c864c5abSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Setup2D_3Dcoords, Setup2D_3Dcoords_loc, &qf_setup));
154c864c5abSJames Wright break;
155c864c5abSJames Wright }
156c864c5abSJames Wright break;
157c864c5abSJames Wright case 3:
158c864c5abSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Setup, Setup_loc, &qf_setup));
159c864c5abSJames Wright break;
160c864c5abSJames Wright }
161da8b59d6SJames Wright PetscCheck(qf_setup, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "%s not valid for DM of dimension %" PetscInt_FMT, __func__, dim);
162c864c5abSJames Wright
163c864c5abSJames Wright // -- Create QFunction for quadrature data
164c864c5abSJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_setup, 0));
165c864c5abSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup, "dx", num_comp_x * (dim - height), CEED_EVAL_GRAD));
166c864c5abSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup, "weight", 1, CEED_EVAL_WEIGHT));
167c864c5abSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_setup, "surface qdata", *q_data_size, CEED_EVAL_NONE));
168c864c5abSJames Wright
169e816a7e4SJames Wright PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, *q_data_size, &elem_restr_qd_created));
170e816a7e4SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_qd_created, &q_data_created, NULL));
171c864c5abSJames Wright
172c864c5abSJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_setup, NULL, NULL, &op_setup));
173c864c5abSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "dx", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE));
174c864c5abSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "weight", CEED_ELEMRESTRICTION_NONE, basis_x, CEED_VECTOR_NONE));
175e816a7e4SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "surface qdata", elem_restr_qd_created, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
176c864c5abSJames Wright
177e816a7e4SJames Wright PetscCallCeed(ceed, CeedOperatorApply(op_setup, x_coord, q_data_created, CEED_REQUEST_IMMEDIATE));
178c864c5abSJames Wright
179e816a7e4SJames Wright PetscCall(QDataGetStored(q_data_created, elem_restr_qd_created, q_data, elem_restr_qd));
180e816a7e4SJames Wright
181e816a7e4SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data_created));
182e816a7e4SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd_created));
1839018c49aSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x));
1849018c49aSJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_x));
1859018c49aSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&x_coord));
186c864c5abSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup));
187c864c5abSJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_setup));
188c864c5abSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
189c864c5abSJames Wright }
190c864c5abSJames Wright
191c864c5abSJames Wright /**
192c864c5abSJames Wright * @brief Get number of components of quadrature data for boundary of domain
193c864c5abSJames Wright *
194c864c5abSJames Wright * @param[in] dm DM where quadrature data would be used
195c864c5abSJames Wright * @param[out] q_data_size Number of components of quadrature data
196c864c5abSJames Wright */
QDataBoundaryGetNumComponents(DM dm,CeedInt * q_data_size)197c864c5abSJames Wright PetscErrorCode QDataBoundaryGetNumComponents(DM dm, CeedInt *q_data_size) {
198c864c5abSJames Wright PetscInt dim;
199c864c5abSJames Wright
200c864c5abSJames Wright PetscFunctionBeginUser;
201c864c5abSJames Wright PetscCall(DMGetDimension(dm, &dim));
202da8b59d6SJames Wright *q_data_size = 0;
203c864c5abSJames Wright switch (dim) {
204c864c5abSJames Wright case 2:
205c864c5abSJames Wright *q_data_size = 3;
206c864c5abSJames Wright break;
207c864c5abSJames Wright case 3:
208c864c5abSJames Wright *q_data_size = 10;
209c864c5abSJames Wright break;
210c864c5abSJames Wright }
211da8b59d6SJames Wright PetscCheck(*q_data_size, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "%s not valid for DM of dimension %" PetscInt_FMT, __func__, dim);
212c864c5abSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
213c864c5abSJames Wright }
214c864c5abSJames Wright
215c864c5abSJames Wright /**
216c864c5abSJames Wright * @brief Create quadrature data for boundary of domain
217c864c5abSJames Wright *
218c864c5abSJames Wright * @param[in] ceed Ceed object quadrature data will be used with
219c864c5abSJames Wright * @param[in] dm DM where quadrature data would be used
220c864c5abSJames Wright * @param[in] domain_label DMLabel that quadrature data would be used one
221c864c5abSJames Wright * @param[in] label_value Value of label
222c864c5abSJames Wright * @param[out] elem_restr_qd CeedElemRestriction of the quadrature data
223c864c5abSJames Wright * @param[out] q_data CeedVector of the quadrature data
2244aea4664SJames Wright * @param[out] q_data_size Number of components of quadrature data
225c864c5abSJames Wright */
QDataBoundaryGet(Ceed ceed,DM dm,DMLabel domain_label,PetscInt label_value,CeedElemRestriction * elem_restr_qd,CeedVector * q_data,CeedInt * q_data_size)22607126c9aSJames Wright PetscErrorCode QDataBoundaryGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedElemRestriction *elem_restr_qd, CeedVector *q_data,
22707126c9aSJames Wright CeedInt *q_data_size) {
22887d3884fSJeremy L Thompson CeedQFunction qf_setup_sur = NULL;
229c864c5abSJames Wright CeedOperator op_setup_sur;
23007126c9aSJames Wright CeedVector q_data_created, x_coord;
23107126c9aSJames Wright CeedElemRestriction elem_restr_qd_created, elem_restr_x;
23207126c9aSJames Wright CeedBasis basis_x;
23307126c9aSJames Wright PetscInt dim, height = 1, num_comp_x;
234c864c5abSJames Wright
235c864c5abSJames Wright PetscFunctionBeginUser;
23607126c9aSJames Wright PetscCall(DMPlexCeedCoordinateCreateField(ceed, dm, domain_label, label_value, height, &elem_restr_x, &basis_x, &x_coord));
237c864c5abSJames Wright PetscCall(QDataBoundaryGetNumComponents(dm, q_data_size));
23807126c9aSJames Wright PetscCall(DMGetCoordinateNumComps(dm, &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));
27007126c9aSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x));
27107126c9aSJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_x));
27207126c9aSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&x_coord));
273c864c5abSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup_sur));
274c864c5abSJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_setup_sur));
275c864c5abSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
276c864c5abSJames Wright }
2778c85b835SJames Wright
2788c85b835SJames Wright /**
2798c85b835SJames Wright * @brief Get number of components of quadrature data for boundary of domain
2808c85b835SJames Wright *
2818c85b835SJames Wright * @param[in] dm DM where quadrature data would be used
2828c85b835SJames Wright * @param[out] q_data_size Number of components of quadrature data
2838c85b835SJames Wright */
QDataBoundaryGradientGetNumComponents(DM dm,CeedInt * q_data_size)2848c85b835SJames Wright PetscErrorCode QDataBoundaryGradientGetNumComponents(DM dm, CeedInt *q_data_size) {
2858c85b835SJames Wright PetscInt dim;
2868c85b835SJames Wright
2878c85b835SJames Wright PetscFunctionBeginUser;
2888c85b835SJames Wright PetscCall(DMGetDimension(dm, &dim));
289da8b59d6SJames Wright *q_data_size = 0;
2908c85b835SJames Wright switch (dim) {
291da8b59d6SJames Wright case 2:
292da8b59d6SJames Wright *q_data_size = 7;
293da8b59d6SJames Wright break;
2948c85b835SJames Wright case 3:
2958c85b835SJames Wright *q_data_size = 13;
2968c85b835SJames Wright break;
2978c85b835SJames Wright }
298da8b59d6SJames Wright PetscCheck(*q_data_size, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "%s not valid for DM of dimension %" PetscInt_FMT, __func__, dim);
2998c85b835SJames Wright PetscFunctionReturn(PETSC_SUCCESS);
3008c85b835SJames Wright }
3018c85b835SJames Wright
3028c85b835SJames Wright /**
3038c85b835SJames Wright @brief Compute `CeedOperator` surface gradient QData
3048c85b835SJames Wright
3058c85b835SJames Wright Collective across MPI processes.
3068c85b835SJames Wright
3078c85b835SJames Wright @param[in] ceed `Ceed` object
3088c85b835SJames Wright @param[in] dm `DMPlex` grid
3098c85b835SJames Wright @param[in] domain_label `DMLabel` for surface
3108c85b835SJames Wright @param[in] label_value `DMPlex` label value for surface
31100dbc7b1SJames Wright @param[out] elem_restr_qd `CeedElemRestriction` for QData
3128c85b835SJames Wright @param[out] q_data `CeedVector` holding QData
3138c85b835SJames Wright @param[out] q_data_size Number of QData components per quadrature point
3148c85b835SJames Wright **/
QDataBoundaryGradientGet(Ceed ceed,DM dm,DMLabel domain_label,PetscInt label_value,CeedElemRestriction * elem_restr_qd,CeedVector * q_data,CeedInt * q_data_size)315*d510ce3cSJames Wright PetscErrorCode QDataBoundaryGradientGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedElemRestriction *elem_restr_qd,
316*d510ce3cSJames Wright CeedVector *q_data, CeedInt *q_data_size) {
3178c85b835SJames Wright PetscInt dim;
3188c85b835SJames Wright const PetscInt height_cell = 0, height_face = 1;
3198c85b835SJames Wright CeedInt num_comp_x;
320f17df9b6SJames Wright CeedElemRestriction elem_restr_x_cell, elem_restr_x_face, elem_restr_qd_created;
321*d510ce3cSJames Wright CeedVector q_data_created, x_coord;
3228c85b835SJames Wright CeedBasis basis_x_cell_to_face, basis_x_face;
323da8b59d6SJames Wright CeedQFunction qf_setup_sur = NULL;
3248c85b835SJames Wright CeedOperator op_setup_sur;
3258c85b835SJames Wright
3268c85b835SJames Wright PetscFunctionBeginUser;
327*d510ce3cSJames Wright PetscCall(DMPlexCeedCoordinateCreateField(ceed, dm, domain_label, label_value, height_face, &elem_restr_x_face, &basis_x_face, &x_coord));
3288c85b835SJames Wright PetscCall(DMPlexCeedBasisCellToFaceCoordinateCreate(ceed, dm, domain_label, label_value, label_value, &basis_x_cell_to_face));
32900dbc7b1SJames Wright PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height_cell, &elem_restr_x_cell));
3308c85b835SJames Wright
3318c85b835SJames Wright PetscCall(QDataBoundaryGradientGetNumComponents(dm, q_data_size));
332f17df9b6SJames Wright PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height_face, *q_data_size, &elem_restr_qd_created));
333f17df9b6SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_qd_created, &q_data_created, NULL));
3348c85b835SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_x_face, &num_comp_x));
3358c85b835SJames Wright
336da8b59d6SJames Wright PetscCall(DMGetDimension(dm, &dim));
337da8b59d6SJames Wright switch (dim) {
338da8b59d6SJames Wright case 2:
339da8b59d6SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Setup2DBoundaryGradient, Setup2DBoundaryGradient_loc, &qf_setup_sur));
340da8b59d6SJames Wright break;
341da8b59d6SJames Wright case 3:
3428c85b835SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, SetupBoundaryGradient, SetupBoundaryGradient_loc, &qf_setup_sur));
343da8b59d6SJames Wright break;
344da8b59d6SJames Wright }
345da8b59d6SJames Wright PetscCheck(qf_setup_sur, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "%s not valid for DM of dimension %" PetscInt_FMT, __func__, dim);
346da8b59d6SJames Wright
3478c85b835SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup_sur, "dx/dX cell", num_comp_x * (dim - height_cell), CEED_EVAL_GRAD));
3488c85b835SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup_sur, "dx/dX face", num_comp_x * (dim - height_face), CEED_EVAL_GRAD));
3498c85b835SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup_sur, "weight", 1, CEED_EVAL_WEIGHT));
3508c85b835SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_setup_sur, "q data", *q_data_size, CEED_EVAL_NONE));
3518c85b835SJames Wright
3528c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_setup_sur, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_setup_sur));
35300dbc7b1SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "dx/dX cell", elem_restr_x_cell, basis_x_cell_to_face, CEED_VECTOR_ACTIVE));
3548c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "dx/dX face", elem_restr_x_face, basis_x_face, CEED_VECTOR_ACTIVE));
3558c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, basis_x_face, CEED_VECTOR_NONE));
356f17df9b6SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "q data", elem_restr_qd_created, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
3578c85b835SJames Wright
358f17df9b6SJames Wright PetscCallCeed(ceed, CeedOperatorApply(op_setup_sur, x_coord, q_data_created, CEED_REQUEST_IMMEDIATE));
359f17df9b6SJames Wright
360f17df9b6SJames Wright PetscCall(QDataGetStored(q_data_created, elem_restr_qd_created, q_data, elem_restr_qd));
3618c85b835SJames Wright
3628c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup_sur));
3638c85b835SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_setup_sur));
364f17df9b6SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data_created));
365*d510ce3cSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&x_coord));
366f17df9b6SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd_created));
3678c85b835SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x_face));
36800dbc7b1SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x_cell));
3698c85b835SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_face));
3708c85b835SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_cell_to_face));
3718c85b835SJames Wright PetscFunctionReturn(PETSC_SUCCESS);
3728c85b835SJames Wright }
373