Serac  0.1
Serac is an implicit thermal strucural mechanics simulation code.
domain.hpp
Go to the documentation of this file.
1 // Copyright (c) 2019-2024, Lawrence Livermore National Security, LLC and
2 // other Serac 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 <vector>
10 #include "mfem.hpp"
11 
13 
14 namespace serac {
15 
21 struct Domain {
23  enum Type
24  {
25  Elements,
26  BoundaryElements
27  };
28 
29  static constexpr int num_types = 2;
30 
32  const mfem::Mesh& mesh_;
33 
35  int dim_;
36 
39 
47  std::vector<int> vertex_ids_;
48  std::vector<int> edge_ids_;
49  std::vector<int> tri_ids_;
50  std::vector<int> quad_ids_;
51  std::vector<int> tet_ids_;
52  std::vector<int> hex_ids_;
54 
55  Domain(const mfem::Mesh& m, int d, Type type = Domain::Type::Elements) : mesh_(m), dim_(d), type_(type) {}
56 
63  static Domain ofVertices(const mfem::Mesh& mesh, std::function<bool(vec2)> func);
64 
66  static Domain ofVertices(const mfem::Mesh& mesh, std::function<bool(vec3)> func);
67 
75  static Domain ofEdges(const mfem::Mesh& mesh, std::function<bool(std::vector<vec2>, int)> func);
76 
78  static Domain ofEdges(const mfem::Mesh& mesh, std::function<bool(std::vector<vec3>)> func);
79 
87  static Domain ofFaces(const mfem::Mesh& mesh, std::function<bool(std::vector<vec2>, int)> func);
88 
90  static Domain ofFaces(const mfem::Mesh& mesh, std::function<bool(std::vector<vec3>, int)> func);
91 
99  static Domain ofElements(const mfem::Mesh& mesh, std::function<bool(std::vector<vec2>, int)> func);
100 
102  static Domain ofElements(const mfem::Mesh& mesh, std::function<bool(std::vector<vec3>, int)> func);
103 
109  static Domain ofBoundaryElements(const mfem::Mesh& mesh, std::function<bool(std::vector<vec2>, int)> func);
110 
112  static Domain ofBoundaryElements(const mfem::Mesh& mesh, std::function<bool(std::vector<vec3>, int)> func);
113 
115  const std::vector<int>& get(mfem::Geometry::Type geom) const
116  {
117  if (geom == mfem::Geometry::POINT) return vertex_ids_;
118  if (geom == mfem::Geometry::SEGMENT) return edge_ids_;
119  if (geom == mfem::Geometry::TRIANGLE) return tri_ids_;
120  if (geom == mfem::Geometry::SQUARE) return quad_ids_;
121  if (geom == mfem::Geometry::TETRAHEDRON) return tet_ids_;
122  if (geom == mfem::Geometry::CUBE) return hex_ids_;
123 
124  exit(1);
125  }
126 };
127 
129 Domain EntireDomain(const mfem::Mesh& mesh);
130 
132 Domain EntireBoundary(const mfem::Mesh& mesh);
133 
135 Domain operator|(const Domain& a, const Domain& b);
136 
138 Domain operator&(const Domain& a, const Domain& b);
139 
141 Domain operator-(const Domain& a, const Domain& b);
142 
144 template <int dim>
145 inline auto by_attr(int value)
146 {
147  return [value](std::vector<tensor<double, dim> >, int attr) { return attr == value; };
148 }
149 
150 } // namespace serac
Accelerator functionality.
Definition: serac.cpp:38
Domain EntireBoundary(const mfem::Mesh &mesh)
constructs a domain from all the boundary elements in a mesh
Definition: domain.cpp:418
auto by_attr(int value)
convenience predicate for creating domains by attribute
Definition: domain.hpp:145
Domain operator-(const Domain &a, const Domain &b)
create a new domain that is the set difference of a and b
Definition: domain.cpp:500
Domain operator&(const Domain &a, const Domain &b)
create a new domain that is the intersection of a and b
Definition: domain.cpp:499
Domain operator|(const Domain &a, const Domain &b)
create a new domain that is the union of a and b
Definition: domain.cpp:498
constexpr SERAC_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:274
Domain EntireDomain(const mfem::Mesh &mesh)
constructs a domain from all the elements in a mesh
Definition: domain.cpp:382
a class for representing a geometric region that can be used for integration
Definition: domain.hpp:21
static Domain ofFaces(const mfem::Mesh &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:213
static Domain ofBoundaryElements(const mfem::Mesh &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:367
static constexpr int num_types
the number of entries in the Type enum
Definition: domain.hpp:29
static Domain ofVertices(const mfem::Mesh &mesh, std::function< bool(vec2)> func)
create a domain from some subset of the vertices in an mfem::Mesh
Definition: domain.cpp:79
Domain(const mfem::Mesh &m, int d, Type type=Domain::Type::Elements)
Definition: domain.hpp:55
Type type_
whether the elements in this domain are on the boundary or not
Definition: domain.hpp:38
Type
enum describing what kind of elements are included in a Domain
Definition: domain.hpp:24
int dim_
the geometric dimension of the domain
Definition: domain.hpp:35
const mfem::Mesh & mesh_
the underyling mesh for this domain
Definition: domain.hpp:32
const std::vector< int > & get(mfem::Geometry::Type geom) const
get elements by geometry type
Definition: domain.hpp:115
static Domain ofElements(const mfem::Mesh &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:290
static Domain ofEdges(const mfem::Mesh &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:132
Arbitrary-rank tensor class.
Definition: tensor.hpp:29
Implementation of the tensor class used by Functional.