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