xref: /libCEED/examples/fluids/src/dm_utils.c (revision f17d818ddb5f0fdb5c67183d33cc31a0016906ee)
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 // Utility function to create local CEED restriction
18d68a66c4SJames Wright PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, CeedInt label_value, PetscInt dm_field,
19d68a66c4SJames Wright                                          CeedElemRestriction *elem_restr) {
20d68a66c4SJames Wright   PetscInt num_elem, elem_size, num_dof, num_comp, *elem_restr_offsets_petsc;
21d68a66c4SJames Wright   CeedInt *elem_restr_offsets_ceed;
22d68a66c4SJames Wright 
23d68a66c4SJames Wright   PetscFunctionBeginUser;
24d68a66c4SJames Wright   PetscCall(
25d68a66c4SJames Wright       DMPlexGetLocalOffsets(dm, domain_label, label_value, height, dm_field, &num_elem, &elem_size, &num_comp, &num_dof, &elem_restr_offsets_petsc));
26d68a66c4SJames Wright 
27d68a66c4SJames Wright   PetscCall(IntArrayP2C(num_elem * elem_size, &elem_restr_offsets_petsc, &elem_restr_offsets_ceed));
28d68a66c4SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreate(ceed, num_elem, elem_size, num_comp, 1, num_dof, CEED_MEM_HOST, CEED_COPY_VALUES,
29d68a66c4SJames Wright                                                 elem_restr_offsets_ceed, elem_restr));
30d68a66c4SJames Wright   PetscCall(PetscFree(elem_restr_offsets_ceed));
31d68a66c4SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
32d68a66c4SJames Wright }
33d68a66c4SJames Wright 
34d68a66c4SJames Wright // Utility function to get Ceed Restriction for each domain
35d68a66c4SJames Wright PetscErrorCode GetRestrictionForDomain(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, PetscInt label_value, PetscInt dm_field, CeedInt Q,
36d68a66c4SJames Wright                                        CeedInt q_data_size, CeedElemRestriction *elem_restr_q, CeedElemRestriction *elem_restr_x,
37d68a66c4SJames Wright                                        CeedElemRestriction *elem_restr_qd_i) {
38d68a66c4SJames Wright   DM                  dm_coord;
39d68a66c4SJames Wright   CeedInt             loc_num_elem;
40d68a66c4SJames Wright   PetscInt            dim;
41d68a66c4SJames Wright   CeedElemRestriction elem_restr_tmp;
42d68a66c4SJames Wright 
43*f17d818dSJames Wright   PetscFunctionBeginUser;
44d68a66c4SJames Wright   PetscCall(DMGetDimension(dm, &dim));
45d68a66c4SJames Wright   dim -= height;
46d68a66c4SJames Wright   PetscCall(CreateRestrictionFromPlex(ceed, dm, height, domain_label, label_value, dm_field, &elem_restr_tmp));
47d68a66c4SJames Wright   if (elem_restr_q) *elem_restr_q = elem_restr_tmp;
48d68a66c4SJames Wright   if (elem_restr_x) {
49d68a66c4SJames Wright     PetscCall(DMGetCellCoordinateDM(dm, &dm_coord));
50d68a66c4SJames Wright     if (!dm_coord) {
51d68a66c4SJames Wright       PetscCall(DMGetCoordinateDM(dm, &dm_coord));
52d68a66c4SJames Wright     }
53d68a66c4SJames Wright     PetscCall(DMPlexSetClosurePermutationTensor(dm_coord, PETSC_DETERMINE, NULL));
54d68a66c4SJames Wright     PetscCall(CreateRestrictionFromPlex(ceed, dm_coord, height, domain_label, label_value, dm_field, elem_restr_x));
55d68a66c4SJames Wright   }
56d68a66c4SJames Wright   if (elem_restr_qd_i) {
57d68a66c4SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionGetNumElements(elem_restr_tmp, &loc_num_elem));
58d68a66c4SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateStrided(ceed, loc_num_elem, Q, q_data_size, q_data_size * loc_num_elem * Q, CEED_STRIDES_BACKEND,
59d68a66c4SJames Wright                                                          elem_restr_qd_i));
60d68a66c4SJames Wright   }
61d68a66c4SJames Wright   if (!elem_restr_q) PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_tmp));
62d68a66c4SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
63d68a66c4SJames Wright }
64d68a66c4SJames Wright 
65d68a66c4SJames Wright /**
66d68a66c4SJames Wright   @brief Convert `DM` field to `DS` field.
67d68a66c4SJames Wright 
68d68a66c4SJames Wright   @param[in]   dm            `DM` holding mesh
69d68a66c4SJames Wright   @param[in]   domain_label  Label for `DM` domain
70d68a66c4SJames Wright   @param[in]   dm_field      Index of `DM` field
71d68a66c4SJames Wright   @param[out]  ds_field      Index of `DS` field
72d68a66c4SJames Wright 
73d68a66c4SJames Wright   @return An error code: 0 - success, otherwise - failure
74d68a66c4SJames Wright **/
75d68a66c4SJames Wright PetscErrorCode DMFieldToDSField(DM dm, DMLabel domain_label, PetscInt dm_field, PetscInt *ds_field) {
76d68a66c4SJames Wright   PetscDS         ds;
77d68a66c4SJames Wright   IS              field_is;
78d68a66c4SJames Wright   const PetscInt *fields;
79d68a66c4SJames Wright   PetscInt        num_fields;
80d68a66c4SJames Wright 
81d68a66c4SJames Wright   PetscFunctionBeginUser;
82d68a66c4SJames Wright   PetscCall(DMGetRegionDS(dm, domain_label, &field_is, &ds, NULL));
83d68a66c4SJames Wright   PetscCall(ISGetIndices(field_is, &fields));
84d68a66c4SJames Wright   PetscCall(ISGetSize(field_is, &num_fields));
85d68a66c4SJames Wright   for (PetscInt i = 0; i < num_fields; i++) {
86d68a66c4SJames Wright     if (dm_field == fields[i]) {
87d68a66c4SJames Wright       *ds_field = i;
88d68a66c4SJames Wright       break;
89d68a66c4SJames Wright     }
90d68a66c4SJames Wright   }
91d68a66c4SJames Wright   PetscCall(ISRestoreIndices(field_is, &fields));
92d68a66c4SJames Wright 
93d68a66c4SJames Wright   if (*ds_field == -1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Could not find dm_field %" PetscInt_FMT " in DS", dm_field);
94d68a66c4SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
95d68a66c4SJames Wright }
96d68a66c4SJames Wright 
97d68a66c4SJames Wright // -----------------------------------------------------------------------------
98d68a66c4SJames Wright // Utility function - convert from DMPolytopeType to CeedElemTopology
99d68a66c4SJames Wright // -----------------------------------------------------------------------------
100d68a66c4SJames Wright static CeedElemTopology ElemTopologyP2C(DMPolytopeType cell_type) {
101d68a66c4SJames Wright   switch (cell_type) {
102d68a66c4SJames Wright     case DM_POLYTOPE_TRIANGLE:
103d68a66c4SJames Wright       return CEED_TOPOLOGY_TRIANGLE;
104d68a66c4SJames Wright     case DM_POLYTOPE_QUADRILATERAL:
105d68a66c4SJames Wright       return CEED_TOPOLOGY_QUAD;
106d68a66c4SJames Wright     case DM_POLYTOPE_TETRAHEDRON:
107d68a66c4SJames Wright       return CEED_TOPOLOGY_TET;
108d68a66c4SJames Wright     case DM_POLYTOPE_HEXAHEDRON:
109d68a66c4SJames Wright       return CEED_TOPOLOGY_HEX;
110d68a66c4SJames Wright     default:
111d68a66c4SJames Wright       return 0;
112d68a66c4SJames Wright   }
113d68a66c4SJames Wright }
114d68a66c4SJames Wright 
115d68a66c4SJames Wright // -----------------------------------------------------------------------------
116d68a66c4SJames Wright // Create libCEED Basis from PetscTabulation
117d68a66c4SJames Wright // -----------------------------------------------------------------------------
118d68a66c4SJames Wright PetscErrorCode BasisCreateFromTabulation(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt face, PetscFE fe,
119d68a66c4SJames Wright                                          PetscTabulation basis_tabulation, PetscQuadrature quadrature, CeedBasis *basis) {
120d68a66c4SJames Wright   PetscInt           first_point;
121d68a66c4SJames Wright   PetscInt           ids[1] = {label_value};
122d68a66c4SJames Wright   DMLabel            depth_label;
123d68a66c4SJames Wright   DMPolytopeType     cell_type;
124d68a66c4SJames Wright   CeedElemTopology   elem_topo;
125d68a66c4SJames Wright   PetscScalar       *q_points, *interp, *grad;
126d68a66c4SJames Wright   const PetscScalar *q_weights;
127d68a66c4SJames Wright   PetscDualSpace     dual_space;
128d68a66c4SJames Wright   PetscInt           num_dual_basis_vectors;
129d68a66c4SJames Wright   PetscInt           dim, num_comp, P, Q;
130d68a66c4SJames Wright 
131d68a66c4SJames Wright   PetscFunctionBeginUser;
132d68a66c4SJames Wright   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
133d68a66c4SJames Wright   PetscCall(PetscFEGetNumComponents(fe, &num_comp));
134d68a66c4SJames Wright   PetscCall(PetscFEGetDualSpace(fe, &dual_space));
135d68a66c4SJames Wright   PetscCall(PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors));
136d68a66c4SJames Wright   P = num_dual_basis_vectors / num_comp;
137d68a66c4SJames Wright 
138d68a66c4SJames Wright   // Use depth label if no domain label present
139d68a66c4SJames Wright   if (!domain_label) {
140d68a66c4SJames Wright     PetscInt depth;
141d68a66c4SJames Wright 
142d68a66c4SJames Wright     PetscCall(DMPlexGetDepth(dm, &depth));
143d68a66c4SJames Wright     PetscCall(DMPlexGetDepthLabel(dm, &depth_label));
144d68a66c4SJames Wright     ids[0] = depth - height;
145d68a66c4SJames Wright   }
146d68a66c4SJames Wright 
147d68a66c4SJames Wright   // Get cell interp, grad, and quadrature data
148d68a66c4SJames Wright   PetscCall(DMGetFirstLabeledPoint(dm, dm, domain_label ? domain_label : depth_label, 1, ids, height, &first_point, NULL));
149d68a66c4SJames Wright   PetscCall(DMPlexGetCellType(dm, first_point, &cell_type));
150d68a66c4SJames Wright   elem_topo = ElemTopologyP2C(cell_type);
151d68a66c4SJames Wright   if (!elem_topo) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "DMPlex topology not supported");
152d68a66c4SJames Wright   {
153d68a66c4SJames Wright     size_t             q_points_size;
154d68a66c4SJames Wright     const PetscScalar *q_points_petsc;
155d68a66c4SJames Wright     PetscInt           q_dim;
156d68a66c4SJames Wright 
157d68a66c4SJames Wright     PetscCall(PetscQuadratureGetData(quadrature, &q_dim, NULL, &Q, &q_points_petsc, &q_weights));
158d68a66c4SJames Wright     q_points_size = Q * dim * sizeof(CeedScalar);
159d68a66c4SJames Wright     PetscCall(PetscCalloc(q_points_size, &q_points));
160d68a66c4SJames Wright     for (PetscInt q = 0; q < Q; q++) {
161d68a66c4SJames Wright       for (PetscInt d = 0; d < q_dim; d++) q_points[q * dim + d] = q_points_petsc[q * q_dim + d];
162d68a66c4SJames Wright     }
163d68a66c4SJames Wright   }
164d68a66c4SJames Wright 
165d68a66c4SJames Wright   {  // Convert to libCEED orientation
166d68a66c4SJames Wright     PetscBool       is_simplex  = PETSC_FALSE;
167d68a66c4SJames Wright     IS              permutation = NULL;
168d68a66c4SJames Wright     const PetscInt *permutation_indices;
169d68a66c4SJames Wright 
170d68a66c4SJames Wright     PetscCall(DMPlexIsSimplex(dm, &is_simplex));
171d68a66c4SJames Wright     if (!is_simplex) {
172d68a66c4SJames Wright       PetscSection section;
173d68a66c4SJames Wright 
174d68a66c4SJames Wright       // -- Get permutation
175d68a66c4SJames Wright       PetscCall(DMGetLocalSection(dm, &section));
176d68a66c4SJames Wright       PetscCall(PetscSectionGetClosurePermutation(section, (PetscObject)dm, dim, num_comp * P, &permutation));
177d68a66c4SJames Wright       PetscCall(ISGetIndices(permutation, &permutation_indices));
178d68a66c4SJames Wright     }
179d68a66c4SJames Wright 
180d68a66c4SJames Wright     // -- Copy interp, grad matrices
181d68a66c4SJames Wright     PetscCall(PetscCalloc(P * Q * sizeof(CeedScalar), &interp));
182d68a66c4SJames Wright     PetscCall(PetscCalloc(P * Q * dim * sizeof(CeedScalar), &grad));
183d68a66c4SJames Wright     const CeedInt c = 0;
184d68a66c4SJames Wright     for (CeedInt q = 0; q < Q; q++) {
185d68a66c4SJames Wright       for (CeedInt p_ceed = 0; p_ceed < P; p_ceed++) {
186d68a66c4SJames Wright         CeedInt p_petsc = is_simplex ? (p_ceed * num_comp) : permutation_indices[p_ceed * num_comp];
187d68a66c4SJames Wright 
188d68a66c4SJames Wright         interp[q * P + p_ceed] = basis_tabulation->T[0][((face * Q + q) * P * num_comp + p_petsc) * num_comp + c];
189d68a66c4SJames Wright         for (CeedInt d = 0; d < dim; d++) {
190d68a66c4SJames 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];
191d68a66c4SJames Wright         }
192d68a66c4SJames Wright       }
193d68a66c4SJames Wright     }
194d68a66c4SJames Wright 
195d68a66c4SJames Wright     // -- Cleanup
196d68a66c4SJames Wright     if (permutation) PetscCall(ISRestoreIndices(permutation, &permutation_indices));
197d68a66c4SJames Wright     PetscCall(ISDestroy(&permutation));
198d68a66c4SJames Wright   }
199d68a66c4SJames Wright 
200d68a66c4SJames Wright   // Finally, create libCEED basis
201d68a66c4SJames Wright   PetscCallCeed(ceed, CeedBasisCreateH1(ceed, elem_topo, num_comp, P, Q, interp, grad, q_points, q_weights, basis));
202d68a66c4SJames Wright   PetscCall(PetscFree(q_points));
203d68a66c4SJames Wright   PetscCall(PetscFree(interp));
204d68a66c4SJames Wright   PetscCall(PetscFree(grad));
205d68a66c4SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
206d68a66c4SJames Wright }
207d68a66c4SJames Wright 
208d68a66c4SJames Wright // -----------------------------------------------------------------------------
209d68a66c4SJames Wright // Get CEED Basis from DMPlex
210d68a66c4SJames Wright // -----------------------------------------------------------------------------
211d68a66c4SJames Wright PetscErrorCode CreateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label, CeedInt label_value, CeedInt height, CeedInt dm_field, CeedBasis *basis) {
212d68a66c4SJames Wright   PetscDS         ds;
213d68a66c4SJames Wright   PetscFE         fe;
214d68a66c4SJames Wright   PetscQuadrature quadrature;
215d68a66c4SJames Wright   PetscBool       is_simplex = PETSC_TRUE;
216d68a66c4SJames Wright   PetscInt        ds_field   = -1;
217d68a66c4SJames Wright 
218d68a66c4SJames Wright   PetscFunctionBeginUser;
219d68a66c4SJames Wright   // Get element information
220d68a66c4SJames Wright   PetscCall(DMGetRegionDS(dm, domain_label, NULL, &ds, NULL));
221d68a66c4SJames Wright   PetscCall(DMFieldToDSField(dm, domain_label, dm_field, &ds_field));
222d68a66c4SJames Wright   PetscCall(PetscDSGetDiscretization(ds, ds_field, (PetscObject *)&fe));
223d68a66c4SJames Wright   PetscCall(PetscFEGetHeightSubspace(fe, height, &fe));
224d68a66c4SJames Wright   PetscCall(PetscFEGetQuadrature(fe, &quadrature));
225d68a66c4SJames Wright 
226d68a66c4SJames Wright   // Check if simplex or tensor-product mesh
227d68a66c4SJames Wright   PetscCall(DMPlexIsSimplex(dm, &is_simplex));
228d68a66c4SJames Wright 
229d68a66c4SJames Wright   // Build libCEED basis
230d68a66c4SJames Wright   if (is_simplex) {
231d68a66c4SJames Wright     PetscTabulation basis_tabulation;
232d68a66c4SJames Wright     PetscInt        num_derivatives = 1, face = 0;
233d68a66c4SJames Wright 
234d68a66c4SJames Wright     PetscCall(PetscFEGetCellTabulation(fe, num_derivatives, &basis_tabulation));
235d68a66c4SJames Wright     PetscCall(BasisCreateFromTabulation(ceed, dm, domain_label, label_value, height, face, fe, basis_tabulation, quadrature, basis));
236d68a66c4SJames Wright   } else {
237d68a66c4SJames Wright     PetscDualSpace dual_space;
238d68a66c4SJames Wright     PetscInt       num_dual_basis_vectors;
239d68a66c4SJames Wright     PetscInt       dim, num_comp, P, Q;
240d68a66c4SJames Wright 
241d68a66c4SJames Wright     PetscCall(PetscFEGetSpatialDimension(fe, &dim));
242d68a66c4SJames Wright     PetscCall(PetscFEGetNumComponents(fe, &num_comp));
243d68a66c4SJames Wright     PetscCall(PetscFEGetDualSpace(fe, &dual_space));
244d68a66c4SJames Wright     PetscCall(PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors));
245d68a66c4SJames Wright     P = num_dual_basis_vectors / num_comp;
246d68a66c4SJames Wright     PetscCall(PetscQuadratureGetData(quadrature, NULL, NULL, &Q, NULL, NULL));
247d68a66c4SJames Wright 
248d68a66c4SJames Wright     CeedInt P_1d = (CeedInt)round(pow(P, 1.0 / dim));
249d68a66c4SJames Wright     CeedInt Q_1d = (CeedInt)round(pow(Q, 1.0 / dim));
250d68a66c4SJames Wright 
251d68a66c4SJames Wright     PetscCallCeed(ceed, CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp, P_1d, Q_1d, CEED_GAUSS, basis));
252d68a66c4SJames Wright   }
253d68a66c4SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
254d68a66c4SJames Wright }
255d68a66c4SJames Wright 
256d68a66c4SJames Wright /**
257d68a66c4SJames Wright   @brief Setup `DM` with FE space of appropriate degree
258d68a66c4SJames Wright 
259d68a66c4SJames Wright   Must be followed by `DMSetupByOrderEnd_FEM`
260d68a66c4SJames Wright 
261d68a66c4SJames Wright   @param[in]   setup_faces     Flag to setup face geometry
262d68a66c4SJames Wright   @param[in]   setup_coords    Flag to setup coordinate spaces
263d68a66c4SJames Wright   @param[in]   degree          Polynomial orders of field
264d68a66c4SJames Wright   @param[in]   coord_order     Polynomial order of coordinate basis, or `RATEL_DECIDE` for default
265d68a66c4SJames Wright   @param[in]   q_extra         Additional quadrature order, or `RATEL_DECIDE` for default
266d68a66c4SJames Wright   @param[in]   num_fields      Number of fields in solution vector
267d68a66c4SJames Wright   @param[in]   field_sizes     Array of number of components for each field
268d68a66c4SJames Wright   @param[out]  dm              `DM` to setup
269d68a66c4SJames Wright 
270d68a66c4SJames Wright   @return An error code: 0 - success, otherwise - failure
271d68a66c4SJames Wright **/
272d68a66c4SJames Wright PetscErrorCode DMSetupByOrderBegin_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra,
273d68a66c4SJames Wright                                        CeedInt num_fields, const CeedInt *field_sizes, DM dm) {
274d68a66c4SJames Wright   PetscInt  dim, q_order = degree + q_extra;
275d68a66c4SJames Wright   PetscBool is_simplex = PETSC_TRUE;
276d68a66c4SJames Wright   PetscFE   fe;
277d68a66c4SJames Wright   MPI_Comm  comm = PetscObjectComm((PetscObject)dm);
278d68a66c4SJames Wright 
279d68a66c4SJames Wright   PetscFunctionBeginUser;
280d68a66c4SJames Wright   PetscCall(DMPlexIsSimplex(dm, &is_simplex));
281d68a66c4SJames Wright 
282d68a66c4SJames Wright   // Setup DM
283d68a66c4SJames Wright   PetscCall(DMGetDimension(dm, &dim));
284d68a66c4SJames Wright   for (PetscInt i = 0; i < num_fields; i++) {
285d68a66c4SJames Wright     PetscFE  fe_face;
286d68a66c4SJames Wright     PetscInt q_order = degree + q_extra;
287d68a66c4SJames Wright 
288d68a66c4SJames Wright     PetscCall(PetscFECreateLagrange(comm, dim, field_sizes[i], is_simplex, degree, q_order, &fe));
289d68a66c4SJames Wright     if (setup_faces) PetscCall(PetscFEGetHeightSubspace(fe, 1, &fe_face));
290d68a66c4SJames Wright     PetscCall(DMAddField(dm, NULL, (PetscObject)fe));
291d68a66c4SJames Wright     PetscCall(PetscFEDestroy(&fe));
292d68a66c4SJames Wright   }
293d68a66c4SJames Wright   PetscCall(DMCreateDS(dm));
294d68a66c4SJames Wright 
295d68a66c4SJames Wright   // Project coordinates to enrich quadrature space
296d68a66c4SJames Wright   if (setup_coords) {
297d68a66c4SJames Wright     DM             dm_coord;
298d68a66c4SJames Wright     PetscDS        ds_coord;
299d68a66c4SJames Wright     PetscFE        fe_coord_current, fe_coord_new, fe_coord_face_new;
300d68a66c4SJames Wright     PetscDualSpace fe_coord_dual_space;
301d68a66c4SJames Wright     PetscInt       fe_coord_order, num_comp_coord;
302d68a66c4SJames Wright 
303d68a66c4SJames Wright     PetscCall(DMGetCoordinateDM(dm, &dm_coord));
304d68a66c4SJames Wright     PetscCall(DMGetCoordinateDim(dm, &num_comp_coord));
305d68a66c4SJames Wright     PetscCall(DMGetRegionDS(dm_coord, NULL, NULL, &ds_coord, NULL));
306d68a66c4SJames Wright     PetscCall(PetscDSGetDiscretization(ds_coord, 0, (PetscObject *)&fe_coord_current));
307d68a66c4SJames Wright     PetscCall(PetscFEGetDualSpace(fe_coord_current, &fe_coord_dual_space));
308d68a66c4SJames Wright     PetscCall(PetscDualSpaceGetOrder(fe_coord_dual_space, &fe_coord_order));
309d68a66c4SJames Wright 
310d68a66c4SJames Wright     // Create FE for coordinates
311d68a66c4SJames Wright     PetscCheck(fe_coord_order == 1 || coord_order == 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER_INPUT,
312d68a66c4SJames Wright                "Only linear mesh geometry supported. Recieved %d order", fe_coord_order);
313d68a66c4SJames Wright     PetscCall(PetscFECreateLagrange(comm, dim, num_comp_coord, is_simplex, fe_coord_order, q_order, &fe_coord_new));
314d68a66c4SJames Wright     if (setup_faces) PetscCall(PetscFEGetHeightSubspace(fe_coord_new, 1, &fe_coord_face_new));
315d68a66c4SJames Wright     PetscCall(DMProjectCoordinates(dm, fe_coord_new));
316d68a66c4SJames Wright     PetscCall(PetscFEDestroy(&fe_coord_new));
317d68a66c4SJames Wright   }
318d68a66c4SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
319d68a66c4SJames Wright }
320d68a66c4SJames Wright 
321d68a66c4SJames Wright /**
322d68a66c4SJames Wright   @brief Finish setting up `DM`
323d68a66c4SJames Wright 
324d68a66c4SJames Wright   Must be called after `DMSetupByOrderBegin_FEM` and all strong BCs have be declared via `DMAddBoundaries`.
325d68a66c4SJames Wright 
326d68a66c4SJames Wright   @param[in]   setup_coords    Flag to setup coordinate spaces
327d68a66c4SJames Wright   @param[out]  dm              `DM` to setup
328d68a66c4SJames Wright 
329d68a66c4SJames Wright   @return An error code: 0 - success, otherwise - failure
330d68a66c4SJames Wright **/
331d68a66c4SJames Wright PetscErrorCode DMSetupByOrderEnd_FEM(PetscBool setup_coords, DM dm) {
332d68a66c4SJames Wright   PetscBool is_simplex;
333d68a66c4SJames Wright 
334*f17d818dSJames Wright   PetscFunctionBeginUser;
335d68a66c4SJames Wright   PetscCall(DMPlexIsSimplex(dm, &is_simplex));
336d68a66c4SJames Wright   // Set tensor permutation if needed
337d68a66c4SJames Wright   if (!is_simplex) {
338d68a66c4SJames Wright     PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL));
339d68a66c4SJames Wright     if (setup_coords) {
340d68a66c4SJames Wright       DM dm_coord;
341d68a66c4SJames Wright 
342d68a66c4SJames Wright       PetscCall(DMGetCoordinateDM(dm, &dm_coord));
343d68a66c4SJames Wright       PetscCall(DMPlexSetClosurePermutationTensor(dm_coord, PETSC_DETERMINE, NULL));
344d68a66c4SJames Wright     }
345d68a66c4SJames Wright   }
346d68a66c4SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
347d68a66c4SJames Wright }
348d68a66c4SJames Wright 
349d68a66c4SJames Wright /**
350d68a66c4SJames Wright   @brief Setup `DM` with FE space of appropriate degree with no boundary conditions
351d68a66c4SJames Wright 
352d68a66c4SJames Wright   Calls `DMSetupByOrderBegin_FEM` and `DMSetupByOrderEnd_FEM` successively
353d68a66c4SJames Wright 
354d68a66c4SJames Wright   @param[in]   setup_faces     Flag to setup face geometry
355d68a66c4SJames Wright   @param[in]   setup_coords    Flag to setup coordinate spaces
356d68a66c4SJames Wright   @param[in]   degree          Polynomial orders of field
357d68a66c4SJames Wright   @param[in]   coord_order     Polynomial order of coordinate basis, or `RATEL_DECIDE` for default
358d68a66c4SJames Wright   @param[in]   q_extra         Additional quadrature order, or `RATEL_DECIDE` for default
359d68a66c4SJames Wright   @param[in]   num_fields      Number of fields in solution vector
360d68a66c4SJames Wright   @param[in]   field_sizes     Array of number of components for each field
361d68a66c4SJames Wright   @param[out]  dm              `DM` to setup
362d68a66c4SJames Wright 
363d68a66c4SJames Wright   @return An error code: 0 - success, otherwise - failure
364d68a66c4SJames Wright **/
365d68a66c4SJames Wright PetscErrorCode DMSetupByOrder_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra,
366d68a66c4SJames Wright                                   CeedInt num_fields, const CeedInt *field_sizes, DM dm) {
367d68a66c4SJames Wright   PetscFunctionBeginUser;
368d68a66c4SJames Wright   PetscCall(DMSetupByOrderBegin_FEM(setup_faces, setup_coords, degree, coord_order, q_extra, num_fields, field_sizes, dm));
369d68a66c4SJames Wright   PetscCall(DMSetupByOrderEnd_FEM(setup_coords, dm));
370d68a66c4SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
371d68a66c4SJames Wright }
372