7 #include "serac/numerics/functional/element_restriction.hpp"
13 #include "serac/numerics/functional/geometry.hpp"
16 std::vector<std::vector<int> > lexicographic_permutations(
int p)
24 std::vector<std::vector<int> > output(mfem::Geometry::Type::NUM_GEOMETRIES);
27 auto P = mfem::H1_SegmentElement(p).GetLexicographicOrdering();
28 std::vector<int> native_to_lex(uint32_t(P.Size()));
29 for (
int i = 0; i < P.Size(); i++) {
30 native_to_lex[uint32_t(i)] = P[i];
32 output[mfem::Geometry::Type::SEGMENT] = native_to_lex;
36 auto P = mfem::H1_TriangleElement(p).GetLexicographicOrdering();
37 std::vector<int> native_to_lex(uint32_t(P.Size()));
38 for (
int i = 0; i < P.Size(); i++) {
39 native_to_lex[uint32_t(i)] = P[i];
41 output[mfem::Geometry::Type::TRIANGLE] = native_to_lex;
45 auto P = mfem::H1_QuadrilateralElement(p).GetLexicographicOrdering();
46 std::vector<int> native_to_lex(uint32_t(P.Size()));
47 for (
int i = 0; i < P.Size(); i++) {
48 native_to_lex[uint32_t(i)] = P[i];
50 output[mfem::Geometry::Type::SQUARE] = native_to_lex;
54 auto P = mfem::H1_TetrahedronElement(p).GetLexicographicOrdering();
55 std::vector<int> native_to_lex(uint32_t(P.Size()));
56 for (
int i = 0; i < P.Size(); i++) {
57 native_to_lex[uint32_t(i)] = P[i];
59 output[mfem::Geometry::Type::TETRAHEDRON] = native_to_lex;
63 auto P = mfem::H1_HexahedronElement(p).GetLexicographicOrdering();
64 std::vector<int> native_to_lex(uint32_t(P.Size()));
65 for (
int i = 0; i < P.Size(); i++) {
66 native_to_lex[uint32_t(i)] = P[i];
68 output[mfem::Geometry::Type::CUBE] = native_to_lex;
76 Array2D<int> face_permutations(mfem::Geometry::Type geom,
int p)
78 if (geom == mfem::Geometry::Type::SEGMENT) {
80 for (
int i = 0; i <= p; i++) {
87 if (geom == mfem::Geometry::Type::TRIANGLE) {
94 auto tri_id = [p](
int x,
int y) {
return x + ((3 + 2 * p - y) * y) / 2; };
95 for (
int j = 0; j <= p; j++) {
96 for (
int i = 0; i <= p - j; i++) {
97 int id = tri_id(i, j);
98 output(0, tri_id(i, j)) = id;
99 output(1, tri_id(p - i - j, j)) = id;
100 output(2, tri_id(j, p - i - j)) = id;
101 output(3, tri_id(i, p - i - j)) = id;
102 output(4, tri_id(p - i - j, i)) = id;
103 output(5, tri_id(j, i)) = id;
109 if (geom == mfem::Geometry::Type::SQUARE) {
117 auto quad_id = [p](
int x,
int y) {
return ((p + 1) * y) + x; };
118 for (
int j = 0; j <= p; j++) {
119 for (
int i = 0; i <= p; i++) {
120 int id = quad_id(i, j);
121 output(0, quad_id(i, j)) = id;
122 output(1, quad_id(j, i)) = id;
123 output(2, quad_id(p - j, i)) = id;
124 output(3, quad_id(p - i, j)) = id;
125 output(4, quad_id(p - i, p - j)) = id;
126 output(5, quad_id(p - j, p - i)) = id;
127 output(6, quad_id(j, p - i)) = id;
128 output(7, quad_id(i, p - j)) = id;
134 std::cout <<
"face_permutation(): unsupported geometry type" << std::endl;
138 std::vector<Array2D<int> > geom_local_face_dofs(
int p)
147 auto tri_id = [p](
int x,
int y) {
return x + ((3 + 2 * p - y) * y) / 2; };
157 auto tet_id = [p](
int x,
int y,
int z) {
158 return x + (z * (11 + p * (12 + 3 * p) - 6 * y + z * (z - 6 - 3 * p)) - 3 * y * (y - 2 * p - 3)) / 6;
161 auto quad_id = [p](
int x,
int y) {
return ((p + 1) * y) + x; };
163 auto hex_id = [p](
int x,
int y,
int z) {
return (p + 1) * ((p + 1) * z + y) + x; };
165 std::vector<Array2D<int> > output(mfem::Geometry::Type::NUM_GEOMETRIES);
168 for (
int k = 0; k <= p; k++) {
169 tris(0, k) = tri_id(k, 0);
170 tris(1, k) = tri_id(p - k, k);
171 tris(2, k) = tri_id(0, p - k);
173 output[mfem::Geometry::Type::TRIANGLE] = tris;
176 for (
int k = 0; k <= p; k++) {
177 quads(0, k) = quad_id(k, 0);
178 quads(1, k) = quad_id(p, k);
179 quads(2, k) = quad_id(p - k, p);
180 quads(3, k) = quad_id(0, p - k);
182 output[mfem::Geometry::Type::SQUARE] = quads;
190 for (
int k = 0; k <= p; k++) {
191 for (
int j = 0; j <= p - k; j++) {
192 int id = tri_id(j, k);
193 tets(0,
id) = tet_id(p - j - k, j, k);
194 tets(1,
id) = tet_id(0, k, j);
195 tets(2,
id) = tet_id(j, 0, k);
196 tets(3,
id) = tet_id(k, j, 0);
199 output[mfem::Geometry::Type::TETRAHEDRON] = tets;
209 for (
int k = 0; k <= p; k++) {
210 for (
int j = 0; j <= p; j++) {
211 int id = quad_id(j, k);
212 hexes(0,
id) = hex_id(j, p - k, 0);
213 hexes(1,
id) = hex_id(j, 0, k);
214 hexes(2,
id) = hex_id(p, j, k);
215 hexes(3,
id) = hex_id(p - j, p, k);
216 hexes(4,
id) = hex_id(0, p - j, k);
217 hexes(5,
id) = hex_id(j, k, p);
220 output[mfem::Geometry::Type::CUBE] = hexes;
225 axom::Array<DoF, 2, serac::detail::host_memory_space> GetElementRestriction(
const serac::fes_t* fes,
226 mfem::Geometry::Type geom)
228 std::vector<DoF> elem_dofs{};
229 mfem::Mesh* mesh = fes->GetMesh();
232 int p = fes->GetElementOrder(0);
233 std::vector<std::vector<int> > lex_perm = lexicographic_permutations(p);
237 for (
int elem = 0; elem < fes->GetNE(); elem++) {
239 if (mesh->GetElementGeometry(elem) != geom)
continue;
241 mfem::Array<int> dofs;
243 [[maybe_unused]]
auto* dof_transformation = fes->GetElementDofs(elem, dofs);
248 for (
int k = 0; k < dofs.Size(); k++) {
249 elem_dofs.push_back({uint64_t(dofs[lex_perm[uint32_t(geom)][uint32_t(k)]])});
261 uint64_t orientation = 0;
262 for (
int k = 0; k < dofs.Size(); k++) {
263 elem_dofs.push_back({uint64_t(dofs[k]), sign, orientation});
270 for (
int k = 0; k < dofs.Size(); k++) {
271 elem_dofs.push_back({uint64_t(dofs[k])});
279 return axom::Array<DoF, 2, serac::detail::host_memory_space>{};
281 uint64_t dofs_per_elem = elem_dofs.size() / n;
282 axom::Array<DoF, 2, serac::detail::host_memory_space> output(n, dofs_per_elem);
283 std::memcpy(output.data(), elem_dofs.data(),
sizeof(
DoF) * n * dofs_per_elem);
288 axom::Array<DoF, 2, serac::detail::host_memory_space> GetElementDofs(
const serac::fes_t* fes, mfem::Geometry::Type geom,
289 const std::vector<int>& mfem_elem_ids)
292 std::vector<DoF> elem_dofs{};
293 mfem::Mesh* mesh = fes->GetMesh();
296 int p = fes->GetElementOrder(0);
297 std::vector<std::vector<int> > lex_perm = lexicographic_permutations(p);
301 for (
auto elem : mfem_elem_ids) {
303 if (mesh->GetElementGeometry(elem) != geom) {
304 SLIC_ERROR(
"encountered incorrect element geometry type");
307 mfem::Array<int> dofs;
309 [[maybe_unused]]
auto* dof_transformation = fes->GetElementDofs(elem, dofs);
314 for (
int k = 0; k < dofs.Size(); k++) {
315 elem_dofs.push_back({uint64_t(dofs[lex_perm[uint32_t(geom)][uint32_t(k)]])});
327 uint64_t orientation = 0;
328 for (
int k = 0; k < dofs.Size(); k++) {
329 elem_dofs.push_back({uint64_t(dofs[k]), sign, orientation});
336 for (
int k = 0; k < dofs.Size(); k++) {
337 elem_dofs.push_back({uint64_t(dofs[k])});
345 return axom::Array<DoF, 2, serac::detail::host_memory_space>{};
347 uint64_t dofs_per_elem = elem_dofs.size() / n;
348 axom::Array<DoF, 2, serac::detail::host_memory_space> output(n, dofs_per_elem);
349 std::memcpy(output.data(), elem_dofs.data(),
sizeof(
DoF) * n * dofs_per_elem);
354 axom::Array<DoF, 2, serac::detail::host_memory_space> GetFaceDofs(
const serac::fes_t* fes,
355 mfem::Geometry::Type face_geom, FaceType
type)
357 std::vector<DoF> face_dofs;
358 mfem::Mesh* mesh = fes->GetMesh();
359 mfem::Table* face_to_elem = mesh->GetFaceToElementTable();
362 int p = fes->GetElementOrder(0);
363 Array2D<int> face_perm = face_permutations(face_geom, p);
364 std::vector<Array2D<int> > local_face_dofs = geom_local_face_dofs(p);
365 std::vector<std::vector<int> > lex_perm = lexicographic_permutations(p);
369 for (
int f = 0; f < fes->GetNF(); f++) {
370 auto faceinfo = mesh->GetFaceInformation(f);
373 if (mesh->GetFaceGeometry(f) != face_geom)
continue;
374 if (faceinfo.IsInterior() &&
type == FaceType::BOUNDARY)
continue;
375 if (faceinfo.IsBoundary() &&
type == FaceType::INTERIOR)
continue;
381 mfem::Array<int> elem_ids;
382 face_to_elem->GetRow(f, elem_ids);
384 for (
auto elem : elem_ids) {
386 mfem::Array<int> elem_side_ids, orientations;
387 if (mesh->Dimension() == 2) {
388 mesh->GetElementEdges(elem, elem_side_ids, orientations);
395 for (
auto& o : orientations) {
396 o = (o == -1) ? 1 : 0;
400 mesh->GetElementFaces(elem, elem_side_ids, orientations);
405 for (i = 0; i < elem_side_ids.Size(); i++) {
406 if (elem_side_ids[i] == f)
break;
410 mfem::Array<int> elem_dof_ids;
411 fes->GetElementDofs(elem, elem_dof_ids);
413 mfem::Geometry::Type elem_geom = mesh->GetElementGeometry(elem);
421 int orientation = (mesh->Dimension() == 2 &&
type == FaceType::BOUNDARY) ? 0 : orientations[i];
424 for (
auto k : face_perm(orientation)) {
425 face_dofs.push_back(uint64_t(elem_dof_ids[local_face_dofs[uint32_t(elem_geom)](i, k)]));
429 if (
type == FaceType::BOUNDARY)
break;
435 mfem::Array<int> dofs;
437 fes->GetFaceDofs(f, dofs);
440 for (
int k = 0; k < dofs.Size(); k++) {
442 face_dofs.push_back(
DoF{uint64_t(dofs[k]), 0});
444 face_dofs.push_back(
DoF{uint64_t(-1 - dofs[k]), 1});
448 for (
int k = 0; k < dofs.Size(); k++) {
449 face_dofs.push_back(uint64_t(dofs[lex_perm[uint32_t(face_geom)][uint32_t(k)]]));
460 return axom::Array<DoF, 2, serac::detail::host_memory_space>{};
462 uint64_t dofs_per_face = face_dofs.size() / n;
463 axom::Array<DoF, 2, serac::detail::host_memory_space> output(n, dofs_per_face);
464 std::memcpy(output.data(), face_dofs.data(),
sizeof(
DoF) * n * dofs_per_face);
469 axom::Array<DoF, 2, serac::detail::host_memory_space> GetFaceDofs(
const serac::fes_t* fes,
470 mfem::Geometry::Type face_geom,
471 const std::vector<int>& mfem_face_ids)
473 std::vector<DoF> face_dofs;
474 mfem::Mesh* mesh = fes->GetMesh();
475 mfem::Table* face_to_elem = mesh->GetFaceToElementTable();
478 int p = fes->GetElementOrder(0);
479 Array2D<int> face_perm = face_permutations(face_geom, p);
480 std::vector<Array2D<int> > local_face_dofs = geom_local_face_dofs(p);
481 std::vector<std::vector<int> > lex_perm = lexicographic_permutations(p);
485 for (
int f : mfem_face_ids) {
486 mfem::Mesh::FaceInformation info = mesh->GetFaceInformation(f);
488 if (mesh->GetFaceGeometry(f) != face_geom) {
489 SLIC_ERROR(
"encountered incorrect face geometry type");
496 mfem::Array<int> elem_ids;
497 face_to_elem->GetRow(f, elem_ids);
499 for (
auto elem : elem_ids) {
501 mfem::Array<int> elem_side_ids, orientations;
502 if (mesh->Dimension() == 2) {
503 mesh->GetElementEdges(elem, elem_side_ids, orientations);
510 for (
auto& o : orientations) {
511 o = (o == -1) ? 1 : 0;
515 mesh->GetElementFaces(elem, elem_side_ids, orientations);
521 for (i = 0; i < elem_side_ids.Size(); i++) {
522 if (elem_side_ids[i] == f)
break;
526 mfem::Array<int> elem_dof_ids;
527 fes->GetElementDofs(elem, elem_dof_ids);
529 mfem::Geometry::Type elem_geom = mesh->GetElementGeometry(elem);
537 int orientation = (mesh->Dimension() == 2 && info.IsBoundary()) ? 0 : orientations[i];
540 for (
auto k : face_perm(orientation)) {
541 face_dofs.push_back(uint64_t(elem_dof_ids[local_face_dofs[uint32_t(elem_geom)](i, k)]));
548 if (info.IsShared()) {
551 int components_per_node = fes->GetVDim();
552 bool by_vdim = fes->GetOrdering() == mfem::Ordering::byVDIM;
553 int LSize = fes->GetProlongationMatrix()->Height();
556 SLIC_ERROR_IF(!by_vdim && components_per_node > 1,
"Unsupported: L2 vector field ordered by nodes");
558 mfem::Array<int> shared_elem_vdof_ids;
559 mfem::Array<int> shared_elem_dof_ids;
566 int other_element_id = info.element[1].index;
567 fes->GetFaceNbrElementVDofs(other_element_id, shared_elem_vdof_ids);
571 int dofs_per_element = shared_elem_vdof_ids.Size() / components_per_node;
572 int stride = (by_vdim) ? components_per_node : 1;
574 shared_elem_dof_ids.SetSize(dofs_per_element);
575 for (
int k = 0; k < dofs_per_element; ++k) {
576 shared_elem_dof_ids[k] = fes->VDofToDof(shared_elem_vdof_ids[k * stride]);
581 mfem::Geometry::Type ghost_geom = fes->GetParMesh()->face_nbr_elements[other_element_id]->GetGeometryType();
582 int j = info.element[1].local_face_id;
583 int ghost_orientation;
584 if (mesh->Dimension() == 2) {
590 ghost_orientation = (orientation == 0) ? 1 : 0;
594 ghost_orientation = info.element[1].orientation;
598 for (
auto k : face_perm(ghost_orientation)) {
599 face_dofs.push_back(uint64_t(shared_elem_dof_ids[local_face_dofs[uint32_t(ghost_geom)](j, k)] +
600 LSize / components_per_node));
608 mfem::Array<int> dofs;
610 fes->GetFaceDofs(f, dofs);
613 for (
int k = 0; k < dofs.Size(); k++) {
615 face_dofs.push_back(
DoF{uint64_t(dofs[k]), 0});
617 face_dofs.push_back(
DoF{uint64_t(-1 - dofs[k]), 1});
621 for (
int k = 0; k < dofs.Size(); k++) {
622 face_dofs.push_back(uint64_t(dofs[lex_perm[uint32_t(face_geom)][uint32_t(k)]]));
633 return axom::Array<DoF, 2, serac::detail::host_memory_space>{};
635 uint64_t dofs_per_face = face_dofs.size() / n;
636 axom::Array<DoF, 2, serac::detail::host_memory_space> output(n, dofs_per_face);
637 std::memcpy(output.data(), face_dofs.data(),
sizeof(
DoF) * n * dofs_per_face);
645 const std::vector<int>& elem_ids)
647 int sdim = fes->GetMesh()->Dimension();
651 dof_info = GetElementDofs(fes, elem_geom, elem_ids);
653 if (gdim + 1 == sdim) {
654 dof_info = GetFaceDofs(fes, elem_geom, elem_ids);
659 lsize = uint64_t(fes->GetVSize());
673 if (
ordering == mfem::Ordering::Type::byNODES) {
696 E_vector[int(E_id)] = L_vector[int(L_id)];
709 L_vector[int(L_id)] += E_vector[int(E_id)];
721 if (domain.mfem_edge_ids_.size() > 0) {
725 if (domain.mfem_tri_ids_.size() > 0) {
729 if (domain.mfem_quad_ids_.size() > 0) {
733 if (domain.mfem_tet_ids_.size() > 0) {
738 if (domain.mfem_hex_ids_.size() > 0) {
747 total += restriction.ESize();
756 mfem::Array<int> offsets(mfem::Geometry::NUM_GEOMETRIES + 1);
759 for (
int i = 0; i < mfem::Geometry::NUM_GEOMETRIES; i++) {
760 auto g = mfem::Geometry::Type(i);
762 offsets[g + 1] = offsets[g] + int(
restrictions.at(g).ESize());
764 offsets[g + 1] = offsets[g];
773 restriction.Gather(L_vector, E_block_vector.GetBlock(geom));
780 restriction.ScatterAdd(E_block_vector.GetBlock(geom), L_vector);
many of the functions in this file amount to extracting element indices from an mesh_t like
Accelerator functionality.
constexpr SERAC_HOST_DEVICE auto type(const tuple< T... > &values)
a function intended to be used for extracting the ith type from a tuple.
constexpr int dimension_of(mfem::Geometry::Type g)
Returns the dimension of an element geometry.
a struct of metadata (index, sign, orientation) associated with a degree of freedom
int sign() const
get the sign field of this DoF
uint64_t index() const
get the index field of this DoF
uint64_t orientation() const
get the orientation field of this DoF
BlockElementRestriction()
default ctor leaves this object uninitialized
mfem::Array< int > bOffsets() const
block offsets used when constructing mfem::HypreParVectors
void ScatterAdd(const mfem::BlockVector &E_block_vector, mfem::Vector &L_vector) const
"E->L" in mfem parlance, each element scatter-adds its local vector into the appropriate place in the...
uint64_t ESize() const
the size of the "E-vector" associated with this restriction operator
std::map< mfem::Geometry::Type, ElementRestriction > restrictions
the individual ElementRestriction operators for each element geometry
void Gather(const mfem::Vector &L_vector, mfem::BlockVector &E_block_vector) const
"L->E" in mfem parlance, each element gathers the values that belong to it, and stores them in the "E...
uint64_t LSize() const
the size of the "L-vector" associated with this restriction operator
a class for representing a geometric region that can be used for integration
mfem::Ordering::Type ordering
whether the underlying dofs are arranged "byNodes" or "byVDim"
axom::Array< DoF, 2, serac::detail::host_memory_space > dof_info
a 2D array (num_elements-by-nodes_per_elem) holding the dof info extracted from the finite element sp...
uint64_t num_nodes
the total number of nodes in the mesh
ElementRestriction()
default ctor leaves this object uninitialized
uint64_t nodes_per_elem
the number of nodes in each element
uint64_t lsize
the size of the "L-vector"
void GetElementVDofs(int i, std::vector< DoF > &dofs) const
Get a list of DoFs for element i
void ScatterAdd(const mfem::Vector &E_vector, mfem::Vector &L_vector) const
"E->L" in mfem parlance, each element scatter-adds its local vector into the appropriate place in the...
uint64_t esize
the size of the "E-vector"
uint64_t num_elements
the number of elements of the given geometry
uint64_t components
the number of components at each node
void Gather(const mfem::Vector &L_vector, mfem::Vector &E_vector) const
"L->E" in mfem parlance, each element gathers the values that belong to it, and stores them in the "E...
uint64_t LSize() const
the size of the "L-vector" associated with this restriction operator
uint64_t ESize() const
the size of the "E-vector" associated with this restriction operator
DoF GetVDof(DoF node, uint64_t component) const
get the dof information for a given node / component