Smith  0.1
Smith is an implicit thermal structural mechanics simulation code.
domain.hpp
Go to the documentation of this file.
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 #pragma once
8 
9 #include <cstdlib>
10 #include <vector>
11 #include <array>
12 #include <cstdint>
13 #include <functional>
14 #include <map>
15 #include <set>
16 
17 #include "mfem.hpp"
18 
21 #include "smith/numerics/functional/element_restriction.hpp"
22 #include "smith/numerics/functional/typedefs.hpp"
23 
24 namespace smith {
25 
26 struct BlockElementRestriction;
27 
33 struct Domain {
35  enum Type
36  {
37  Elements,
38  BoundaryElements,
39  InteriorFaces
40  };
41 
42  static constexpr int num_types = 3;
43 
45  const mesh_t& mesh_;
46 
48  int dim_;
49 
52 
83 
85  std::vector<int> edge_ids_;
86  std::vector<int> tri_ids_;
87  std::vector<int> quad_ids_;
88  std::vector<int> tet_ids_;
89  std::vector<int> hex_ids_;
90 
91  std::vector<int> mfem_edge_ids_;
92  std::vector<int> mfem_tri_ids_;
93  std::vector<int> mfem_quad_ids_;
94  std::vector<int> mfem_tet_ids_;
95  std::vector<int> mfem_hex_ids_;
97 
103  std::map<FunctionSpace, BlockElementRestriction> restriction_operators;
104 
106  std::vector<int> shared_interior_face_ids_;
107 
111  Domain(const mesh_t& m, int d, Type type) : mesh_(m), dim_(d), type_(type) {}
112 
119  static Domain ofVertices(const mesh_t& mesh, std::function<bool(vec2)> func);
120 
122  static Domain ofVertices(const mesh_t& mesh, std::function<bool(vec3)> func);
123 
131  static Domain ofEdges(const mesh_t& mesh, std::function<bool(std::vector<vec2>, int)> func);
132 
134  static Domain ofEdges(const mesh_t& mesh, std::function<bool(std::vector<vec3>)> func);
135 
143  static Domain ofFaces(const mesh_t& mesh, std::function<bool(std::vector<vec2>, int)> func);
144 
146  static Domain ofFaces(const mesh_t& mesh, std::function<bool(std::vector<vec3>, int)> func);
147 
155  static Domain ofElements(const mesh_t& mesh, std::function<bool(std::vector<vec2>, int)> func);
156 
158  static Domain ofElements(const mesh_t& mesh, std::function<bool(std::vector<vec3>, int)> func);
159 
165  static Domain ofBoundaryElements(const mesh_t& mesh, std::function<bool(std::vector<vec2>, int)> func);
166 
168  static Domain ofBoundaryElements(const mesh_t& mesh, std::function<bool(std::vector<vec3>, int)> func);
169 
171  const std::vector<int>& get(mfem::Geometry::Type geom) const
172  {
173  if (geom == mfem::Geometry::SEGMENT) return edge_ids_;
174  if (geom == mfem::Geometry::TRIANGLE) return tri_ids_;
175  if (geom == mfem::Geometry::SQUARE) return quad_ids_;
176  if (geom == mfem::Geometry::TETRAHEDRON) return tet_ids_;
177  if (geom == mfem::Geometry::CUBE) return hex_ids_;
178 
179  exit(1);
180  }
181 
183  const std::vector<int>& get_mfem_ids(mfem::Geometry::Type geom) const
184  {
185  if (geom == mfem::Geometry::SEGMENT) return mfem_edge_ids_;
186  if (geom == mfem::Geometry::TRIANGLE) return mfem_tri_ids_;
187  if (geom == mfem::Geometry::SQUARE) return mfem_quad_ids_;
188  if (geom == mfem::Geometry::TETRAHEDRON) return mfem_tet_ids_;
189  if (geom == mfem::Geometry::CUBE) return mfem_hex_ids_;
190 
191  exit(1);
192  }
193 
197  int total_elements() const
198  {
199  return int(edge_ids_.size() + tri_ids_.size() + quad_ids_.size() + tet_ids_.size() + hex_ids_.size());
200  }
201 
206  mfem::Array<int> bOffsets() const
207  {
208  mfem::Array<int> offsets(mfem::Geometry::NUM_GEOMETRIES + 1);
209 
210  int total = 0;
211  offsets[mfem::Geometry::POINT] = total;
212  total += 0; // vertices;
213  offsets[mfem::Geometry::SEGMENT] = total;
214  total += int(edge_ids_.size());
215  offsets[mfem::Geometry::TRIANGLE] = total;
216  total += int(tri_ids_.size());
217  offsets[mfem::Geometry::SQUARE] = total;
218  total += int(quad_ids_.size());
219  offsets[mfem::Geometry::TETRAHEDRON] = total;
220  total += int(tet_ids_.size());
221  offsets[mfem::Geometry::CUBE] = total;
222  total += int(hex_ids_.size());
223  offsets[mfem::Geometry::PRISM] = total;
224  offsets[mfem::Geometry::PYRAMID] = total;
225  offsets[mfem::Geometry::NUM_GEOMETRIES] = total;
226 
227  return offsets;
228  }
229 
231  mfem::Array<int> dof_list(const fes_t* fes) const;
232 
237  void insert_restriction(const fes_t* fes, FunctionSpace space);
238 
241 
247  void addElement(int geom_id, int elem_id, mfem::Geometry::Type element_geometry);
248 
254  void addElements(const std::vector<int>& geom_id, const std::vector<int>& elem_id,
255  mfem::Geometry::Type element_geometry);
256 
262 };
263 
265 Domain EntireDomain(const mesh_t& mesh);
266 
268 Domain EntireBoundary(const mesh_t& mesh);
269 
271 Domain InteriorFaces(const mesh_t& mesh);
272 
274 Domain operator|(const Domain& a, const Domain& b);
275 
277 Domain operator&(const Domain& a, const Domain& b);
278 
280 Domain operator-(const Domain& a, const Domain& b);
281 
283 template <int dim>
284 inline auto by_attr(int value)
285 {
286  return [value](std::vector<tensor<double, dim>>, int attr) { return value == attr; };
287 }
288 
290 template <int dim>
291 inline auto by_attr(std::set<int> values)
292 {
293  return [values](std::vector<tensor<double, dim>>, int attr) { return values.find(attr) != values.end(); };
294 }
295 
300 inline std::array<uint32_t, mfem::Geometry::NUM_GEOMETRIES> geometry_counts(const Domain& domain)
301 {
302  std::array<uint32_t, mfem::Geometry::NUM_GEOMETRIES> counts{};
303 
304  constexpr std::array<mfem::Geometry::Type, 5> geometries = {mfem::Geometry::SEGMENT, mfem::Geometry::TRIANGLE,
305  mfem::Geometry::SQUARE, mfem::Geometry::TETRAHEDRON,
306  mfem::Geometry::CUBE};
307  for (auto geom : geometries) {
308  counts[uint32_t(geom)] = uint32_t(domain.get(geom).size());
309  }
310  return counts;
311 }
312 
316 template <int dim>
317 inline tensor<double, dim> average(std::vector<tensor<double, dim>>& positions)
318 {
319  tensor<double, dim> total{};
320  for (auto x : positions) {
321  total += x;
322  }
323  return total / double(positions.size());
324 }
325 
326 } // namespace smith
This file contains helper traits and enumerations for classifying finite elements.
Accelerator functionality.
Definition: smith.cpp:36
tensor< double, dim > average(std::vector< tensor< double, dim >> &positions)
convenience function for computing the arithmetic mean of some list of vectors
Definition: domain.hpp:317
std::array< uint32_t, mfem::Geometry::NUM_GEOMETRIES > geometry_counts(const Domain &domain)
count the number of elements of each geometry in a domain
Definition: domain.hpp:300
Domain operator&(const Domain &a, const Domain &b)
create a new domain that is the intersection of a and b
Definition: domain.cpp:760
Domain operator|(const Domain &a, const Domain &b)
create a new domain that is the union of a and b
Definition: domain.cpp:759
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
Domain EntireBoundary(const mesh_t &mesh)
constructs a domain from all the boundary elements in a mesh
Definition: domain.cpp:609
Domain InteriorFaces(const mesh_t &mesh)
constructs a domain from all the interior face elements in a mesh
Definition: domain.cpp:625
auto by_attr(int value)
convenience predicate for creating domains by attribute
Definition: domain.hpp:284
Domain EntireDomain(const mesh_t &mesh)
constructs a domain from all the elements in a mesh
Definition: domain.cpp:594
Domain operator-(const Domain &a, const Domain &b)
create a new domain that is the set difference of a and b
Definition: domain.cpp:761
a generalization of mfem::ElementRestriction that works with multiple kinds of element geometries....
a class for representing a geometric region that can be used for integration
Definition: domain.hpp:33
int dim_
the geometric dimension of the domain
Definition: domain.hpp:48
std::vector< int > shared_interior_face_ids_
Ids of interior faces that lie on the boundary shared by two processors.
Definition: domain.hpp:106
const BlockElementRestriction & get_restriction(FunctionSpace space)
getter for accessing a restriction operator by its function space
Definition: domain.cpp:544
void insert_restriction(const fes_t *fes, FunctionSpace space)
create a restriction operator over this domain, using its FunctionSpace as a key
Definition: domain.cpp:536
static Domain ofElements(const mesh_t &mesh, std::function< bool(std::vector< vec2 >, int)> func)
create a domain from some subset of the elements (spatial dim == geometry dim) in an mfem::Mesh
Definition: domain.cpp:258
void addElement(int geom_id, int elem_id, mfem::Geometry::Type element_geometry)
Add an element to the domain.
Definition: domain.cpp:268
mfem::Array< int > dof_list(const fes_t *fes) const
get mfem degree of freedom list for a given FiniteElementSpace
Definition: domain.cpp:466
const std::vector< int > & get(mfem::Geometry::Type geom) const
get elements by geometry type
Definition: domain.hpp:171
std::map< FunctionSpace, BlockElementRestriction > restriction_operators
a collection of restriction operators for the different test/trial spaces appearing in integrals eval...
Definition: domain.hpp:103
int total_elements() const
returns how many elements of any type belong to this domain
Definition: domain.hpp:197
Type
enum describing what kind of elements are included in a Domain
Definition: domain.hpp:36
static Domain ofBoundaryElements(const mesh_t &mesh, std::function< bool(std::vector< vec2 >, int)> func)
create a domain from some subset of the boundary elements (spatial dim == geometry dim + 1) in an mfe...
Definition: domain.cpp:368
Type type_
whether the elements in this domain are on the boundary or not
Definition: domain.hpp:51
void addElements(const std::vector< int > &geom_id, const std::vector< int > &elem_id, mfem::Geometry::Type element_geometry)
Add a batch of elements to the domain.
Definition: domain.cpp:290
static Domain ofFaces(const mesh_t &mesh, std::function< bool(std::vector< vec2 >, int)> func)
create a domain from some subset of the faces in an mfem::Mesh
Definition: domain.cpp:182
static Domain ofVertices(const mesh_t &mesh, std::function< bool(vec3)> func)
This is an overloaded member function, provided for convenience. It differs from the above function o...
void insert_shared_interior_face_list()
Find the list of interior faces shared by two processors to make sure these faces are only integrated...
Definition: domain.cpp:546
Domain(const mesh_t &m, int d, Type type)
empty Domain constructor, with connectivity info to be populated later
Definition: domain.hpp:111
static Domain ofEdges(const mesh_t &mesh, std::function< bool(std::vector< vec2 >, int)> func)
create a domain from some subset of the edges in an mfem::Mesh
Definition: domain.cpp:102
const std::vector< int > & get_mfem_ids(mfem::Geometry::Type geom) const
get elements by geometry type
Definition: domain.hpp:183
static constexpr int num_types
the number of entries in the Type enum
Definition: domain.hpp:42
mfem::Array< int > bOffsets() const
returns an array of the prefix sum of element counts belonging to this domain. Primarily intended to ...
Definition: domain.hpp:206
static Domain ofVertices(const mesh_t &mesh, std::function< bool(vec2)> func)
create a domain from some subset of the vertices in an mfem::Mesh
const mesh_t & mesh_
the underyling mesh for this domain
Definition: domain.hpp:45
a small POD class for tracking function space metadata
Arbitrary-rank tensor class.
Definition: tensor.hpp:28
Implementation of the tensor class used by Functional.