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;
   115                     std::size_t sub_domain,
   116                     bool check_midpoint=
true) 
const;
   127                      std::size_t sub_domain,
   128                      bool check_midpoint=
true) 
const;
   143               std::size_t sub_domain,
   144               bool check_midpoint=
true) 
const;
   157               std::size_t sub_domain,
   158               bool check_midpoint=
true) 
const;
   170               bool check_midpoint=
true) 
const;
   182               bool check_midpoint=
true) 
const;
   194               bool check_midpoint=
true) 
const;
   209               std::size_t sub_domain,
   211               bool check_midpoint=
true) 
const;
   226               bool check_midpoint=
true) 
const;
   241               bool check_midpoint=
true) 
const;
   256               bool check_midpoint=
true) 
const;
   268     virtual void set_property(std::string name, 
double value);
   286     template<
typename S, 
typename T>
   287     void apply_markers(S& sub_domains,
   290                        bool check_midpoint) 
const;
   293       void apply_markers(std::map<std::size_t, std::size_t>& sub_domains,
   297                          bool check_midpoint) 
const;
   301     friend class PeriodicBC;
   305     mutable std::size_t _geometric_dimension;
 
virtual ~SubDomain()
Destructor. 
Definition: SubDomain.cpp:46
 
Definition: SubDomain.h:42
 
void mark(Mesh &mesh, std::size_t dim, std::size_t sub_domain, bool check_midpoint=true) const
Definition: SubDomain.cpp:94
 
virtual bool inside(const Array< double > &x, bool on_boundary) const
Definition: SubDomain.cpp:51
 
void mark_facets(Mesh &mesh, std::size_t sub_domain, bool check_midpoint=true) const
Definition: SubDomain.cpp:87
 
const double map_tolerance
Definition: SubDomain.h:280
 
virtual void map(const Array< double > &x, Array< double > &y) const
Definition: SubDomain.cpp:65
 
SubDomain(const double map_tol=1.0e-10)
Definition: SubDomain.cpp:40
 
virtual void snap(Array< double > &x) const
Definition: SubDomain.h:102
 
Interface for setting (strong) Dirichlet boundary conditions. 
Definition: DirichletBC.h:124
 
Definition: GenericFile.h:38
 
virtual void set_property(std::string name, double value)
Definition: SubDomain.cpp:381
 
std::size_t geometric_dimension() const
Definition: SubDomain.cpp:165
 
void mark_cells(Mesh &mesh, std::size_t sub_domain, bool check_midpoint=true) const
Definition: SubDomain.cpp:80
 
virtual double get_property(std::string name) const
Definition: SubDomain.cpp:388