Smith  0.1
Smith is an implicit thermal structural mechanics simulation code.
element_restriction.cpp
1 // Copyright (c) Lawrence Livermore National Security, LLC and
2 // other Smith Project Developers. See the top-level LICENSE file for
3 // details.
4 //
5 // SPDX-License-Identifier: (BSD-3-Clause)
6 
7 #include "smith/numerics/functional/element_restriction.hpp"
8 
9 #include <cstdlib>
10 #include <cstring>
11 #include <iostream>
12 
13 #include "smith/numerics/functional/geometry.hpp"
15 
16 std::vector<std::vector<int> > lexicographic_permutations(int p)
17 {
18  // p == 0 is admissible for L2 spaces, but lexicographic permutations
19  // aren't needed in that corner case
20  if (p == 0) {
21  return {};
22  }
23 
24  std::vector<std::vector<int> > output(mfem::Geometry::Type::NUM_GEOMETRIES);
25 
26  {
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];
31  }
32  output[mfem::Geometry::Type::SEGMENT] = native_to_lex;
33  }
34 
35  {
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];
40  }
41  output[mfem::Geometry::Type::TRIANGLE] = native_to_lex;
42  }
43 
44  {
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];
49  }
50  output[mfem::Geometry::Type::SQUARE] = native_to_lex;
51  }
52 
53  {
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];
58  }
59  output[mfem::Geometry::Type::TETRAHEDRON] = native_to_lex;
60  }
61 
62  {
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];
67  }
68  output[mfem::Geometry::Type::CUBE] = native_to_lex;
69  }
70 
71  // other geometries are not defined, as they are not currently used
72 
73  return output;
74 }
75 
76 Array2D<int> face_permutations(mfem::Geometry::Type geom, int p)
77 {
78  if (geom == mfem::Geometry::Type::SEGMENT) {
79  Array2D<int> output(2, p + 1);
80  for (int i = 0; i <= p; i++) {
81  output(0, i) = i;
82  output(1, i) = p - i;
83  }
84  return output;
85  }
86 
87  if (geom == mfem::Geometry::Type::TRIANGLE) {
88  // v = {{0, 0}, {1, 0}, {0, 1}};
89  // f = Transpose[{{0, 1, 2}, {1, 0, 2}, {2, 0, 1}, {2, 1, 0}, {1, 2, 0}, {0, 2, 1}} + 1];
90  // p v[[f[[1]]]] + (v[[f[[2]]]] - v[[f[[1]]]]) i + (v[[f[[3]]]] - v[[f[[1]]]]) j
91  //
92  // {{i, j}, {p-i-j, j}, {j, p-i-j}, {i, p-i-j}, {p-i-j, i}, {j, i}}
93  Array2D<int> output(6, (p + 1) * (p + 2) / 2);
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;
104  }
105  }
106  return output;
107  }
108 
109  if (geom == mfem::Geometry::Type::SQUARE) {
110  // v = {{0, 0}, {1, 0}, {1, 1}, {0, 1}};
111  // f = Transpose[{{0, 1, 2, 3}, {0, 3, 2, 1}, {1, 2, 3, 0}, {1, 0, 3, 2},
112  // {2, 3, 0, 1}, {2, 1, 0, 3}, {3, 0, 1, 2}, {3, 2, 1, 0}} + 1];
113  // p v[[f[[1]]]] + (v[[f[[2]]]] - v[[f[[1]]]]) i + (v[[f[[4]]]] - v[[f[[1]]]]) j
114  //
115  // {{i,j}, {j,i}, {p-j,i}, {p-i,j}, {p-i, p-j}, {p-j, p-i}, {j, p-i}, {i, p-j}}
116  Array2D<int> output(8, (p + 1) * (p + 1));
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;
129  }
130  }
131  return output;
132  }
133 
134  std::cout << "face_permutation(): unsupported geometry type" << std::endl;
135  exit(1);
136 }
137 
138 std::vector<Array2D<int> > geom_local_face_dofs(int p)
139 {
140  // FullSimplify[InterpolatingPolynomial[{
141  // {{0, 2}, (p + 1) + p},
142  // {{0, 1}, p + 1}, {{1, 1}, p + 2},
143  // {{0, 0}, 0}, {{1, 0}, 1}, {{2, 0}, 2}
144  // }, {x, y}]]
145  //
146  // x + 1/2 (3 + 2 p - y) y
147  auto tri_id = [p](int x, int y) { return x + ((3 + 2 * p - y) * y) / 2; };
148 
149  // FullSimplify[InterpolatingPolynomial[{
150  // {{0, 3}, ((p - 1) (p) + (p) (p + 1) + (p + 1) (p + 2))/2},
151  // {{0, 2}, ((p) (p + 1) + (p + 1) (p + 2))/2}, {{1, 2}, p - 1 + ((p) (p + 1) + (p + 1) (p + 2))/2},
152  // {{0, 1}, (p + 1) (p + 2)/2}, {{1, 1}, p + (p + 1) (p + 2)/2}, {{2, 1}, 2 p - 1 + (p + 1) (p + 2)/2},
153  // {{0, 0}, 0}, {{1, 0}, p + 1}, {{2, 0}, 2 p + 1}, {{3, 0}, 3 p}
154  // }, {y, z}]] + x
155  //
156  // x + (z (11 + p (12 + 3 p) - 6 y + z (z - 6 - 3 p)) - 3 y (y - 2 p - 3))/6
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;
159  };
160 
161  auto quad_id = [p](int x, int y) { return ((p + 1) * y) + x; };
162 
163  auto hex_id = [p](int x, int y, int z) { return (p + 1) * ((p + 1) * z + y) + x; };
164 
165  std::vector<Array2D<int> > output(mfem::Geometry::Type::NUM_GEOMETRIES);
166 
167  Array2D<int> tris(3, p + 1);
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);
172  }
173  output[mfem::Geometry::Type::TRIANGLE] = tris;
174 
175  Array2D<int> quads(4, p + 1);
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);
181  }
182  output[mfem::Geometry::Type::SQUARE] = quads;
183 
184  // v = {{0, 0, 0}, {1, 0, 0}, {0, 1, 0}, {0, 0, 1}};
185  // f = Transpose[{{1, 2, 3}, {0, 3, 2}, {0, 1, 3}, {0, 2, 1}} + 1];
186  // p v[[f[[1]]]] + (v[[f[[2]]]] - v[[f[[1]]]]) j + (v[[f[[3]]]] - v[[f[[1]]]]) k
187  //
188  // {{p-j-k, j, k}, {0, k, j}, {j, 0, k}, {k, j, 0}}
189  Array2D<int> tets(4, (p + 1) * (p + 2) / 2);
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);
197  }
198  }
199  output[mfem::Geometry::Type::TETRAHEDRON] = tets;
200 
201  // v = {{0, 0, 0}, {1, 0, 0}, {1, 1, 0}, {0, 1, 0},
202  // {0, 0, 1}, {1, 0, 1}, {1, 1, 1}, {0, 1, 1}};
203  // f = Transpose[{{3, 2, 1, 0}, {0, 1, 5, 4}, {1, 2, 6, 5},
204  // {2, 3, 7, 6}, {3, 0, 4, 7}, {4, 5, 6, 7}} + 1];
205  // p v[[f[[1]]]] + (v[[f[[2]]]] - v[[f[[1]]]]) j + (v[[f[[4]]]] - v[[f[[1]]]]) k
206  //
207  // {{j, p-k, 0}, {j, 0, k}, {p, j, k}, {p-j, p, k}, {0, p-j, k}, {j, k, p}}
208  Array2D<int> hexes(6, (p + 1) * (p + 1));
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);
218  }
219  }
220  output[mfem::Geometry::Type::CUBE] = hexes;
221 
222  return output;
223 }
224 
225 axom::Array<DoF, 2, smith::detail::host_memory_space> GetElementRestriction(const smith::fes_t* fes,
226  mfem::Geometry::Type geom)
227 {
228  std::vector<DoF> elem_dofs{};
229  mfem::Mesh* mesh = fes->GetMesh();
230 
231  // note: this assumes that all the elements are the same polynomial order
232  int p = fes->GetElementOrder(0);
233  std::vector<std::vector<int> > lex_perm = lexicographic_permutations(p);
234 
235  uint64_t n = 0;
236 
237  for (int elem = 0; elem < fes->GetNE(); elem++) {
238  // discard elements with the wrong geometry
239  if (mesh->GetElementGeometry(elem) != geom) continue;
240 
241  mfem::Array<int> dofs;
242 
243  [[maybe_unused]] auto* dof_transformation = fes->GetElementDofs(elem, dofs);
244 
245  // mfem returns the H1 dofs in "native" order, so we need
246  // to apply the native-to-lexicographic permutation
247  if (isH1(*fes)) {
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)]])});
250  }
251  }
252 
253  // the dofs mfem returns for Hcurl include information about
254  // dof orientation, but not for triangle faces on 3D elements.
255  // So, we need to manually
256  if (isHcurl(*fes)) {
257  // TODO
258  // TODO
259  // TODO
260  uint64_t sign = 1;
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});
264  }
265  }
266 
267  // mfem returns DG dofs in lexicographic order already
268  // so no permutation is required here
269  if (isDG(*fes)) {
270  for (int k = 0; k < dofs.Size(); k++) {
271  elem_dofs.push_back({uint64_t(dofs[k])});
272  }
273  }
274 
275  n++;
276  }
277 
278  if (n == 0) {
279  return axom::Array<DoF, 2, smith::detail::host_memory_space>{};
280  } else {
281  uint64_t dofs_per_elem = elem_dofs.size() / n;
282  axom::Array<DoF, 2, smith::detail::host_memory_space> output(n, dofs_per_elem);
283  std::memcpy(output.data(), elem_dofs.data(), sizeof(DoF) * n * dofs_per_elem);
284  return output;
285  }
286 }
287 
288 axom::Array<DoF, 2, smith::detail::host_memory_space> GetElementDofs(const smith::fes_t* fes, mfem::Geometry::Type geom,
289  const std::vector<int>& mfem_elem_ids)
290 
291 {
292  std::vector<DoF> elem_dofs{};
293  mfem::Mesh* mesh = fes->GetMesh();
294 
295  // note: this assumes that all the elements are the same polynomial order
296  int p = fes->GetElementOrder(0);
297  std::vector<std::vector<int> > lex_perm = lexicographic_permutations(p);
298 
299  uint64_t n = 0;
300 
301  for (auto elem : mfem_elem_ids) {
302  // discard elements with the wrong geometry
303  if (mesh->GetElementGeometry(elem) != geom) {
304  SLIC_ERROR("encountered incorrect element geometry type");
305  }
306 
307  mfem::Array<int> dofs;
308 
309  [[maybe_unused]] auto* dof_transformation = fes->GetElementDofs(elem, dofs);
310 
311  // mfem returns the H1 dofs in "native" order, so we need
312  // to apply the native-to-lexicographic permutation
313  if (isH1(*fes)) {
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)]])});
316  }
317  }
318 
319  // the dofs mfem returns for Hcurl include information about
320  // dof orientation, but not for triangle faces on 3D elements.
321  // So, we need to manually
322  if (isHcurl(*fes)) {
323  // TODO
324  // TODO
325  // TODO
326  uint64_t sign = 1;
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});
330  }
331  }
332 
333  // mfem returns DG dofs in lexicographic order already
334  // so no permutation is required here
335  if (isDG(*fes)) {
336  for (int k = 0; k < dofs.Size(); k++) {
337  elem_dofs.push_back({uint64_t(dofs[k])});
338  }
339  }
340 
341  n++;
342  }
343 
344  if (n == 0) {
345  return axom::Array<DoF, 2, smith::detail::host_memory_space>{};
346  } else {
347  uint64_t dofs_per_elem = elem_dofs.size() / n;
348  axom::Array<DoF, 2, smith::detail::host_memory_space> output(n, dofs_per_elem);
349  std::memcpy(output.data(), elem_dofs.data(), sizeof(DoF) * n * dofs_per_elem);
350  return output;
351  }
352 }
353 
354 axom::Array<DoF, 2, smith::detail::host_memory_space> GetFaceDofs(const smith::fes_t* fes,
355  mfem::Geometry::Type face_geom, FaceType type)
356 {
357  std::vector<DoF> face_dofs;
358  mfem::Mesh* mesh = fes->GetMesh();
359  mfem::Table* face_to_elem = mesh->GetFaceToElementTable();
360 
361  // note: this assumes that all the elements are the same polynomial order
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);
366 
367  uint64_t n = 0;
368 
369  for (int f = 0; f < fes->GetNF(); f++) {
370  auto faceinfo = mesh->GetFaceInformation(f);
371 
372  // discard faces with the wrong geometry or type
373  if (mesh->GetFaceGeometry(f) != face_geom) continue;
374  if (faceinfo.IsInterior() && type == FaceType::BOUNDARY) continue;
375  if (faceinfo.IsBoundary() && type == FaceType::INTERIOR) continue;
376 
377  // mfem doesn't provide this connectivity info for DG spaces directly,
378  // so we have to get at it indirectly in several steps:
379  if (isDG(*fes)) {
380  // 1. find the element(s) that this face belongs to
381  mfem::Array<int> elem_ids;
382  face_to_elem->GetRow(f, elem_ids);
383 
384  for (auto elem : elem_ids) {
385  // 2a. get the list of faces (and their orientations) that belong to that element ...
386  mfem::Array<int> elem_side_ids, orientations;
387  if (mesh->Dimension() == 2) {
388  mesh->GetElementEdges(elem, elem_side_ids, orientations);
389 
390  // mfem returns {-1, 1} for edge orientations,
391  // but {0, 1, ... , n} for face orientations.
392  // Here, we renumber the edge orientations to
393  // {0 (no permutation), 1 (reversed)} so the values can be
394  // consistently used as indices into a permutation table
395  for (auto& o : orientations) {
396  o = (o == -1) ? 1 : 0;
397  }
398 
399  } else {
400  mesh->GetElementFaces(elem, elem_side_ids, orientations);
401  }
402 
403  // 2b. ... and find `i` such that `elem_side_ids[i] == f`
404  int i;
405  for (i = 0; i < elem_side_ids.Size(); i++) {
406  if (elem_side_ids[i] == f) break;
407  }
408 
409  // 3. get the dofs for the entire element
410  mfem::Array<int> elem_dof_ids;
411  fes->GetElementDofs(elem, elem_dof_ids);
412 
413  mfem::Geometry::Type elem_geom = mesh->GetElementGeometry(elem);
414 
415  // mfem uses different conventions for boundary element orientations in 2D and 3D.
416  // In 2D, mfem's official edge orientations on the boundary will always be a mix of
417  // CW and CCW, so we have to discard mfem's orientation information in order
418  // to get a consistent winding.
419  //
420  // In 3D, mfem does use a consistently CCW winding for boundary faces (I think).
421  int orientation = (mesh->Dimension() == 2 && type == FaceType::BOUNDARY) ? 0 : orientations[i];
422 
423  // 4. extract only the dofs that correspond to side `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)]));
426  }
427 
428  // boundary faces only belong to 1 element, so we exit early
429  if (type == FaceType::BOUNDARY) break;
430  }
431 
432  // H1 and Hcurl spaces are more straight-forward, since
433  // we can use FiniteElementSpace::GetFaceDofs() directly
434  } else {
435  mfem::Array<int> dofs;
436 
437  fes->GetFaceDofs(f, dofs);
438 
439  if (isHcurl(*fes)) {
440  for (int k = 0; k < dofs.Size(); k++) {
441  if (dofs[k] >= 0) {
442  face_dofs.push_back(DoF{uint64_t(dofs[k]), 0});
443  } else {
444  face_dofs.push_back(DoF{uint64_t(-1 - dofs[k]), 1});
445  }
446  }
447  } else {
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)]]));
450  }
451  }
452  }
453 
454  n++;
455  }
456 
457  delete face_to_elem;
458 
459  if (n == 0) {
460  return axom::Array<DoF, 2, smith::detail::host_memory_space>{};
461  } else {
462  uint64_t dofs_per_face = face_dofs.size() / n;
463  axom::Array<DoF, 2, smith::detail::host_memory_space> output(n, dofs_per_face);
464  std::memcpy(output.data(), face_dofs.data(), sizeof(DoF) * n * dofs_per_face);
465  return output;
466  }
467 }
468 
469 axom::Array<DoF, 2, smith::detail::host_memory_space> GetFaceDofs(const smith::fes_t* fes,
470  mfem::Geometry::Type face_geom,
471  const std::vector<int>& mfem_face_ids)
472 {
473  std::vector<DoF> face_dofs;
474  mfem::Mesh* mesh = fes->GetMesh();
475  mfem::Table* face_to_elem = mesh->GetFaceToElementTable();
476 
477  // note: this assumes that all the elements are the same polynomial order
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);
482 
483  uint64_t n = 0;
484 
485  for (int f : mfem_face_ids) {
486  mfem::Mesh::FaceInformation info = mesh->GetFaceInformation(f);
487 
488  if (mesh->GetFaceGeometry(f) != face_geom) {
489  SLIC_ERROR("encountered incorrect face geometry type");
490  }
491 
492  // mfem doesn't provide this connectivity info for DG spaces directly (?),
493  // so we have to get at it indirectly in several steps:
494  if (isDG(*fes)) {
495  // 1. find the element(s) that this face belongs to
496  mfem::Array<int> elem_ids;
497  face_to_elem->GetRow(f, elem_ids);
498 
499  for (auto elem : elem_ids) {
500  // 2a. get the list of faces (and their orientations) that belong to that element ...
501  mfem::Array<int> elem_side_ids, orientations;
502  if (mesh->Dimension() == 2) {
503  mesh->GetElementEdges(elem, elem_side_ids, orientations);
504 
505  // mfem returns {-1, 1} for edge orientations,
506  // but {0, 1, ... , n} for face orientations.
507  // Here, we renumber the edge orientations to
508  // {0 (no permutation), 1 (reversed)} so the values can be
509  // consistently used as indices into a permutation table
510  for (auto& o : orientations) {
511  o = (o == -1) ? 1 : 0;
512  }
513 
514  } else {
515  mesh->GetElementFaces(elem, elem_side_ids, orientations);
516  }
517 
518  // 2b. ... and find `i` such that `elem_side_ids[i] == f`
519  // This generates a local_face_index that coincides with ones from GetFaceInformation
520  int i;
521  for (i = 0; i < elem_side_ids.Size(); i++) {
522  if (elem_side_ids[i] == f) break;
523  }
524 
525  // 3. get the dofs for the entire element
526  mfem::Array<int> elem_dof_ids;
527  fes->GetElementDofs(elem, elem_dof_ids);
528 
529  mfem::Geometry::Type elem_geom = mesh->GetElementGeometry(elem);
530 
531  // mfem uses different conventions for boundary element orientations in 2D and 3D.
532  // In 2D, mfem's official edge orientations on the boundary will always be a mix of
533  // CW and CCW, so we have to discard mfem's orientation information in order
534  // to get a consistent winding.
535  //
536  // In 3D, mfem does use a consistently CCW winding for boundary faces (I think).
537  int orientation = (mesh->Dimension() == 2 && info.IsBoundary()) ? 0 : orientations[i];
538 
539  // 4. extract only the dofs that correspond to side `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)]));
542  }
543 
544  // 5. add remaining dofs that were omitted on shared faces
545  // FaceNbrData includes all dofs in the "volumetric" element adjacent to the shared faces
546  // We need to again find our face on the ghost element and extract only the face dofs
547  // This may not yet work for NC faces
548  if (info.IsShared()) {
549  // FaceNbrData should be already exchanged in Functional when invoked by ParGridFunction
550  // which also invokes ExchangeFaceNbrData on subsequent ParFiniteElementSpace and ParMesh
551  int components_per_node = fes->GetVDim();
552  bool by_vdim = fes->GetOrdering() == mfem::Ordering::byVDIM;
553  int LSize = fes->GetProlongationMatrix()->Height();
554 
555  // Additional safeguard to make sure ByNODES ordering in vector L2 field raise an error
556  SLIC_ERROR_IF(!by_vdim && components_per_node > 1, "Unsupported: L2 vector field ordered by nodes");
557 
558  mfem::Array<int> shared_elem_vdof_ids;
559  mfem::Array<int> shared_elem_dof_ids;
560 
561  // on my processor
562  // VVVVVVV
563  // dofs for the face [0 1 3 4 | 5 8 12 13]
564  // ^^^^^^^^^
565  // on the other processor
566  int other_element_id = info.element[1].index;
567  fes->GetFaceNbrElementVDofs(other_element_id, shared_elem_vdof_ids); // indices into vector from FaceNbrData
568 
569  // byVDIM == x y z x y z x y z x y z
570  // byNODES == x x x x y y y y z z z z
571  int dofs_per_element = shared_elem_vdof_ids.Size() / components_per_node;
572  int stride = (by_vdim) ? components_per_node : 1;
573 
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]); // Change vdof to dof indices
577  }
578 
579  // Find the dofs on the shared face on the ghost element side
580  // Neighbor element local_face_index and orientation is already available from GetFaceInformation
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) {
585  // In 2D, orientations[i] sometimes is inverted as compared to the ones from GetFaceInformation
586  // In this case, we stay with the orientations constructed previously by GetElementEdges, which is strictly
587  // CCW So we set ghost orientation to be exactly the opposite of the original orientation (orientations[i])
588  // side note: even if you use info.element[1].orientation, looks like it's still fine.
589  // The only difference is the order those face dof indices gets registered in std::vector face_dofs
590  ghost_orientation = (orientation == 0) ? 1 : 0;
591  } else {
592  // In 3D, I think it's consistently CCW. So we directly use the orientation from GetFaceInformation
593  // This orientation is already consistent with the permutation in {0, 1, ... , n} form
594  ghost_orientation = info.element[1].orientation;
595  }
596 
597  // extract only the dofs that correspond to side `j`
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));
601  }
602  }
603  }
604 
605  // H1 and Hcurl spaces are more straight-forward, since
606  // we can use FiniteElementSpace::GetFaceDofs() directly
607  } else {
608  mfem::Array<int> dofs;
609 
610  fes->GetFaceDofs(f, dofs);
611 
612  if (isHcurl(*fes)) {
613  for (int k = 0; k < dofs.Size(); k++) {
614  if (dofs[k] >= 0) {
615  face_dofs.push_back(DoF{uint64_t(dofs[k]), 0});
616  } else {
617  face_dofs.push_back(DoF{uint64_t(-1 - dofs[k]), 1});
618  }
619  }
620  } else {
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)]]));
623  }
624  }
625  }
626 
627  n++;
628  }
629 
630  delete face_to_elem;
631 
632  if (n == 0) {
633  return axom::Array<DoF, 2, smith::detail::host_memory_space>{};
634  } else {
635  uint64_t dofs_per_face = face_dofs.size() / n;
636  axom::Array<DoF, 2, smith::detail::host_memory_space> output(n, dofs_per_face);
637  std::memcpy(output.data(), face_dofs.data(), sizeof(DoF) * n * dofs_per_face);
638  return output;
639  }
640 }
641 
642 namespace smith {
643 
644 ElementRestriction::ElementRestriction(const fes_t* fes, mfem::Geometry::Type elem_geom,
645  const std::vector<int>& elem_ids)
646 {
647  int sdim = fes->GetMesh()->Dimension();
648  int gdim = dimension_of(elem_geom);
649 
650  if (gdim == sdim) {
651  dof_info = GetElementDofs(fes, elem_geom, elem_ids);
652  }
653  if (gdim + 1 == sdim) {
654  dof_info = GetFaceDofs(fes, elem_geom, elem_ids);
655  }
656 
657  ordering = fes->GetOrdering();
658 
659  lsize = uint64_t(fes->GetVSize());
660  components = uint64_t(fes->GetVDim());
662  num_elements = uint64_t(dof_info.shape()[0]);
663  nodes_per_elem = uint64_t(dof_info.shape()[1]);
665 }
666 
667 uint64_t ElementRestriction::ESize() const { return esize; }
668 
669 uint64_t ElementRestriction::LSize() const { return lsize; }
670 
671 DoF ElementRestriction::GetVDof(DoF node, uint64_t component) const
672 {
673  if (ordering == mfem::Ordering::Type::byNODES) {
674  return DoF{component * num_nodes + node.index(), (node.sign() == 1) ? 0ull : 1ull, node.orientation()};
675  } else {
676  return DoF{node.index() * components + component, (node.sign() == 1) ? 0ull : 1ull, node.orientation()};
677  }
678 }
679 
680 void ElementRestriction::GetElementVDofs(int i, std::vector<DoF>& vdofs) const
681 {
682  for (uint64_t c = 0; c < components; c++) {
683  for (uint64_t j = 0; j < nodes_per_elem; j++) {
684  vdofs[c * nodes_per_elem + j] = GetVDof(dof_info(i, j), c);
685  }
686  }
687 }
688 
689 void ElementRestriction::Gather(const mfem::Vector& L_vector, mfem::Vector& E_vector) const
690 {
691  for (uint64_t i = 0; i < num_elements; i++) {
692  for (uint64_t c = 0; c < components; c++) {
693  for (uint64_t j = 0; j < nodes_per_elem; j++) {
694  uint64_t E_id = (i * components + c) * nodes_per_elem + j;
695  uint64_t L_id = GetVDof(dof_info(i, j), c).index();
696  E_vector[int(E_id)] = L_vector[int(L_id)];
697  }
698  }
699  }
700 }
701 
702 void ElementRestriction::ScatterAdd(const mfem::Vector& E_vector, mfem::Vector& L_vector) const
703 {
704  for (uint64_t i = 0; i < num_elements; i++) {
705  for (uint64_t c = 0; c < components; c++) {
706  for (uint64_t j = 0; j < nodes_per_elem; j++) {
707  uint64_t E_id = (i * components + c) * nodes_per_elem + j;
708  uint64_t L_id = GetVDof(dof_info(i, j), c).index();
709  L_vector[int(L_id)] += E_vector[int(E_id)];
710  }
711  }
712  }
713 }
714 
716 
719 {
720  // TODO: changing the mfem_XXX_ids arrays to mfem_ids[XXX] would simplify this
721  if (domain.mfem_edge_ids_.size() > 0) {
722  restrictions[mfem::Geometry::SEGMENT] = ElementRestriction(fes, mfem::Geometry::SEGMENT, domain.mfem_edge_ids_);
723  }
724 
725  if (domain.mfem_tri_ids_.size() > 0) {
726  restrictions[mfem::Geometry::TRIANGLE] = ElementRestriction(fes, mfem::Geometry::TRIANGLE, domain.mfem_tri_ids_);
727  }
728 
729  if (domain.mfem_quad_ids_.size() > 0) {
730  restrictions[mfem::Geometry::SQUARE] = ElementRestriction(fes, mfem::Geometry::SQUARE, domain.mfem_quad_ids_);
731  }
732 
733  if (domain.mfem_tet_ids_.size() > 0) {
734  restrictions[mfem::Geometry::TETRAHEDRON] =
735  ElementRestriction(fes, mfem::Geometry::TETRAHEDRON, domain.mfem_tet_ids_);
736  }
737 
738  if (domain.mfem_hex_ids_.size() > 0) {
739  restrictions[mfem::Geometry::CUBE] = ElementRestriction(fes, mfem::Geometry::CUBE, domain.mfem_hex_ids_);
740  }
741 }
742 
744 {
745  uint64_t total = 0;
746  for (auto& [geom, restriction] : restrictions) {
747  total += restriction.ESize();
748  }
749  return total;
750 }
751 
752 uint64_t BlockElementRestriction::LSize() const { return (*restrictions.begin()).second.LSize(); }
753 
754 mfem::Array<int> BlockElementRestriction::bOffsets() const
755 {
756  mfem::Array<int> offsets(mfem::Geometry::NUM_GEOMETRIES + 1);
757 
758  offsets[0] = 0;
759  for (int i = 0; i < mfem::Geometry::NUM_GEOMETRIES; i++) {
760  auto g = mfem::Geometry::Type(i);
761  if (restrictions.count(g) > 0) {
762  offsets[g + 1] = offsets[g] + int(restrictions.at(g).ESize());
763  } else {
764  offsets[g + 1] = offsets[g];
765  }
766  }
767  return offsets;
768 };
769 
770 void BlockElementRestriction::Gather(const mfem::Vector& L_vector, mfem::BlockVector& E_block_vector) const
771 {
772  for (auto& [geom, restriction] : restrictions) {
773  restriction.Gather(L_vector, E_block_vector.GetBlock(geom));
774  }
775 }
776 
777 void BlockElementRestriction::ScatterAdd(const mfem::BlockVector& E_block_vector, mfem::Vector& L_vector) const
778 {
779  for (auto& [geom, restriction] : restrictions) {
780  restriction.ScatterAdd(E_block_vector.GetBlock(geom), L_vector);
781  }
782 }
783 
784 } // namespace smith
many of the functions in this file amount to extracting element indices from an mesh_t like
Accelerator functionality.
Definition: smith.cpp:36
constexpr SMITH_HOST_DEVICE auto type(const tuple< T... > &values)
a function intended to be used for extracting the ith type from a tuple.
Definition: tuple.hpp:376
constexpr int dimension_of(mfem::Geometry::Type g)
Returns the dimension of an element geometry.
Definition: geometry.hpp:55
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
uint64_t ESize() const
the size of the "E-vector" associated with this restriction operator
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...
BlockElementRestriction()
default ctor leaves this object uninitialized
uint64_t LSize() const
the size of the "L-vector" associated with this restriction operator
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...
std::map< mfem::Geometry::Type, ElementRestriction > restrictions
the individual ElementRestriction operators for each element geometry
mfem::Array< int > bOffsets() const
block offsets used when constructing mfem::HypreParVectors
a class for representing a geometric region that can be used for integration
Definition: domain.hpp:33
void GetElementVDofs(int i, std::vector< DoF > &dofs) const
Get a list of DoFs for element i
uint64_t esize
the size of the "E-vector"
uint64_t LSize() const
the size of the "L-vector" associated with this restriction operator
mfem::Ordering::Type ordering
whether the underlying dofs are arranged "byNodes" or "byVDim"
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 lsize
the size of the "L-vector"
ElementRestriction()
default ctor leaves this object uninitialized
uint64_t num_nodes
the total number of nodes in the mesh
axom::Array< DoF, 2, smith::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...
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 components
the number of components at each node
uint64_t num_elements
the number of elements of the given geometry
DoF GetVDof(DoF node, uint64_t component) const
get the dof information for a given node / component
uint64_t nodes_per_elem
the number of nodes in each element
uint64_t ESize() const
the size of the "E-vector" associated with this restriction operator