DOLFIN
DOLFIN C++ interface
Classes | Public Types | Public Member Functions | Static Public Member Functions | Public Attributes | List of all members
dolfin::DirichletBC Class Reference

Interface for setting (strong) Dirichlet boundary conditions. More...

#include <DirichletBC.h>

Inheritance diagram for dolfin::DirichletBC:
Inheritance graph
[legend]
Collaboration diagram for dolfin::DirichletBC:
Collaboration graph
[legend]

Public Types

typedef std::unordered_map< std::size_t, double > Map
 map type used by DirichletBC
 

Public Member Functions

 DirichletBC (std::shared_ptr< const FunctionSpace > V, std::shared_ptr< const GenericFunction > g, std::shared_ptr< const SubDomain > sub_domain, std::string method="topological", bool check_midpoint=true)
 
 DirichletBC (std::shared_ptr< const FunctionSpace > V, std::shared_ptr< const GenericFunction > g, std::shared_ptr< const MeshFunction< std::size_t >> sub_domains, std::size_t sub_domain, std::string method="topological")
 
 DirichletBC (std::shared_ptr< const FunctionSpace > V, std::shared_ptr< const GenericFunction > g, std::size_t sub_domain, std::string method="topological")
 
 DirichletBC (std::shared_ptr< const FunctionSpace > V, std::shared_ptr< const GenericFunction > g, const std::vector< std::size_t > &markers, std::string method="topological")
 
 DirichletBC (const DirichletBC &bc)
 
 ~DirichletBC ()
 Destructor.
 
const DirichletBCoperator= (const DirichletBC &bc)
 
void apply (GenericMatrix &A) const
 
void apply (GenericVector &b) const
 
void apply (GenericMatrix &A, GenericVector &b) const
 
void apply (GenericVector &b, const GenericVector &x) const
 
void apply (GenericMatrix &A, GenericVector &b, const GenericVector &x) const
 
void get_boundary_values (Map &boundary_values) const
 
void gather (Map &boundary_values) const
 
void zero (GenericMatrix &A) const
 
void zero_columns (GenericMatrix &A, GenericVector &b, double diag_val=0) const
 
const std::vector< std::size_t > & markers () const
 
std::shared_ptr< const FunctionSpacefunction_space () const
 
std::shared_ptr< const GenericFunctionvalue () const
 
std::shared_ptr< const SubDomainuser_sub_domain () const
 
void set_value (std::shared_ptr< const GenericFunction > g)
 
void homogenize ()
 Set value to 0.0.
 
std::string method () const
 
- Public Member Functions inherited from dolfin::Hierarchical< DirichletBC >
 Hierarchical (DirichletBC &self)
 Constructor.
 
virtual ~Hierarchical ()
 Destructor.
 
std::size_t depth () const
 
bool has_parent () const
 
bool has_child () const
 
DirichletBCparent ()
 
const DirichletBCparent () const
 Return parent in hierarchy (const version).
 
std::shared_ptr< DirichletBCparent_shared_ptr ()
 
std::shared_ptr< const DirichletBCparent_shared_ptr () const
 Return shared pointer to parent (const version).
 
DirichletBCchild ()
 
const DirichletBCchild () const
 Return child in hierarchy (const version).
 
std::shared_ptr< DirichletBCchild_shared_ptr ()
 
std::shared_ptr< const DirichletBCchild_shared_ptr () const
 Return shared pointer to child (const version).
 
DirichletBCroot_node ()
 
const DirichletBCroot_node () const
 Return root node object in hierarchy (const version).
 
std::shared_ptr< DirichletBCroot_node_shared_ptr ()
 
std::shared_ptr< const DirichletBCroot_node_shared_ptr () const
 Return shared pointer to root node object in hierarchy (const version).
 
DirichletBCleaf_node ()
 
const DirichletBCleaf_node () const
 Return leaf node object in hierarchy (const version).
 
std::shared_ptr< DirichletBCleaf_node_shared_ptr ()
 
std::shared_ptr< const DirichletBCleaf_node_shared_ptr () const
 Return shared pointer to leaf node object in hierarchy (const version).
 
void set_parent (std::shared_ptr< DirichletBC > parent)
 Set parent.
 
void clear_child ()
 Clear child.
 
void set_child (std::shared_ptr< DirichletBC > child)
 Set child.
 
const Hierarchicaloperator= (const Hierarchical &hierarchical)
 Assignment operator.
 
void _debug () const
 Function useful for debugging the hierarchy.
 
- Public Member Functions inherited from dolfin::Variable
 Variable ()
 Create unnamed variable.
 
 Variable (const std::string name, const std::string label)
 Create variable with given name and label.
 
 Variable (const Variable &variable)
 Copy constructor.
 
virtual ~Variable ()
 Destructor.
 
const Variableoperator= (const Variable &variable)
 Assignment operator.
 
void rename (const std::string name, const std::string label)
 Rename variable.
 
std::string name () const
 Return name.
 
std::string label () const
 Return label (description)
 
std::size_t id () const
 
virtual std::string str (bool verbose) const
 Return informal string representation (pretty-print)
 

Static Public Member Functions

static Parameters default_parameters ()
 

Public Attributes

std::shared_ptr< const SubDomain_user_sub_domain
 
- Public Attributes inherited from dolfin::Variable
Parameters parameters
 Parameters.
 

Detailed Description

Interface for setting (strong) Dirichlet boundary conditions.

u = g on G,

where u is the solution to be computed, g is a function and G is a sub domain of the mesh.

A DirichletBC is specified by the function g, the function space (trial space) and boundary indicators on (a subset of) the mesh boundary.

The boundary indicators may be specified in a number of different ways.

The simplest approach is to specify the boundary by a SubDomain object, using the inside() function to specify on which facets the boundary conditions should be applied. The boundary facets will then be searched for and marked only on the first call to apply. This means that the mesh could be moved after the first apply and the boundary markers would still remain intact.

Alternatively, the boundary may be specified by a MeshFunction over facets labeling all mesh facets together with a number that specifies which facets should be included in the boundary.

The third option is to attach the boundary information to the mesh. This is handled automatically when exporting a mesh from for example VMTK.

The 'method' variable may be used to specify the type of method used to identify degrees of freedom on the boundary. Available methods are: topological approach (default), geometric approach, and pointwise approach. The topological approach is faster, but will only identify degrees of freedom that are located on a facet that is entirely on the boundary. In particular, the topological approach will not identify degrees of freedom for discontinuous elements (which are all internal to the cell). A remedy for this is to use the geometric approach. In the geometric approach, each dof on each facet that matches the boundary condition will be checked. To apply pointwise boundary conditions e.g. pointloads, one will have to use the pointwise approach. The three possibilities are "topological", "geometric" and "pointwise".

Note: when using "pointwise", the boolean argument on_boundary in SubDomain::inside will always be false.

The 'check_midpoint' variable can be used to decide whether or not the midpoint of each facet should be checked when a user-defined SubDomain is used to define the domain of the boundary condition. By default, midpoints are always checked. Note that this variable may be of importance close to corners, in which case it is sometimes important to check the midpoint to avoid including facets "on the diagonal close" to a corner. This variable is also of importance for curved boundaries (like on a sphere or cylinder), in which case it is important not to check the midpoint which will be located in the interior of a domain defined relative to a radius.

Note that there may be caching employed in BC computation for performance reasons. In particular, applicable DOFs are cached by some methods on a first apply(). This means that changing a supplied object (defining boundary subdomain) after first use may have no effect. But this is implementation and method specific.

Constructor & Destructor Documentation

◆ DirichletBC() [1/5]

DirichletBC::DirichletBC ( std::shared_ptr< const FunctionSpace V,
std::shared_ptr< const GenericFunction g,
std::shared_ptr< const SubDomain sub_domain,
std::string  method = "topological",
bool  check_midpoint = true 
)

Create boundary condition for subdomain

Parameters
[in]V(FunctionSpace) The function space
[in]g(GenericFunction) The value
[in]sub_domain(SubDomain) The subdomain
[in]method(std::string) Optional argument: A string specifying the method to identify dofs
[in]check_midpoint(bool)

◆ DirichletBC() [2/5]

DirichletBC::DirichletBC ( std::shared_ptr< const FunctionSpace V,
std::shared_ptr< const GenericFunction g,
std::shared_ptr< const MeshFunction< std::size_t >>  sub_domains,
std::size_t  sub_domain,
std::string  method = "topological" 
)

Create boundary condition for subdomain specified by index

Parameters
[in]V(FunctionSpace) The function space.
[in]g(GenericFunction) The value.
[in]sub_domains(MeshFnunction<std::size_t>) Subdomain markers
[in]sub_domain(std::size_t) The subdomain index (number)
[in]method(std::string) Optional argument: A string specifying the method to identify dofs.

◆ DirichletBC() [3/5]

DirichletBC::DirichletBC ( std::shared_ptr< const FunctionSpace V,
std::shared_ptr< const GenericFunction g,
std::size_t  sub_domain,
std::string  method = "topological" 
)

Create boundary condition for boundary data included in the mesh

Parameters
[in]V(FunctionSpace) The function space.
[in]g(GenericFunction) The value.
[in]sub_domain(std::size_t) The subdomain index (number)
[in]method(std::string) Optional argument: A string specifying the method to identify dofs.

◆ DirichletBC() [4/5]

DirichletBC::DirichletBC ( std::shared_ptr< const FunctionSpace V,
std::shared_ptr< const GenericFunction g,
const std::vector< std::size_t > &  markers,
std::string  method = "topological" 
)

Create boundary condition for subdomain by boundary markers (cells, local facet numbers)

Parameters
[in]V(FunctionSpace) The function space.
[in]g(GenericFunction) The value.
[in]markers(std::vector<std:size_t>&) Subdomain markers (facet index local to process)
[in]method(std::string) Optional argument: A string specifying the method to identify dofs.

◆ DirichletBC() [5/5]

DirichletBC::DirichletBC ( const DirichletBC bc)

Copy constructor. Either cached DOF data are copied.

Parameters
[in]bc(DirichletBC&) The object to be copied.

Member Function Documentation

◆ apply() [1/5]

void DirichletBC::apply ( GenericMatrix A) const

Apply boundary condition to a matrix

Parameters
[in,out]A(GenericMatrix) The matrix to apply boundary condition to.

◆ apply() [2/5]

void DirichletBC::apply ( GenericVector b) const

Apply boundary condition to a vector

Parameters
[in,out]b(GenericVector) The vector to apply boundary condition to.

◆ apply() [3/5]

void DirichletBC::apply ( GenericMatrix A,
GenericVector b 
) const

Apply boundary condition to a linear system

Parameters
[in,out]A(GenericMatrix) The matrix to apply boundary condition to.
[in,out]b(GenericVector) The vector to apply boundary condition to.

◆ apply() [4/5]

void DirichletBC::apply ( GenericVector b,
const GenericVector x 
) const

Apply boundary condition to vectors for a nonlinear problem

Parameters
[in,out]b(GenericVector) The vector to apply boundary conditions to.
[in]x(GenericVector) Another vector (nonlinear problem).

◆ apply() [5/5]

void DirichletBC::apply ( GenericMatrix A,
GenericVector b,
const GenericVector x 
) const

Apply boundary condition to a linear system for a nonlinear problem

Parameters
[in,out]A(GenericMatrix) The matrix to apply boundary conditions to.
[in,out]b(GenericVector) The vector to apply boundary conditions to.
[in]x(GenericVector) Another vector (nonlinear problem).

◆ default_parameters()

static Parameters dolfin::DirichletBC::default_parameters ( )
inlinestatic

Default parameter values

Returns
Parameters

◆ function_space()

std::shared_ptr<const FunctionSpace> dolfin::DirichletBC::function_space ( ) const
inline

Return function space V

Returns
FunctionSpace The function space to which boundary conditions are applied.

◆ gather()

void DirichletBC::gather ( Map boundary_values) const

Get boundary values from neighbour processes. If a method other than "pointwise" is used, this is necessary to ensure all boundary dofs are marked on all processes.

Parameters
[in,out]boundary_values(Map&) Map from dof to boundary value.

◆ get_boundary_values()

void DirichletBC::get_boundary_values ( Map boundary_values) const

Get Dirichlet dofs and values. If a method other than 'pointwise' is used in parallel, the map may not be complete for local vertices since a vertex can have a bc applied, but the partition might not have a facet on the boundary. To ensure all local boundary dofs are marked, it is necessary to call gather() on the returned boundary values.

Parameters
[in,out]boundary_values(Map&) Map from dof to boundary value.

◆ markers()

const std::vector< std::size_t > & DirichletBC::markers ( ) const

Return boundary markers

Returns
std::vector<std::size_t>& Boundary markers (facets stored as pairs of cells and local facet numbers).

◆ method()

std::string DirichletBC::method ( ) const

Return method used for computing Dirichlet dofs

Returns
std::string Method used for computing Dirichlet dofs ("topological", "geometric" or "pointwise").

◆ operator=()

const DirichletBC & DirichletBC::operator= ( const DirichletBC bc)

Assignment operator. Either cached DOF data are assigned.

Parameters
[in]bc(DirichletBC) Another DirichletBC object.

◆ set_value()

void DirichletBC::set_value ( std::shared_ptr< const GenericFunction g)

Set value g for boundary condition, domain remains unchanged

Parameters
[in]g(GenericFucntion) The value.

◆ user_sub_domain()

std::shared_ptr< const SubDomain > DirichletBC::user_sub_domain ( ) const

Return shared pointer to subdomain

Returns
SubDomain Shared pointer to subdomain.

◆ value()

std::shared_ptr< const GenericFunction > DirichletBC::value ( ) const

Return boundary value g

Returns
GenericFunction The boundary values.

◆ zero()

void DirichletBC::zero ( GenericMatrix A) const

Make rows of matrix associated with boundary condition zero, useful for non-diagonal matrices in a block matrix.

Parameters
[in,out]A(GenericMatrix&) The matrix

◆ zero_columns()

void DirichletBC::zero_columns ( GenericMatrix A,
GenericVector b,
double  diag_val = 0 
) const

Make columns of matrix associated with boundary condition zero, and update a (right-hand side) vector to reflect the changes. Useful for non-diagonals.

Parameters
[in,out]A(GenericMatrix&) The matrix
[in,out]b(GenericVector&) The vector
[in]diag_val(double) This parameter would normally be -1, 0 or 1.

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