xref: /honee/src/qdata.c (revision e816a7e40144ce229b5323f192a798e8c639fa87)
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 
10*e816a7e4SJames Wright static CeedVector          *q_data_vecs;
11*e816a7e4SJames Wright static CeedElemRestriction *q_data_restrictions;
12*e816a7e4SJames Wright static PetscInt             num_q_data_stored;
13*e816a7e4SJames Wright 
14*e816a7e4SJames Wright /**
15*e816a7e4SJames Wright   @brief Get stored QData objects that match created objects, if present
16*e816a7e4SJames Wright 
17*e816a7e4SJames Wright   If created objects are not present, they are added to the storage and returned in the output
18*e816a7e4SJames Wright 
19*e816a7e4SJames Wright   Not Collective
20*e816a7e4SJames Wright 
21*e816a7e4SJames Wright   @param[in]  q_data_created        Vector with created QData
22*e816a7e4SJames Wright   @param[in]  elem_restr_qd_created Restriction for created QData
23*e816a7e4SJames Wright   @param[out] q_data_stored         Vector from storage matching QData
24*e816a7e4SJames Wright   @param[out] elem_restr_qd_stored  Restriction from storage matching QData
25*e816a7e4SJames Wright **/
26*e816a7e4SJames Wright static PetscErrorCode QDataGetStored(CeedVector q_data_created, CeedElemRestriction elem_restr_qd_created, CeedVector *q_data_stored,
27*e816a7e4SJames Wright                                      CeedElemRestriction *elem_restr_qd_stored) {
28*e816a7e4SJames Wright   Ceed     ceed = CeedVectorReturnCeed(q_data_created);
29*e816a7e4SJames Wright   CeedSize created_length, stored_length;
30*e816a7e4SJames Wright   PetscInt q_data_stored_index = -1;
31*e816a7e4SJames Wright 
32*e816a7e4SJames Wright   PetscFunctionBeginUser;
33*e816a7e4SJames Wright   PetscCallCeed(ceed, CeedVectorGetLength(q_data_created, &created_length));
34*e816a7e4SJames Wright   for (PetscInt i = 0; i < num_q_data_stored; i++) {
35*e816a7e4SJames Wright     CeedVector difference_cvec;
36*e816a7e4SJames Wright     CeedScalar max_difference;
37*e816a7e4SJames Wright 
38*e816a7e4SJames Wright     PetscCallCeed(ceed, CeedVectorGetLength(q_data_vecs[0], &stored_length));
39*e816a7e4SJames Wright     if (created_length != stored_length) continue;
40*e816a7e4SJames Wright     PetscCallCeed(ceed, CeedVectorCreate(ceed, stored_length, &difference_cvec));
41*e816a7e4SJames Wright     PetscCallCeed(ceed, CeedVectorCopy(q_data_vecs[i], difference_cvec));
42*e816a7e4SJames Wright     PetscCallCeed(ceed, CeedVectorAXPY(difference_cvec, -1, q_data_created));
43*e816a7e4SJames Wright     PetscCallCeed(ceed, CeedVectorNorm(difference_cvec, CEED_NORM_MAX, &max_difference));
44*e816a7e4SJames Wright     if (max_difference < 100 * CEED_EPSILON) {
45*e816a7e4SJames Wright       q_data_stored_index = i;
46*e816a7e4SJames Wright       break;
47*e816a7e4SJames Wright     }
48*e816a7e4SJames Wright   }
49*e816a7e4SJames Wright 
50*e816a7e4SJames Wright   if (q_data_stored_index == -1) {
51*e816a7e4SJames Wright     q_data_stored_index = num_q_data_stored++;
52*e816a7e4SJames Wright 
53*e816a7e4SJames Wright     PetscCall(PetscRealloc(num_q_data_stored * sizeof(CeedVector), &q_data_vecs));
54*e816a7e4SJames Wright     PetscCall(PetscRealloc(num_q_data_stored * sizeof(CeedElemRestriction), &q_data_restrictions));
55*e816a7e4SJames Wright     q_data_vecs[q_data_stored_index]         = NULL;  // Must set to NULL for ReferenceCopy
56*e816a7e4SJames Wright     q_data_restrictions[q_data_stored_index] = NULL;  // Must set to NULL for ReferenceCopy
57*e816a7e4SJames Wright     PetscCallCeed(ceed, CeedVectorReferenceCopy(q_data_created, &q_data_vecs[q_data_stored_index]));
58*e816a7e4SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionReferenceCopy(elem_restr_qd_created, &q_data_restrictions[q_data_stored_index]));
59*e816a7e4SJames Wright   }
60*e816a7e4SJames Wright   *q_data_stored        = NULL;  // Must set to NULL for ReferenceCopy
61*e816a7e4SJames Wright   *elem_restr_qd_stored = NULL;  // Must set to NULL for ReferenceCopy
62*e816a7e4SJames Wright   PetscCallCeed(ceed, CeedVectorReferenceCopy(q_data_vecs[q_data_stored_index], q_data_stored));
63*e816a7e4SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionReferenceCopy(q_data_restrictions[q_data_stored_index], elem_restr_qd_stored));
64*e816a7e4SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
65*e816a7e4SJames Wright }
66*e816a7e4SJames Wright 
67*e816a7e4SJames Wright /**
68*e816a7e4SJames Wright   @brief Clear stored QData objects
69*e816a7e4SJames Wright **/
70*e816a7e4SJames Wright PetscErrorCode QDataClearStoredData() {
71*e816a7e4SJames Wright   PetscFunctionBeginUser;
72*e816a7e4SJames Wright   for (PetscInt i = 0; i < num_q_data_stored; i++) {
73*e816a7e4SJames Wright     Ceed ceed = CeedVectorReturnCeed(q_data_vecs[i]);
74*e816a7e4SJames Wright 
75*e816a7e4SJames Wright     PetscCallCeed(ceed, CeedVectorDestroy(&q_data_vecs[i]));
76*e816a7e4SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&q_data_restrictions[i]));
77*e816a7e4SJames Wright   }
78*e816a7e4SJames Wright   PetscCall(PetscFree(q_data_vecs));
79*e816a7e4SJames Wright   PetscCall(PetscFree(q_data_restrictions));
80*e816a7e4SJames Wright   q_data_vecs         = NULL;
81*e816a7e4SJames Wright   q_data_restrictions = NULL;
82*e816a7e4SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
83*e816a7e4SJames Wright }
84*e816a7e4SJames Wright 
85c864c5abSJames Wright /**
86c864c5abSJames Wright  * @brief Get number of components of quadrature data for domain
87c864c5abSJames Wright  *
88c864c5abSJames Wright  * @param[in]  dm          DM where quadrature data would be used
89c864c5abSJames Wright  * @param[out] q_data_size Number of components of quadrature data
90c864c5abSJames Wright  */
91c864c5abSJames Wright PetscErrorCode QDataGetNumComponents(DM dm, CeedInt *q_data_size) {
92c864c5abSJames Wright   PetscInt num_comp_x, dim;
93c864c5abSJames Wright 
94c864c5abSJames Wright   PetscFunctionBeginUser;
95c864c5abSJames Wright   PetscCall(DMGetDimension(dm, &dim));
96c864c5abSJames Wright   {  // Get number of coordinate components
97c864c5abSJames Wright     DM           dm_coord;
98c864c5abSJames Wright     PetscSection section_coord;
99c864c5abSJames Wright     PetscInt     field = 0;  // Default field has the coordinates
100c864c5abSJames Wright     PetscCall(DMGetCoordinateDM(dm, &dm_coord));
101c864c5abSJames Wright     PetscCall(DMGetLocalSection(dm_coord, &section_coord));
102c864c5abSJames Wright     PetscCall(PetscSectionGetFieldComponents(section_coord, field, &num_comp_x));
103c864c5abSJames Wright   }
104c864c5abSJames Wright   switch (dim) {
105c864c5abSJames Wright     case 2:
106c864c5abSJames Wright       switch (num_comp_x) {
107c864c5abSJames Wright         case 2:
108c864c5abSJames Wright           *q_data_size = 5;
109c864c5abSJames Wright           break;
110c864c5abSJames Wright         case 3:
111c864c5abSJames Wright           *q_data_size = 7;
112c864c5abSJames Wright           break;
113c864c5abSJames Wright         default:
114c864c5abSJames Wright           SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP,
115c864c5abSJames Wright                   "QData not valid for DM of dimension %" PetscInt_FMT " and coordinates with dimension %" PetscInt_FMT, dim, num_comp_x);
116c864c5abSJames Wright           break;
117c864c5abSJames Wright       }
118c864c5abSJames Wright       break;
119c864c5abSJames Wright     case 3:
120c864c5abSJames Wright       *q_data_size = 10;
121c864c5abSJames Wright       break;
122c864c5abSJames Wright     default:
123c864c5abSJames Wright       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP,
124c864c5abSJames Wright               "QData not valid for DM of dimension %" PetscInt_FMT " and coordinates with dimension %" PetscInt_FMT, dim, num_comp_x);
125c864c5abSJames Wright       break;
126c864c5abSJames Wright   }
127c864c5abSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
128c864c5abSJames Wright }
129c864c5abSJames Wright 
130c864c5abSJames Wright /**
131c864c5abSJames Wright  * @brief Create quadrature data for domain
132c864c5abSJames Wright  *
133c864c5abSJames Wright  * @param[in]  ceed          Ceed object quadrature data will be used with
134c864c5abSJames Wright  * @param[in]  dm            DM where quadrature data would be used
135c864c5abSJames Wright  * @param[in]  domain_label  DMLabel that quadrature data would be used one
136c864c5abSJames Wright  * @param[in]  label_value   Value of label
137c864c5abSJames Wright  * @param[in]  elem_restr_x  CeedElemRestriction of the coordinates (must match `domain_label` and `label_value` selections)
138c864c5abSJames Wright  * @param[in]  basis_x       CeedBasis of the coordinates
139c864c5abSJames Wright  * @param[in]  x_coord       CeedVector of the coordinates
140c864c5abSJames Wright  * @param[out] elem_restr_qd CeedElemRestriction of the quadrature data
141c864c5abSJames Wright  * @param[out] q_data        CeedVector of the quadrature data
142c864c5abSJames Wright  * @param[out] q_data_size   number of components of quadrature data
143c864c5abSJames Wright  */
144c864c5abSJames Wright PetscErrorCode QDataGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedElemRestriction elem_restr_x, CeedBasis basis_x,
145c864c5abSJames Wright                         CeedVector x_coord, CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size) {
14687d3884fSJeremy L Thompson   CeedQFunction       qf_setup = NULL;
147c864c5abSJames Wright   CeedOperator        op_setup;
148*e816a7e4SJames Wright   CeedVector          q_data_created;
149*e816a7e4SJames Wright   CeedElemRestriction elem_restr_qd_created;
150c864c5abSJames Wright   CeedInt             num_comp_x;
151c864c5abSJames Wright   PetscInt            dim, height = 0;
152c864c5abSJames Wright 
153c864c5abSJames Wright   PetscFunctionBeginUser;
154c864c5abSJames Wright   PetscCall(QDataGetNumComponents(dm, q_data_size));
155c864c5abSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_x, &num_comp_x));
156c864c5abSJames Wright   PetscCall(DMGetDimension(dm, &dim));
157c864c5abSJames Wright   switch (dim) {
158c864c5abSJames Wright     case 2:
159c864c5abSJames Wright       switch (num_comp_x) {
160c864c5abSJames Wright         case 2:
161c864c5abSJames Wright           PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Setup2d, Setup2d_loc, &qf_setup));
162c864c5abSJames Wright           break;
163c864c5abSJames Wright         case 3:
164c864c5abSJames Wright           PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Setup2D_3Dcoords, Setup2D_3Dcoords_loc, &qf_setup));
165c864c5abSJames Wright           break;
166c864c5abSJames Wright       }
167c864c5abSJames Wright       break;
168c864c5abSJames Wright     case 3:
169c864c5abSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Setup, Setup_loc, &qf_setup));
170c864c5abSJames Wright       break;
171c864c5abSJames Wright   }
172c864c5abSJames Wright 
173c864c5abSJames Wright   // -- Create QFunction for quadrature data
174c864c5abSJames Wright   PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_setup, 0));
175c864c5abSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup, "dx", num_comp_x * (dim - height), CEED_EVAL_GRAD));
176c864c5abSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup, "weight", 1, CEED_EVAL_WEIGHT));
177c864c5abSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_setup, "surface qdata", *q_data_size, CEED_EVAL_NONE));
178c864c5abSJames Wright 
179*e816a7e4SJames Wright   PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, *q_data_size, &elem_restr_qd_created));
180*e816a7e4SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_qd_created, &q_data_created, NULL));
181c864c5abSJames Wright 
182c864c5abSJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_setup, NULL, NULL, &op_setup));
183c864c5abSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "dx", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE));
184c864c5abSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "weight", CEED_ELEMRESTRICTION_NONE, basis_x, CEED_VECTOR_NONE));
185*e816a7e4SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "surface qdata", elem_restr_qd_created, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
186c864c5abSJames Wright 
187*e816a7e4SJames Wright   PetscCallCeed(ceed, CeedOperatorApply(op_setup, x_coord, q_data_created, CEED_REQUEST_IMMEDIATE));
188c864c5abSJames Wright 
189*e816a7e4SJames Wright   PetscCall(QDataGetStored(q_data_created, elem_restr_qd_created, q_data, elem_restr_qd));
190*e816a7e4SJames Wright 
191*e816a7e4SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&q_data_created));
192*e816a7e4SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd_created));
193c864c5abSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup));
194c864c5abSJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_setup));
195c864c5abSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
196c864c5abSJames Wright }
197c864c5abSJames Wright 
198c864c5abSJames Wright /**
199c864c5abSJames Wright  * @brief Get number of components of quadrature data for boundary of domain
200c864c5abSJames Wright  *
201c864c5abSJames Wright  * @param[in]  dm          DM where quadrature data would be used
202c864c5abSJames Wright  * @param[out] q_data_size Number of components of quadrature data
203c864c5abSJames Wright  */
204c864c5abSJames Wright PetscErrorCode QDataBoundaryGetNumComponents(DM dm, CeedInt *q_data_size) {
205c864c5abSJames Wright   PetscInt dim;
206c864c5abSJames Wright 
207c864c5abSJames Wright   PetscFunctionBeginUser;
208c864c5abSJames Wright   PetscCall(DMGetDimension(dm, &dim));
209c864c5abSJames Wright   switch (dim) {
210c864c5abSJames Wright     case 2:
211c864c5abSJames Wright       *q_data_size = 3;
212c864c5abSJames Wright       break;
213c864c5abSJames Wright     case 3:
214c864c5abSJames Wright       *q_data_size = 10;
215c864c5abSJames Wright       break;
216c864c5abSJames Wright     default:
217c864c5abSJames Wright       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "QDataBoundary not valid for DM of dimension %" PetscInt_FMT, dim);
218c864c5abSJames Wright       break;
219c864c5abSJames Wright   }
220c864c5abSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
221c864c5abSJames Wright }
222c864c5abSJames Wright 
223c864c5abSJames Wright /**
224c864c5abSJames Wright  * @brief Create quadrature data for boundary of domain
225c864c5abSJames Wright  *
226c864c5abSJames Wright  * @param[in]  ceed          Ceed object quadrature data will be used with
227c864c5abSJames Wright  * @param[in]  dm            DM where quadrature data would be used
228c864c5abSJames Wright  * @param[in]  domain_label  DMLabel that quadrature data would be used one
229c864c5abSJames Wright  * @param[in]  label_value   Value of label
230c864c5abSJames Wright  * @param[in]  elem_restr_x  CeedElemRestriction of the coordinates (must match `domain_label` and `label_value` selections)
231c864c5abSJames Wright  * @param[in]  basis_x       CeedBasis of the coordinates
232c864c5abSJames Wright  * @param[in]  x_coord       CeedVector of the coordinates
233c864c5abSJames Wright  * @param[out] elem_restr_qd CeedElemRestriction of the quadrature data
234c864c5abSJames Wright  * @param[out] q_data        CeedVector of the quadrature data
235c864c5abSJames Wright  * @param[out] q_data_size   number of components of quadrature data
236c864c5abSJames Wright  */
237c864c5abSJames Wright PetscErrorCode QDataBoundaryGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedElemRestriction elem_restr_x, CeedBasis basis_x,
238c864c5abSJames Wright                                 CeedVector x_coord, CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size) {
23987d3884fSJeremy L Thompson   CeedQFunction       qf_setup_sur = NULL;
240c864c5abSJames Wright   CeedOperator        op_setup_sur;
241*e816a7e4SJames Wright   CeedVector          q_data_created;
242*e816a7e4SJames Wright   CeedElemRestriction elem_restr_qd_created;
243c864c5abSJames Wright   CeedInt             num_comp_x;
244c864c5abSJames Wright   PetscInt            dim, height = 1;
245c864c5abSJames Wright 
246c864c5abSJames Wright   PetscFunctionBeginUser;
247c864c5abSJames Wright   PetscCall(QDataBoundaryGetNumComponents(dm, q_data_size));
248c864c5abSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_x, &num_comp_x));
249c864c5abSJames Wright   PetscCall(DMGetDimension(dm, &dim));
250c864c5abSJames Wright   switch (dim) {
251c864c5abSJames Wright     case 2:
252c864c5abSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, SetupBoundary2d, SetupBoundary2d_loc, &qf_setup_sur));
253c864c5abSJames Wright       break;
254c864c5abSJames Wright     case 3:
255c864c5abSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, SetupBoundary, SetupBoundary_loc, &qf_setup_sur));
256c864c5abSJames Wright       break;
257c864c5abSJames Wright   }
258c864c5abSJames Wright 
259c864c5abSJames Wright   // -- Create QFunction for quadrature data
260c864c5abSJames Wright   PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_setup_sur, 0));
261c864c5abSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup_sur, "dx", num_comp_x * (dim - height), CEED_EVAL_GRAD));
262c864c5abSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup_sur, "weight", 1, CEED_EVAL_WEIGHT));
263c864c5abSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_setup_sur, "surface qdata", *q_data_size, CEED_EVAL_NONE));
264c864c5abSJames Wright 
265*e816a7e4SJames Wright   PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, *q_data_size, &elem_restr_qd_created));
266*e816a7e4SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_qd_created, &q_data_created, NULL));
267c864c5abSJames Wright 
268c864c5abSJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_setup_sur, NULL, NULL, &op_setup_sur));
269c864c5abSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "dx", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE));
270c864c5abSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, basis_x, CEED_VECTOR_NONE));
271*e816a7e4SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "surface qdata", elem_restr_qd_created, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
272c864c5abSJames Wright 
273*e816a7e4SJames Wright   PetscCallCeed(ceed, CeedOperatorApply(op_setup_sur, x_coord, q_data_created, CEED_REQUEST_IMMEDIATE));
274c864c5abSJames Wright 
275*e816a7e4SJames Wright   PetscCall(QDataGetStored(q_data_created, elem_restr_qd_created, q_data, elem_restr_qd));
276*e816a7e4SJames Wright 
277*e816a7e4SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&q_data_created));
278*e816a7e4SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd_created));
279c864c5abSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup_sur));
280c864c5abSJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_setup_sur));
281c864c5abSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
282c864c5abSJames Wright }
283