| 
    DOLFIN
    
   DOLFIN C++ interface 
   | 
 
#include <SubDomain.h>

Public Member Functions | |
| SubDomain (const double map_tol=1.0e-10) | |
| virtual | ~SubDomain () | 
| Destructor.  | |
| virtual bool | inside (const Array< double > &x, bool on_boundary) const | 
| virtual bool | inside (Eigen::Ref< const Eigen::VectorXd > x, bool on_boundary) const | 
| virtual void | map (const Array< double > &x, Array< double > &y) const | 
| virtual void | map (Eigen::Ref< const Eigen::VectorXd > x, Eigen::Ref< Eigen::VectorXd > y) const | 
| virtual void | snap (Array< double > &x) const | 
| void | mark_cells (Mesh &mesh, std::size_t sub_domain, bool check_midpoint=true) const | 
| void | mark_facets (Mesh &mesh, std::size_t sub_domain, bool check_midpoint=true) const | 
| void | mark (Mesh &mesh, std::size_t dim, std::size_t sub_domain, bool check_midpoint=true) const | 
| void | mark (MeshFunction< std::size_t > &sub_domains, std::size_t sub_domain, bool check_midpoint=true) const | 
| void | mark (MeshFunction< int > &sub_domains, int sub_domain, bool check_midpoint=true) const | 
| void | mark (MeshFunction< double > &sub_domains, double sub_domain, bool check_midpoint=true) const | 
| void | mark (MeshFunction< bool > &sub_domains, bool sub_domain, bool check_midpoint=true) const | 
| void | mark (MeshValueCollection< std::size_t > &sub_domains, std::size_t sub_domain, const Mesh &mesh, bool check_midpoint=true) const | 
| void | mark (MeshValueCollection< int > &sub_domains, int sub_domain, const Mesh &mesh, bool check_midpoint=true) const | 
| void | mark (MeshValueCollection< double > &sub_domains, double sub_domain, const Mesh &mesh, bool check_midpoint=true) const | 
| void | mark (MeshValueCollection< bool > &sub_domains, bool sub_domain, const Mesh &mesh, bool check_midpoint=true) const | 
| std::size_t | geometric_dimension () const | 
| virtual void | set_property (std::string name, double value) | 
| virtual double | get_property (std::string name) const | 
Public Attributes | |
| const double | map_tolerance | 
Friends | |
| class | DirichletBC | 
| class | PeriodicBC | 
This class defines the interface for definition of subdomains. Alternatively, subdomains may be defined by a Mesh and a MeshFunction<std::size_t> over the mesh.
| SubDomain::SubDomain | ( | const double | map_tol = 1.0e-10 | ) | 
Constructor
| map_tol | (double) The tolerance used when identifying mapped points using the function SubDomain::map. | 
| std::size_t SubDomain::geometric_dimension | ( | ) | const | 
Return geometric dimension
      
  | 
  virtual | 
Property getter
| name | 
      
  | 
  virtual | 
Return true for points inside the subdomain
| x | (Array<double>) The coordinates of the point. | 
| on_boundary | (bool) True for points on the boundary. | 
Reimplemented in dolfin::DomainBoundary.
      
  | 
  virtual | 
Return true for points inside the subdomain
| x | (Eigen::Ref<const Eigen::VectorXd>) The coordinates of the point. | 
| on_boundary | (bool) True for points on the boundary. | 
Reimplemented in EmptySubDomain.
Map coordinate x in domain H to coordinate y in domain G (used for periodic boundary conditions)
| x | (Array<double>) The coordinates in domain H. | 
| y | (Array<double>) The coordinates in domain G. | 
      
  | 
  virtual | 
Map coordinate x in domain H to coordinate y in domain G (used for periodic boundary conditions)
| x | (Eigen::Ref<const Eigen::VectorXd>) The coordinates in domain H. | 
| y | (Eigen::Ref<Eigen::VectorXd>) The coordinates in domain G. | 
| void SubDomain::mark | ( | Mesh & | mesh, | 
| std::size_t | dim, | ||
| std::size_t | sub_domain, | ||
| bool | check_midpoint = true  | 
        ||
| ) | const | 
Set subdomain markers (std::size_t) for given topological dimension and subdomain number
| mesh | (Mesh) The mesh to be marked. | 
| dim | (std::size_t) The topological dimension of entities to be marked. | 
| sub_domain | (std::size_t) The subdomain number. | 
| check_midpoint | (bool) Flag for whether midpoint of cell should be checked (default). | 
| void SubDomain::mark | ( | MeshFunction< std::size_t > & | sub_domains, | 
| std::size_t | sub_domain, | ||
| bool | check_midpoint = true  | 
        ||
| ) | const | 
Set subdomain markers (std::size_t) for given subdomain number
| sub_domains | (MeshFunction<std::size_t>) The subdomain markers. | 
| sub_domain | (std::size_t) The subdomain number. | 
| check_midpoint | (bool) Flag for whether midpoint of cell should be checked (default). | 
| void SubDomain::mark | ( | MeshFunction< int > & | sub_domains, | 
| int | sub_domain, | ||
| bool | check_midpoint = true  | 
        ||
| ) | const | 
Set subdomain markers (int) for given subdomain number
| sub_domains | (MeshFunction<int>) The subdomain markers. | 
| sub_domain | (int) The subdomain number. | 
| check_midpoint | (bool) Flag for whether midpoint of cell should be checked (default). | 
| void SubDomain::mark | ( | MeshFunction< double > & | sub_domains, | 
| double | sub_domain, | ||
| bool | check_midpoint = true  | 
        ||
| ) | const | 
Set subdomain markers (double) for given subdomain number
| sub_domains | (MeshFunction<double>) The subdomain markers. | 
| sub_domain | (double) The subdomain number. | 
| check_midpoint | (bool) Flag for whether midpoint of cell should be checked (default). | 
| void SubDomain::mark | ( | MeshFunction< bool > & | sub_domains, | 
| bool | sub_domain, | ||
| bool | check_midpoint = true  | 
        ||
| ) | const | 
Set subdomain markers (bool) for given subdomain
| sub_domains | (MeshFunction<bool>) The subdomain markers. | 
| sub_domain | (bool) The subdomain number. | 
| check_midpoint | (bool) Flag for whether midpoint of cell should be checked (default). | 
| void SubDomain::mark | ( | MeshValueCollection< std::size_t > & | sub_domains, | 
| std::size_t | sub_domain, | ||
| const Mesh & | mesh, | ||
| bool | check_midpoint = true  | 
        ||
| ) | const | 
| void SubDomain::mark | ( | MeshValueCollection< int > & | sub_domains, | 
| int | sub_domain, | ||
| const Mesh & | mesh, | ||
| bool | check_midpoint = true  | 
        ||
| ) | const | 
| void SubDomain::mark | ( | MeshValueCollection< double > & | sub_domains, | 
| double | sub_domain, | ||
| const Mesh & | mesh, | ||
| bool | check_midpoint = true  | 
        ||
| ) | const | 
| void SubDomain::mark | ( | MeshValueCollection< bool > & | sub_domains, | 
| bool | sub_domain, | ||
| const Mesh & | mesh, | ||
| bool | check_midpoint = true  | 
        ||
| ) | const | 
| void SubDomain::mark_cells | ( | Mesh & | mesh, | 
| std::size_t | sub_domain, | ||
| bool | check_midpoint = true  | 
        ||
| ) | const | 
| void SubDomain::mark_facets | ( | Mesh & | mesh, | 
| std::size_t | sub_domain, | ||
| bool | check_midpoint = true  | 
        ||
| ) | const | 
      
  | 
  virtual | 
Property setter
| name | |
| value | 
      
  | 
  inlinevirtual | 
Snap coordinate to boundary of subdomain
| x | (Array<double>) The coordinates. | 
| const double dolfin::SubDomain::map_tolerance | 
Return tolerance uses to find matching point via map function
 1.8.13