21 #ifndef __SUB_DOMAIN_H
22 #define __SUB_DOMAIN_H
26 #include <dolfin/common/constants.h>
27 #include <Eigen/Dense>
34 template <
typename T>
class MeshFunction;
35 template <
typename T>
class MeshValueCollection;
36 template <
typename T>
class Array;
76 virtual bool inside(Eigen::Ref<const Eigen::VectorXd> x,
bool on_boundary)
const;
95 virtual void map(Eigen::Ref<const Eigen::VectorXd> x, Eigen::Ref<Eigen::VectorXd> y)
const;
108 virtual void snap(Eigen::Ref<Eigen::VectorXd> x)
const;
121 std::size_t sub_domain,
122 bool check_midpoint=
true)
const;
133 std::size_t sub_domain,
134 bool check_midpoint=
true)
const;
149 std::size_t sub_domain,
150 bool check_midpoint=
true)
const;
163 std::size_t sub_domain,
164 bool check_midpoint=
true)
const;
176 bool check_midpoint=
true)
const;
188 bool check_midpoint=
true)
const;
200 bool check_midpoint=
true)
const;
215 std::size_t sub_domain,
217 bool check_midpoint=
true)
const;
232 bool check_midpoint=
true)
const;
247 bool check_midpoint=
true)
const;
262 bool check_midpoint=
true)
const;
274 virtual void set_property(std::string name,
double value);
292 template<
typename S,
typename T>
293 void apply_markers(S& sub_domains,
296 bool check_midpoint)
const;
299 void apply_markers(std::map<std::size_t, std::size_t>& sub_domains,
303 bool check_midpoint)
const;
307 friend class PeriodicBC;
311 mutable std::size_t _geometric_dimension;
Interface for setting (strong) Dirichlet boundary conditions.
Definition: DirichletBC.h:125
Definition: MeshValueCollection.h:51
Definition: SubDomain.h:43
virtual double get_property(std::string name) const
Definition: SubDomain.cpp:401
SubDomain(const double map_tol=1.0e-10)
Definition: SubDomain.cpp:40
void mark(Mesh &mesh, std::size_t dim, std::size_t sub_domain, bool check_midpoint=true) const
Definition: SubDomain.cpp:107
virtual ~SubDomain()
Destructor.
Definition: SubDomain.cpp:46
virtual bool inside(const Array< double > &x, bool on_boundary) const
Definition: SubDomain.cpp:51
virtual void set_property(std::string name, double value)
Definition: SubDomain.cpp:394
std::size_t geometric_dimension() const
Definition: SubDomain.cpp:178
const double map_tolerance
Definition: SubDomain.h:286
void mark_facets(Mesh &mesh, std::size_t sub_domain, bool check_midpoint=true) const
Definition: SubDomain.cpp:100
virtual void snap(Array< double > &x) const
Definition: SubDomain.cpp:80
void mark_cells(Mesh &mesh, std::size_t sub_domain, bool check_midpoint=true) const
Definition: SubDomain.cpp:93
virtual void map(const Array< double > &x, Array< double > &y) const
Definition: SubDomain.cpp:65