xref: /libCEED/examples/fluids/src/dm_utils.c (revision fe1e732ecc3813db3c31d69a221ef0a8f0b3ac14)
1d68a66c4SJames Wright // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2d68a66c4SJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3d68a66c4SJames Wright //
4d68a66c4SJames Wright // SPDX-License-Identifier: BSD-2-Clause
5d68a66c4SJames Wright //
6d68a66c4SJames Wright // This file is part of CEED:  http://github.com/ceed
7d68a66c4SJames Wright 
8d68a66c4SJames Wright /// @file
9d68a66c4SJames Wright /// Utilities for setting up DM and PetscFE
10d68a66c4SJames Wright 
11d68a66c4SJames Wright #include <ceed.h>
12d68a66c4SJames Wright #include <petscdmplex.h>
13d68a66c4SJames Wright #include <petscds.h>
14d68a66c4SJames Wright 
15d68a66c4SJames Wright #include "../navierstokes.h"
16d68a66c4SJames Wright 
17d68a66c4SJames Wright /**
18d68a66c4SJames Wright   @brief Convert `DM` field to `DS` field.
19d68a66c4SJames Wright 
20d68a66c4SJames Wright   @param[in]   dm            `DM` holding mesh
21d68a66c4SJames Wright   @param[in]   domain_label  Label for `DM` domain
22d68a66c4SJames Wright   @param[in]   dm_field      Index of `DM` field
23d68a66c4SJames Wright   @param[out]  ds_field      Index of `DS` field
24d68a66c4SJames Wright 
25d68a66c4SJames Wright   @return An error code: 0 - success, otherwise - failure
26d68a66c4SJames Wright **/
27d68a66c4SJames Wright PetscErrorCode DMFieldToDSField(DM dm, DMLabel domain_label, PetscInt dm_field, PetscInt *ds_field) {
28d68a66c4SJames Wright   PetscDS         ds;
29d68a66c4SJames Wright   IS              field_is;
30d68a66c4SJames Wright   const PetscInt *fields;
31d68a66c4SJames Wright   PetscInt        num_fields;
32d68a66c4SJames Wright 
33d68a66c4SJames Wright   PetscFunctionBeginUser;
34d68a66c4SJames Wright   PetscCall(DMGetRegionDS(dm, domain_label, &field_is, &ds, NULL));
35d68a66c4SJames Wright   PetscCall(ISGetIndices(field_is, &fields));
36d68a66c4SJames Wright   PetscCall(ISGetSize(field_is, &num_fields));
37d68a66c4SJames Wright   for (PetscInt i = 0; i < num_fields; i++) {
38d68a66c4SJames Wright     if (dm_field == fields[i]) {
39d68a66c4SJames Wright       *ds_field = i;
40d68a66c4SJames Wright       break;
41d68a66c4SJames Wright     }
42d68a66c4SJames Wright   }
43d68a66c4SJames Wright   PetscCall(ISRestoreIndices(field_is, &fields));
44d68a66c4SJames Wright 
45*fe1e732eSJames Wright   PetscCheck(*ds_field != -1, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Could not find dm_field %" PetscInt_FMT " in DS", dm_field);
46d68a66c4SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
47d68a66c4SJames Wright }
48d68a66c4SJames Wright 
49bb85d312SJames Wright /**
50bb85d312SJames Wright   @brief Create `CeedElemRestriction` from `DMPlex`.
51bb85d312SJames Wright 
52bb85d312SJames Wright   Not collective across MPI processes.
53bb85d312SJames Wright 
54bb85d312SJames Wright   @param[in]   ceed          `Ceed` context
55bb85d312SJames Wright   @param[in]   dm            `DMPlex` holding mesh
56bb85d312SJames Wright   @param[in]   domain_label  `DMLabel` for `DMPlex` domain
57bb85d312SJames Wright   @param[in]   label_value   Stratum value
58bb85d312SJames Wright   @param[in]   height        Height of `DMPlex` topology
59bb85d312SJames Wright   @param[in]   dm_field      Index of `DMPlex` field
60bb85d312SJames Wright   @param[out]  restriction   `CeedElemRestriction` for `DMPlex`
61bb85d312SJames Wright 
62bb85d312SJames Wright   @return An error code: 0 - success, otherwise - failure
63bb85d312SJames Wright **/
64bb85d312SJames Wright PetscErrorCode DMPlexCeedElemRestrictionCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field,
65bb85d312SJames Wright                                                CeedElemRestriction *restriction) {
66bb85d312SJames Wright   PetscInt num_elem, elem_size, num_dof, num_comp, *restriction_offsets_petsc;
67bb85d312SJames Wright   CeedInt *restriction_offsets_ceed = NULL;
68bb85d312SJames Wright 
69bb85d312SJames Wright   PetscFunctionBeginUser;
70bb85d312SJames Wright   PetscCall(
71bb85d312SJames Wright       DMPlexGetLocalOffsets(dm, domain_label, label_value, height, dm_field, &num_elem, &elem_size, &num_comp, &num_dof, &restriction_offsets_petsc));
72bb85d312SJames Wright   PetscCall(IntArrayP2C(num_elem * elem_size, &restriction_offsets_petsc, &restriction_offsets_ceed));
73bb85d312SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreate(ceed, num_elem, elem_size, num_comp, 1, num_dof, CEED_MEM_HOST, CEED_COPY_VALUES,
74bb85d312SJames Wright                                                 restriction_offsets_ceed, restriction));
75bb85d312SJames Wright   PetscCall(PetscFree(restriction_offsets_ceed));
76bb85d312SJames Wright 
77bb85d312SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
78bb85d312SJames Wright }
79bb85d312SJames Wright 
80bb85d312SJames Wright /**
81bb85d312SJames Wright   @brief Create `CeedElemRestriction` from `DMPlex` domain for mesh coordinates.
82bb85d312SJames Wright 
83bb85d312SJames Wright   Not collective across MPI processes.
84bb85d312SJames Wright 
85bb85d312SJames Wright   @param[in]   ceed          `Ceed` context
86bb85d312SJames Wright   @param[in]   dm            `DMPlex` holding mesh
87bb85d312SJames Wright   @param[in]   domain_label  Label for `DMPlex` domain
88bb85d312SJames Wright   @param[in]   label_value   Stratum value
89bb85d312SJames Wright   @param[in]   height        Height of `DMPlex` topology
90bb85d312SJames Wright   @param[out]  restriction   `CeedElemRestriction` for mesh
91bb85d312SJames Wright 
92bb85d312SJames Wright   @return An error code: 0 - success, otherwise - failure
93bb85d312SJames Wright **/
94bb85d312SJames Wright PetscErrorCode DMPlexCeedElemRestrictionCoordinateCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height,
95bb85d312SJames Wright                                                          CeedElemRestriction *restriction) {
96bb85d312SJames Wright   DM dm_coord;
97bb85d312SJames Wright 
98bb85d312SJames Wright   PetscFunctionBeginUser;
99bb85d312SJames Wright   PetscCall(DMGetCellCoordinateDM(dm, &dm_coord));
100bb85d312SJames Wright   if (!dm_coord) {
101bb85d312SJames Wright     PetscCall(DMGetCoordinateDM(dm, &dm_coord));
102bb85d312SJames Wright   }
103bb85d312SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm_coord, domain_label, label_value, height, 0, restriction));
104bb85d312SJames Wright 
105bb85d312SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
106bb85d312SJames Wright }
107bb85d312SJames Wright 
108bb85d312SJames Wright /**
109bb85d312SJames Wright   @brief Create `CeedElemRestriction` from `DMPlex` domain for auxilury `QFunction` data.
110bb85d312SJames Wright 
111bb85d312SJames Wright   Not collective across MPI processes.
112bb85d312SJames Wright 
113bb85d312SJames Wright   @param[in]   ceed           `Ceed` context
114bb85d312SJames Wright   @param[in]   dm             `DMPlex` holding mesh
115bb85d312SJames Wright   @param[in]   domain_label   Label for `DMPlex` domain
116bb85d312SJames Wright   @param[in]   label_value    Stratum value
117bb85d312SJames Wright   @param[in]   height         Height of `DMPlex` topology
118bb85d312SJames Wright   @param[in]   q_data_size    Number of components for `QFunction` data
119bb85d312SJames Wright   @param[in]   is_collocated  Boolean flag indicating if the data is collocated on the nodes (`PETSC_TRUE`) on on quadrature points (`PETSC_FALSE`)
120bb85d312SJames Wright   @param[out]  restriction    Strided `CeedElemRestriction` for `QFunction` data
121bb85d312SJames Wright 
122bb85d312SJames Wright   @return An error code: 0 - success, otherwise - failure
123bb85d312SJames Wright **/
124bb85d312SJames Wright static PetscErrorCode DMPlexCeedElemRestrictionStridedCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height,
125bb85d312SJames Wright                                                              PetscInt q_data_size, PetscBool is_collocated, CeedElemRestriction *restriction) {
126bb85d312SJames Wright   PetscInt num_elem, num_qpts, dm_field = 0;
127bb85d312SJames Wright 
128bb85d312SJames Wright   PetscFunctionBeginUser;
129bb85d312SJames Wright   {  // Get number of elements
130bb85d312SJames Wright     PetscInt depth;
131bb85d312SJames Wright     DMLabel  depth_label;
132bb85d312SJames Wright     IS       point_is, depth_is;
133bb85d312SJames Wright 
134bb85d312SJames Wright     PetscCall(DMPlexGetDepth(dm, &depth));
135bb85d312SJames Wright     PetscCall(DMPlexGetDepthLabel(dm, &depth_label));
136bb85d312SJames Wright     PetscCall(DMLabelGetStratumIS(depth_label, depth - height, &depth_is));
137bb85d312SJames Wright     if (domain_label) {
138bb85d312SJames Wright       IS domain_is;
139bb85d312SJames Wright 
140bb85d312SJames Wright       PetscCall(DMLabelGetStratumIS(domain_label, label_value, &domain_is));
141bb85d312SJames Wright       if (domain_is) {
142bb85d312SJames Wright         PetscCall(ISIntersect(depth_is, domain_is, &point_is));
143bb85d312SJames Wright         PetscCall(ISDestroy(&domain_is));
144bb85d312SJames Wright       } else {
145bb85d312SJames Wright         point_is = NULL;
146bb85d312SJames Wright       }
147bb85d312SJames Wright       PetscCall(ISDestroy(&depth_is));
148bb85d312SJames Wright     } else {
149bb85d312SJames Wright       point_is = depth_is;
150bb85d312SJames Wright     }
151bb85d312SJames Wright     if (point_is) {
152bb85d312SJames Wright       PetscCall(ISGetLocalSize(point_is, &num_elem));
153bb85d312SJames Wright     } else {
154bb85d312SJames Wright       num_elem = 0;
155bb85d312SJames Wright     }
156bb85d312SJames Wright     PetscCall(ISDestroy(&point_is));
157bb85d312SJames Wright   }
158bb85d312SJames Wright 
159bb85d312SJames Wright   {  // Get number of quadrature points
160bb85d312SJames Wright     PetscDS  ds;
161bb85d312SJames Wright     PetscFE  fe;
162bb85d312SJames Wright     PetscInt ds_field = -1;
163bb85d312SJames Wright 
164bb85d312SJames Wright     PetscCall(DMGetRegionDS(dm, domain_label, NULL, &ds, NULL));
165bb85d312SJames Wright     PetscCall(DMFieldToDSField(dm, domain_label, dm_field, &ds_field));
166bb85d312SJames Wright     PetscCall(PetscDSGetDiscretization(ds, ds_field, (PetscObject *)&fe));
167bb85d312SJames Wright     PetscCall(PetscFEGetHeightSubspace(fe, height, &fe));
168bb85d312SJames Wright     if (is_collocated) {
169bb85d312SJames Wright       PetscDualSpace dual_space;
170bb85d312SJames Wright       PetscInt       num_dual_basis_vectors, dim, num_comp;
171bb85d312SJames Wright 
172bb85d312SJames Wright       PetscCall(PetscFEGetSpatialDimension(fe, &dim));
173bb85d312SJames Wright       PetscCall(PetscFEGetNumComponents(fe, &num_comp));
174bb85d312SJames Wright       PetscCall(PetscFEGetDualSpace(fe, &dual_space));
175bb85d312SJames Wright       PetscCall(PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors));
176bb85d312SJames Wright       num_qpts = num_dual_basis_vectors / num_comp;
177bb85d312SJames Wright     } else {
178bb85d312SJames Wright       PetscQuadrature quadrature;
179bb85d312SJames Wright 
180bb85d312SJames Wright       PetscCall(DMGetRegionDS(dm, domain_label, NULL, &ds, NULL));
181bb85d312SJames Wright       PetscCall(DMFieldToDSField(dm, domain_label, dm_field, &ds_field));
182bb85d312SJames Wright       PetscCall(PetscDSGetDiscretization(ds, ds_field, (PetscObject *)&fe));
183bb85d312SJames Wright       PetscCall(PetscFEGetHeightSubspace(fe, height, &fe));
184bb85d312SJames Wright       PetscCall(PetscFEGetQuadrature(fe, &quadrature));
185bb85d312SJames Wright       PetscCall(PetscQuadratureGetData(quadrature, NULL, NULL, &num_qpts, NULL, NULL));
186bb85d312SJames Wright     }
187bb85d312SJames Wright   }
188bb85d312SJames Wright 
189bb85d312SJames Wright   // Create the restriction
190bb85d312SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateStrided(ceed, num_elem, num_qpts, q_data_size, num_elem * num_qpts * q_data_size, CEED_STRIDES_BACKEND,
191bb85d312SJames Wright                                                        restriction));
192bb85d312SJames Wright 
193bb85d312SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
194bb85d312SJames Wright }
195bb85d312SJames Wright 
196bb85d312SJames Wright /**
197bb85d312SJames Wright   @brief Create `CeedElemRestriction` from `DMPlex` domain for auxilury `QFunction` data.
198bb85d312SJames Wright 
199bb85d312SJames Wright   Not collective across MPI processes.
200bb85d312SJames Wright 
201bb85d312SJames Wright   @param[in]   ceed           `Ceed` context
202bb85d312SJames Wright   @param[in]   dm             `DMPlex` holding mesh
203bb85d312SJames Wright   @param[in]   domain_label   Label for `DMPlex` domain
204bb85d312SJames Wright   @param[in]   label_value    Stratum value
205bb85d312SJames Wright   @param[in]   height         Height of `DMPlex` topology
206bb85d312SJames Wright   @param[in]   q_data_size    Number of components for `QFunction` data
207bb85d312SJames Wright   @param[out]  restriction    Strided `CeedElemRestriction` for `QFunction` data
208bb85d312SJames Wright 
209bb85d312SJames Wright   @return An error code: 0 - success, otherwise - failure
210bb85d312SJames Wright **/
211bb85d312SJames Wright PetscErrorCode DMPlexCeedElemRestrictionQDataCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height,
212bb85d312SJames Wright                                                     PetscInt q_data_size, CeedElemRestriction *restriction) {
213bb85d312SJames Wright   PetscFunctionBeginUser;
214bb85d312SJames Wright   PetscCall(DMPlexCeedElemRestrictionStridedCreate(ceed, dm, domain_label, label_value, height, q_data_size, PETSC_FALSE, restriction));
215bb85d312SJames Wright 
216bb85d312SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
217bb85d312SJames Wright }
218bb85d312SJames Wright 
219bb85d312SJames Wright /**
220bb85d312SJames Wright   @brief Create `CeedElemRestriction` from `DMPlex` domain for nodally collocated auxilury `QFunction` data.
221bb85d312SJames Wright 
222bb85d312SJames Wright   Not collective across MPI processes.
223bb85d312SJames Wright 
224bb85d312SJames Wright   @param[in]   ceed           `Ceed` context
225bb85d312SJames Wright   @param[in]   dm             `DMPlex` holding mesh
226bb85d312SJames Wright   @param[in]   domain_label   Label for `DMPlex` domain
227bb85d312SJames Wright   @param[in]   label_value    Stratum value
228bb85d312SJames Wright   @param[in]   height         Height of `DMPlex` topology
229bb85d312SJames Wright   @param[in]   q_data_size    Number of components for `QFunction` data
230bb85d312SJames Wright   @param[out]  restriction    Strided `CeedElemRestriction` for `QFunction` data
231bb85d312SJames Wright 
232bb85d312SJames Wright   @return An error code: 0 - success, otherwise - failure
233bb85d312SJames Wright **/
234bb85d312SJames Wright PetscErrorCode DMPlexCeedElemRestrictionCollocatedCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height,
235bb85d312SJames Wright                                                          PetscInt q_data_size, CeedElemRestriction *restriction) {
236bb85d312SJames Wright   PetscFunctionBeginUser;
237bb85d312SJames Wright   PetscCall(DMPlexCeedElemRestrictionStridedCreate(ceed, dm, domain_label, label_value, height, q_data_size, PETSC_TRUE, restriction));
238bb85d312SJames Wright 
239bb85d312SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
240bb85d312SJames Wright }
241bb85d312SJames Wright 
242d68a66c4SJames Wright // -----------------------------------------------------------------------------
243d68a66c4SJames Wright // Utility function - convert from DMPolytopeType to CeedElemTopology
244d68a66c4SJames Wright // -----------------------------------------------------------------------------
245d68a66c4SJames Wright static CeedElemTopology ElemTopologyP2C(DMPolytopeType cell_type) {
246d68a66c4SJames Wright   switch (cell_type) {
247d68a66c4SJames Wright     case DM_POLYTOPE_TRIANGLE:
248d68a66c4SJames Wright       return CEED_TOPOLOGY_TRIANGLE;
249d68a66c4SJames Wright     case DM_POLYTOPE_QUADRILATERAL:
250d68a66c4SJames Wright       return CEED_TOPOLOGY_QUAD;
251d68a66c4SJames Wright     case DM_POLYTOPE_TETRAHEDRON:
252d68a66c4SJames Wright       return CEED_TOPOLOGY_TET;
253d68a66c4SJames Wright     case DM_POLYTOPE_HEXAHEDRON:
254d68a66c4SJames Wright       return CEED_TOPOLOGY_HEX;
255d68a66c4SJames Wright     default:
256d68a66c4SJames Wright       return 0;
257d68a66c4SJames Wright   }
258d68a66c4SJames Wright }
259d68a66c4SJames Wright 
260d68a66c4SJames Wright // -----------------------------------------------------------------------------
261d68a66c4SJames Wright // Create libCEED Basis from PetscTabulation
262d68a66c4SJames Wright // -----------------------------------------------------------------------------
263d68a66c4SJames Wright PetscErrorCode BasisCreateFromTabulation(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt face, PetscFE fe,
264d68a66c4SJames Wright                                          PetscTabulation basis_tabulation, PetscQuadrature quadrature, CeedBasis *basis) {
265d68a66c4SJames Wright   PetscInt           first_point;
266d68a66c4SJames Wright   PetscInt           ids[1] = {label_value};
267d68a66c4SJames Wright   DMLabel            depth_label;
268d68a66c4SJames Wright   DMPolytopeType     cell_type;
269d68a66c4SJames Wright   CeedElemTopology   elem_topo;
270d68a66c4SJames Wright   PetscScalar       *q_points, *interp, *grad;
271d68a66c4SJames Wright   const PetscScalar *q_weights;
272d68a66c4SJames Wright   PetscDualSpace     dual_space;
273d68a66c4SJames Wright   PetscInt           num_dual_basis_vectors;
274d68a66c4SJames Wright   PetscInt           dim, num_comp, P, Q;
275d68a66c4SJames Wright 
276d68a66c4SJames Wright   PetscFunctionBeginUser;
277d68a66c4SJames Wright   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
278d68a66c4SJames Wright   PetscCall(PetscFEGetNumComponents(fe, &num_comp));
279d68a66c4SJames Wright   PetscCall(PetscFEGetDualSpace(fe, &dual_space));
280d68a66c4SJames Wright   PetscCall(PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors));
281d68a66c4SJames Wright   P = num_dual_basis_vectors / num_comp;
282d68a66c4SJames Wright 
283d68a66c4SJames Wright   // Use depth label if no domain label present
284d68a66c4SJames Wright   if (!domain_label) {
285d68a66c4SJames Wright     PetscInt depth;
286d68a66c4SJames Wright 
287d68a66c4SJames Wright     PetscCall(DMPlexGetDepth(dm, &depth));
288d68a66c4SJames Wright     PetscCall(DMPlexGetDepthLabel(dm, &depth_label));
289d68a66c4SJames Wright     ids[0] = depth - height;
290d68a66c4SJames Wright   }
291d68a66c4SJames Wright 
292d68a66c4SJames Wright   // Get cell interp, grad, and quadrature data
293d68a66c4SJames Wright   PetscCall(DMGetFirstLabeledPoint(dm, dm, domain_label ? domain_label : depth_label, 1, ids, height, &first_point, NULL));
294d68a66c4SJames Wright   PetscCall(DMPlexGetCellType(dm, first_point, &cell_type));
295d68a66c4SJames Wright   elem_topo = ElemTopologyP2C(cell_type);
296*fe1e732eSJames Wright   PetscCheck(elem_topo, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "DMPlex topology not supported");
297d68a66c4SJames Wright   {
298d68a66c4SJames Wright     size_t             q_points_size;
299d68a66c4SJames Wright     const PetscScalar *q_points_petsc;
300d68a66c4SJames Wright     PetscInt           q_dim;
301d68a66c4SJames Wright 
302d68a66c4SJames Wright     PetscCall(PetscQuadratureGetData(quadrature, &q_dim, NULL, &Q, &q_points_petsc, &q_weights));
303d68a66c4SJames Wright     q_points_size = Q * dim * sizeof(CeedScalar);
304d68a66c4SJames Wright     PetscCall(PetscCalloc(q_points_size, &q_points));
305d68a66c4SJames Wright     for (PetscInt q = 0; q < Q; q++) {
306d68a66c4SJames Wright       for (PetscInt d = 0; d < q_dim; d++) q_points[q * dim + d] = q_points_petsc[q * q_dim + d];
307d68a66c4SJames Wright     }
308d68a66c4SJames Wright   }
309d68a66c4SJames Wright 
310d68a66c4SJames Wright   {  // Convert to libCEED orientation
311d68a66c4SJames Wright     PetscBool       is_simplex  = PETSC_FALSE;
312d68a66c4SJames Wright     IS              permutation = NULL;
313d68a66c4SJames Wright     const PetscInt *permutation_indices;
314d68a66c4SJames Wright 
315d68a66c4SJames Wright     PetscCall(DMPlexIsSimplex(dm, &is_simplex));
316d68a66c4SJames Wright     if (!is_simplex) {
317d68a66c4SJames Wright       PetscSection section;
318d68a66c4SJames Wright 
319d68a66c4SJames Wright       // -- Get permutation
320d68a66c4SJames Wright       PetscCall(DMGetLocalSection(dm, &section));
321d68a66c4SJames Wright       PetscCall(PetscSectionGetClosurePermutation(section, (PetscObject)dm, dim, num_comp * P, &permutation));
322d68a66c4SJames Wright       PetscCall(ISGetIndices(permutation, &permutation_indices));
323d68a66c4SJames Wright     }
324d68a66c4SJames Wright 
325d68a66c4SJames Wright     // -- Copy interp, grad matrices
326d68a66c4SJames Wright     PetscCall(PetscCalloc(P * Q * sizeof(CeedScalar), &interp));
327d68a66c4SJames Wright     PetscCall(PetscCalloc(P * Q * dim * sizeof(CeedScalar), &grad));
328d68a66c4SJames Wright     const CeedInt c = 0;
329d68a66c4SJames Wright     for (CeedInt q = 0; q < Q; q++) {
330d68a66c4SJames Wright       for (CeedInt p_ceed = 0; p_ceed < P; p_ceed++) {
331d68a66c4SJames Wright         CeedInt p_petsc = is_simplex ? (p_ceed * num_comp) : permutation_indices[p_ceed * num_comp];
332d68a66c4SJames Wright 
333d68a66c4SJames Wright         interp[q * P + p_ceed] = basis_tabulation->T[0][((face * Q + q) * P * num_comp + p_petsc) * num_comp + c];
334d68a66c4SJames Wright         for (CeedInt d = 0; d < dim; d++) {
335d68a66c4SJames Wright           grad[(d * Q + q) * P + p_ceed] = basis_tabulation->T[1][(((face * Q + q) * P * num_comp + p_petsc) * num_comp + c) * dim + d];
336d68a66c4SJames Wright         }
337d68a66c4SJames Wright       }
338d68a66c4SJames Wright     }
339d68a66c4SJames Wright 
340d68a66c4SJames Wright     // -- Cleanup
341d68a66c4SJames Wright     if (permutation) PetscCall(ISRestoreIndices(permutation, &permutation_indices));
342d68a66c4SJames Wright     PetscCall(ISDestroy(&permutation));
343d68a66c4SJames Wright   }
344d68a66c4SJames Wright 
345d68a66c4SJames Wright   // Finally, create libCEED basis
346d68a66c4SJames Wright   PetscCallCeed(ceed, CeedBasisCreateH1(ceed, elem_topo, num_comp, P, Q, interp, grad, q_points, q_weights, basis));
347d68a66c4SJames Wright   PetscCall(PetscFree(q_points));
348d68a66c4SJames Wright   PetscCall(PetscFree(interp));
349d68a66c4SJames Wright   PetscCall(PetscFree(grad));
350d68a66c4SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
351d68a66c4SJames Wright }
352d68a66c4SJames Wright 
353d68a66c4SJames Wright // -----------------------------------------------------------------------------
354d68a66c4SJames Wright // Get CEED Basis from DMPlex
355d68a66c4SJames Wright // -----------------------------------------------------------------------------
356d68a66c4SJames Wright PetscErrorCode CreateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label, CeedInt label_value, CeedInt height, CeedInt dm_field, CeedBasis *basis) {
357d68a66c4SJames Wright   PetscDS         ds;
358d68a66c4SJames Wright   PetscFE         fe;
359d68a66c4SJames Wright   PetscQuadrature quadrature;
360d68a66c4SJames Wright   PetscBool       is_simplex = PETSC_TRUE;
361d68a66c4SJames Wright   PetscInt        ds_field   = -1;
362d68a66c4SJames Wright 
363d68a66c4SJames Wright   PetscFunctionBeginUser;
364d68a66c4SJames Wright   // Get element information
365d68a66c4SJames Wright   PetscCall(DMGetRegionDS(dm, domain_label, NULL, &ds, NULL));
366d68a66c4SJames Wright   PetscCall(DMFieldToDSField(dm, domain_label, dm_field, &ds_field));
367d68a66c4SJames Wright   PetscCall(PetscDSGetDiscretization(ds, ds_field, (PetscObject *)&fe));
368d68a66c4SJames Wright   PetscCall(PetscFEGetHeightSubspace(fe, height, &fe));
369d68a66c4SJames Wright   PetscCall(PetscFEGetQuadrature(fe, &quadrature));
370d68a66c4SJames Wright 
371d68a66c4SJames Wright   // Check if simplex or tensor-product mesh
372d68a66c4SJames Wright   PetscCall(DMPlexIsSimplex(dm, &is_simplex));
373d68a66c4SJames Wright 
374d68a66c4SJames Wright   // Build libCEED basis
375d68a66c4SJames Wright   if (is_simplex) {
376d68a66c4SJames Wright     PetscTabulation basis_tabulation;
377d68a66c4SJames Wright     PetscInt        num_derivatives = 1, face = 0;
378d68a66c4SJames Wright 
379d68a66c4SJames Wright     PetscCall(PetscFEGetCellTabulation(fe, num_derivatives, &basis_tabulation));
380d68a66c4SJames Wright     PetscCall(BasisCreateFromTabulation(ceed, dm, domain_label, label_value, height, face, fe, basis_tabulation, quadrature, basis));
381d68a66c4SJames Wright   } else {
382d68a66c4SJames Wright     PetscDualSpace dual_space;
383d68a66c4SJames Wright     PetscInt       num_dual_basis_vectors;
384d68a66c4SJames Wright     PetscInt       dim, num_comp, P, Q;
385d68a66c4SJames Wright 
386d68a66c4SJames Wright     PetscCall(PetscFEGetSpatialDimension(fe, &dim));
387d68a66c4SJames Wright     PetscCall(PetscFEGetNumComponents(fe, &num_comp));
388d68a66c4SJames Wright     PetscCall(PetscFEGetDualSpace(fe, &dual_space));
389d68a66c4SJames Wright     PetscCall(PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors));
390d68a66c4SJames Wright     P = num_dual_basis_vectors / num_comp;
391d68a66c4SJames Wright     PetscCall(PetscQuadratureGetData(quadrature, NULL, NULL, &Q, NULL, NULL));
392d68a66c4SJames Wright 
393d68a66c4SJames Wright     CeedInt P_1d = (CeedInt)round(pow(P, 1.0 / dim));
394d68a66c4SJames Wright     CeedInt Q_1d = (CeedInt)round(pow(Q, 1.0 / dim));
395d68a66c4SJames Wright 
396d68a66c4SJames Wright     PetscCallCeed(ceed, CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp, P_1d, Q_1d, CEED_GAUSS, basis));
397d68a66c4SJames Wright   }
398d68a66c4SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
399d68a66c4SJames Wright }
400d68a66c4SJames Wright 
401d68a66c4SJames Wright /**
402d68a66c4SJames Wright   @brief Setup `DM` with FE space of appropriate degree
403d68a66c4SJames Wright 
404d68a66c4SJames Wright   Must be followed by `DMSetupByOrderEnd_FEM`
405d68a66c4SJames Wright 
406d68a66c4SJames Wright   @param[in]   setup_faces     Flag to setup face geometry
407d68a66c4SJames Wright   @param[in]   setup_coords    Flag to setup coordinate spaces
408d68a66c4SJames Wright   @param[in]   degree          Polynomial orders of field
409a0b9cdb5SJames Wright   @param[in]   coord_order     Polynomial order of coordinate basis, or `PETSC_DECIDE` for default
410a0b9cdb5SJames Wright   @param[in]   q_extra         Additional quadrature order
411d68a66c4SJames Wright   @param[in]   num_fields      Number of fields in solution vector
412d68a66c4SJames Wright   @param[in]   field_sizes     Array of number of components for each field
413d68a66c4SJames Wright   @param[out]  dm              `DM` to setup
414d68a66c4SJames Wright 
415d68a66c4SJames Wright   @return An error code: 0 - success, otherwise - failure
416d68a66c4SJames Wright **/
417d68a66c4SJames Wright PetscErrorCode DMSetupByOrderBegin_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra,
41885df0564SJames Wright                                        PetscInt num_fields, const PetscInt *field_sizes, DM dm) {
419d68a66c4SJames Wright   PetscInt  dim, q_order = degree + q_extra;
420d68a66c4SJames Wright   PetscBool is_simplex = PETSC_TRUE;
421d68a66c4SJames Wright   PetscFE   fe;
422d68a66c4SJames Wright   MPI_Comm  comm = PetscObjectComm((PetscObject)dm);
423d68a66c4SJames Wright 
424d68a66c4SJames Wright   PetscFunctionBeginUser;
425d68a66c4SJames Wright   PetscCall(DMPlexIsSimplex(dm, &is_simplex));
426d68a66c4SJames Wright 
427d68a66c4SJames Wright   // Setup DM
428d68a66c4SJames Wright   PetscCall(DMGetDimension(dm, &dim));
429d68a66c4SJames Wright   for (PetscInt i = 0; i < num_fields; i++) {
430d68a66c4SJames Wright     PetscFE  fe_face;
431d68a66c4SJames Wright     PetscInt q_order = degree + q_extra;
432d68a66c4SJames Wright 
433d68a66c4SJames Wright     PetscCall(PetscFECreateLagrange(comm, dim, field_sizes[i], is_simplex, degree, q_order, &fe));
434d68a66c4SJames Wright     if (setup_faces) PetscCall(PetscFEGetHeightSubspace(fe, 1, &fe_face));
435d68a66c4SJames Wright     PetscCall(DMAddField(dm, NULL, (PetscObject)fe));
436d68a66c4SJames Wright     PetscCall(PetscFEDestroy(&fe));
437d68a66c4SJames Wright   }
438d68a66c4SJames Wright   PetscCall(DMCreateDS(dm));
439d68a66c4SJames Wright 
440d68a66c4SJames Wright   // Project coordinates to enrich quadrature space
441d68a66c4SJames Wright   if (setup_coords) {
442d68a66c4SJames Wright     DM             dm_coord;
443d68a66c4SJames Wright     PetscDS        ds_coord;
444d68a66c4SJames Wright     PetscFE        fe_coord_current, fe_coord_new, fe_coord_face_new;
445d68a66c4SJames Wright     PetscDualSpace fe_coord_dual_space;
446d68a66c4SJames Wright     PetscInt       fe_coord_order, num_comp_coord;
447d68a66c4SJames Wright 
448d68a66c4SJames Wright     PetscCall(DMGetCoordinateDM(dm, &dm_coord));
449d68a66c4SJames Wright     PetscCall(DMGetCoordinateDim(dm, &num_comp_coord));
450d68a66c4SJames Wright     PetscCall(DMGetRegionDS(dm_coord, NULL, NULL, &ds_coord, NULL));
451d68a66c4SJames Wright     PetscCall(PetscDSGetDiscretization(ds_coord, 0, (PetscObject *)&fe_coord_current));
452d68a66c4SJames Wright     PetscCall(PetscFEGetDualSpace(fe_coord_current, &fe_coord_dual_space));
453d68a66c4SJames Wright     PetscCall(PetscDualSpaceGetOrder(fe_coord_dual_space, &fe_coord_order));
454d68a66c4SJames Wright 
455d68a66c4SJames Wright     // Create FE for coordinates
456a0b9cdb5SJames Wright     if (coord_order != PETSC_DECIDE) fe_coord_order = coord_order;
457d68a66c4SJames Wright     PetscCall(PetscFECreateLagrange(comm, dim, num_comp_coord, is_simplex, fe_coord_order, q_order, &fe_coord_new));
458d68a66c4SJames Wright     if (setup_faces) PetscCall(PetscFEGetHeightSubspace(fe_coord_new, 1, &fe_coord_face_new));
45949a40c8aSKenneth E. Jansen     PetscCall(DMSetCoordinateDisc(dm, fe_coord_new, PETSC_TRUE));
460736e557eSJames Wright     PetscCall(DMLocalizeCoordinates(dm));  // Update CellCoordinateDM with projected coordinates
461d68a66c4SJames Wright     PetscCall(PetscFEDestroy(&fe_coord_new));
462d68a66c4SJames Wright   }
463d68a66c4SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
464d68a66c4SJames Wright }
465d68a66c4SJames Wright 
466d68a66c4SJames Wright /**
467d68a66c4SJames Wright   @brief Finish setting up `DM`
468d68a66c4SJames Wright 
469d68a66c4SJames Wright   Must be called after `DMSetupByOrderBegin_FEM` and all strong BCs have be declared via `DMAddBoundaries`.
470d68a66c4SJames Wright 
471d68a66c4SJames Wright   @param[in]   setup_coords    Flag to setup coordinate spaces
472d68a66c4SJames Wright   @param[out]  dm              `DM` to setup
473d68a66c4SJames Wright 
474d68a66c4SJames Wright   @return An error code: 0 - success, otherwise - failure
475d68a66c4SJames Wright **/
476d68a66c4SJames Wright PetscErrorCode DMSetupByOrderEnd_FEM(PetscBool setup_coords, DM dm) {
477d68a66c4SJames Wright   PetscBool is_simplex;
478d68a66c4SJames Wright 
479f17d818dSJames Wright   PetscFunctionBeginUser;
480d68a66c4SJames Wright   PetscCall(DMPlexIsSimplex(dm, &is_simplex));
481d68a66c4SJames Wright   // Set tensor permutation if needed
482d68a66c4SJames Wright   if (!is_simplex) {
483d68a66c4SJames Wright     PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL));
484d68a66c4SJames Wright     if (setup_coords) {
485d68a66c4SJames Wright       DM dm_coord;
486d68a66c4SJames Wright 
487d68a66c4SJames Wright       PetscCall(DMGetCoordinateDM(dm, &dm_coord));
488d68a66c4SJames Wright       PetscCall(DMPlexSetClosurePermutationTensor(dm_coord, PETSC_DETERMINE, NULL));
489736e557eSJames Wright       PetscCall(DMGetCellCoordinateDM(dm, &dm_coord));
490736e557eSJames Wright       if (dm_coord) PetscCall(DMPlexSetClosurePermutationTensor(dm_coord, PETSC_DETERMINE, NULL));
491d68a66c4SJames Wright     }
492d68a66c4SJames Wright   }
493d68a66c4SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
494d68a66c4SJames Wright }
495d68a66c4SJames Wright 
496d68a66c4SJames Wright /**
497d68a66c4SJames Wright   @brief Setup `DM` with FE space of appropriate degree with no boundary conditions
498d68a66c4SJames Wright 
499d68a66c4SJames Wright   Calls `DMSetupByOrderBegin_FEM` and `DMSetupByOrderEnd_FEM` successively
500d68a66c4SJames Wright 
501d68a66c4SJames Wright   @param[in]   setup_faces     Flag to setup face geometry
502d68a66c4SJames Wright   @param[in]   setup_coords    Flag to setup coordinate spaces
503d68a66c4SJames Wright   @param[in]   degree          Polynomial orders of field
504a0b9cdb5SJames Wright   @param[in]   coord_order     Polynomial order of coordinate basis, or `PETSC_DECIDE` for default
505a0b9cdb5SJames Wright   @param[in]   q_extra         Additional quadrature order
506d68a66c4SJames Wright   @param[in]   num_fields      Number of fields in solution vector
507d68a66c4SJames Wright   @param[in]   field_sizes     Array of number of components for each field
508d68a66c4SJames Wright   @param[out]  dm              `DM` to setup
509d68a66c4SJames Wright 
510d68a66c4SJames Wright   @return An error code: 0 - success, otherwise - failure
511d68a66c4SJames Wright **/
512d68a66c4SJames Wright PetscErrorCode DMSetupByOrder_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra,
51385df0564SJames Wright                                   PetscInt num_fields, const PetscInt *field_sizes, DM dm) {
514d68a66c4SJames Wright   PetscFunctionBeginUser;
515d68a66c4SJames Wright   PetscCall(DMSetupByOrderBegin_FEM(setup_faces, setup_coords, degree, coord_order, q_extra, num_fields, field_sizes, dm));
516d68a66c4SJames Wright   PetscCall(DMSetupByOrderEnd_FEM(setup_coords, dm));
517d68a66c4SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
518d68a66c4SJames Wright }
519