DOLFIN
DOLFIN C++ interface
Public Member Functions | Public Attributes | Friends | List of all members
dolfin::SubDomain Class Reference

#include <SubDomain.h>

Inheritance diagram for dolfin::SubDomain:
Inheritance graph
[legend]

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
 

Detailed Description

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.

Constructor & Destructor Documentation

◆ SubDomain()

SubDomain::SubDomain ( const double  map_tol = 1.0e-10)

Constructor

Parameters
map_tol(double) The tolerance used when identifying mapped points using the function SubDomain::map.

Member Function Documentation

◆ geometric_dimension()

std::size_t SubDomain::geometric_dimension ( ) const

Return geometric dimension

Returns
std::size_t The geometric dimension.

◆ get_property()

double SubDomain::get_property ( std::string  name) const
virtual

Property getter

Parameters
name
Returns
double

◆ inside() [1/2]

bool SubDomain::inside ( const Array< double > &  x,
bool  on_boundary 
) const
virtual

Return true for points inside the subdomain

Parameters
x(Array<double>) The coordinates of the point.
on_boundary(bool) True for points on the boundary.
Returns
bool True for points inside the subdomain.

Reimplemented in dolfin::DomainBoundary.

◆ inside() [2/2]

bool SubDomain::inside ( Eigen::Ref< const Eigen::VectorXd >  x,
bool  on_boundary 
) const
virtual

Return true for points inside the subdomain

Parameters
x(Eigen::Ref<const Eigen::VectorXd>) The coordinates of the point.
on_boundary(bool) True for points on the boundary.
Returns
bool True for points inside the subdomain.

Reimplemented in EmptySubDomain.

◆ map() [1/2]

void SubDomain::map ( const Array< double > &  x,
Array< double > &  y 
) const
virtual

Map coordinate x in domain H to coordinate y in domain G (used for periodic boundary conditions)

Parameters
x(Array<double>) The coordinates in domain H.
y(Array<double>) The coordinates in domain G.

◆ map() [2/2]

void SubDomain::map ( Eigen::Ref< const Eigen::VectorXd >  x,
Eigen::Ref< Eigen::VectorXd >  y 
) const
virtual

Map coordinate x in domain H to coordinate y in domain G (used for periodic boundary conditions)

Parameters
x(Eigen::Ref<const Eigen::VectorXd>) The coordinates in domain H.
y(Eigen::Ref<Eigen::VectorXd>) The coordinates in domain G.

◆ mark() [1/9]

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

Parameters
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).

◆ mark() [2/9]

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

Parameters
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).

◆ mark() [3/9]

void SubDomain::mark ( MeshFunction< int > &  sub_domains,
int  sub_domain,
bool  check_midpoint = true 
) const

Set subdomain markers (int) for given subdomain number

Parameters
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).

◆ mark() [4/9]

void SubDomain::mark ( MeshFunction< double > &  sub_domains,
double  sub_domain,
bool  check_midpoint = true 
) const

Set subdomain markers (double) for given subdomain number

Parameters
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).

◆ mark() [5/9]

void SubDomain::mark ( MeshFunction< bool > &  sub_domains,
bool  sub_domain,
bool  check_midpoint = true 
) const

Set subdomain markers (bool) for given subdomain

Parameters
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).

◆ mark() [6/9]

void SubDomain::mark ( MeshValueCollection< std::size_t > &  sub_domains,
std::size_t  sub_domain,
const Mesh mesh,
bool  check_midpoint = true 
) const

Set subdomain markers (std::size_t) for given subdomain number

Parameters
sub_domains(MeshValueCollection<std::size_t>) The subdomain markers.
sub_domain(std::size_t) The subdomain number.
mesh(Mesh) The mesh.
check_midpoint(bool) Flag for whether midpoint of cell should be checked (default).

◆ mark() [7/9]

void SubDomain::mark ( MeshValueCollection< int > &  sub_domains,
int  sub_domain,
const Mesh mesh,
bool  check_midpoint = true 
) const

Set subdomain markers (int) for given subdomain number

Parameters
sub_domains(MeshValueCollection<int>) The subdomain markers
sub_domain(int) The subdomain number
mesh(Mesh) The mesh.
check_midpoint(bool) Flag for whether midpoint of cell should be checked (default).

◆ mark() [8/9]

void SubDomain::mark ( MeshValueCollection< double > &  sub_domains,
double  sub_domain,
const Mesh mesh,
bool  check_midpoint = true 
) const

Set subdomain markers (double) for given subdomain number

Parameters
sub_domains(MeshValueCollection<double>) The subdomain markers.
sub_domain(double) The subdomain number
mesh(Mesh) The mesh.
check_midpoint(bool) Flag for whether midpoint of cell should be checked (default).

◆ mark() [9/9]

void SubDomain::mark ( MeshValueCollection< bool > &  sub_domains,
bool  sub_domain,
const Mesh mesh,
bool  check_midpoint = true 
) const

Set subdomain markers (bool) for given subdomain

Parameters
sub_domains(MeshValueCollection<bool>) The subdomain markers
sub_domain(bool) The subdomain number
mesh(Mesh) The mesh.
check_midpoint(bool) Flag for whether midpoint of cell should be checked (default).

◆ mark_cells()

void SubDomain::mark_cells ( Mesh mesh,
std::size_t  sub_domain,
bool  check_midpoint = true 
) const

Set subdomain markers (std::size_t) on cells for given subdomain number

Parameters
mesh(Mesh) The mesh to be marked.
sub_domain(std::size_t) The subdomain number.
check_midpoint(bool) Flag for whether midpoint of cell should be checked (default).

◆ mark_facets()

void SubDomain::mark_facets ( Mesh mesh,
std::size_t  sub_domain,
bool  check_midpoint = true 
) const

Set subdomain markers (std::size_t) on facets for given subdomain number

Parameters
mesh(Mesh) The mesh to be marked.
sub_domain(std::size_t) The subdomain number.
check_midpoint(bool) Flag for whether midpoint of cell should be checked (default).

◆ set_property()

void SubDomain::set_property ( std::string  name,
double  value 
)
virtual

Property setter

Parameters
name
value

◆ snap()

virtual void dolfin::SubDomain::snap ( Array< double > &  x) const
inlinevirtual

Snap coordinate to boundary of subdomain

Parameters
x(Array<double>) The coordinates.

Member Data Documentation

◆ map_tolerance

const double dolfin::SubDomain::map_tolerance

Return tolerance uses to find matching point via map function

Returns
double The tolerance.

The documentation for this class was generated from the following files: