Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/v0.1.0/v0.9.0/cpp
DOLFINx  0.1.0
DOLFINx C++ interface
Public Member Functions | List of all members
dolfinx::fem::DirichletBC< T > Class Template Reference

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

#include <DirichletBC.h>

Public Member Functions

template<typename U >
 DirichletBC (const std::shared_ptr< const fem::Function< T >> &g, U &&dofs)
 Create a representation of a Dirichlet boundary condition where the space being constrained is the same as the function that defines the constraint values, i.e. share the same FunctionSpace. More...
 
 DirichletBC (const std::shared_ptr< const fem::Function< T >> &g, const std::array< std::vector< std::int32_t >, 2 > &V_g_dofs, std::shared_ptr< const fem::FunctionSpace > V)
 Create a representation of a Dirichlet boundary condition where the space being constrained and the function that defines the constraint values do not share the same FunctionSpace. A typical examples is when applying a constraint on a subspace. The (sub)space and the constrain function must have the same finite element. More...
 
 DirichletBC (const DirichletBC &bc)=default
 Copy constructor. More...
 
 DirichletBC (DirichletBC &&bc)=default
 Move constructor. More...
 
 ~DirichletBC ()=default
 Destructor.
 
DirichletBCoperator= (const DirichletBC &bc)=default
 Assignment operator. More...
 
DirichletBCoperator= (DirichletBC &&bc)=default
 Move assignment operator.
 
std::shared_ptr< const fem::FunctionSpacefunction_space () const
 The function space to which boundary conditions are applied. More...
 
std::shared_ptr< const fem::Function< T > > value () const
 Return boundary value function g. More...
 
std::pair< xtl::span< const std::int32_t >, std::int32_t > dof_indices () const
 Access dof indices (local indices, unrolled), including ghosts, to which a Dirichlet condition is applied, and the index to the first non-owned (ghost) index. The array of indices is sorted. More...
 
void set (xtl::span< T > x, double scale=1.0) const
 Set bc entries in x to scale * x_bc More...
 
void set (xtl::span< T > x, const xtl::span< const T > &x0, double scale=1.0) const
 Set bc entries in x to scale * (x0 - x_bc) More...
 
void dof_values (xtl::span< T > values) const
 
void mark_dofs (std::vector< bool > &markers) const
 Set markers[i] = true if dof i has a boundary condition applied. Value of markers[i] is not changed otherwise. More...
 

Detailed Description

template<typename T>
class dolfinx::fem::DirichletBC< T >

Interface for setting (strong) Dirichlet boundary conditions.

\(u = g \ \text{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 degrees of freedom to which the boundary condition applies.

Constructor & Destructor Documentation

◆ DirichletBC() [1/4]

template<typename T >
template<typename U >
dolfinx::fem::DirichletBC< T >::DirichletBC ( const std::shared_ptr< const fem::Function< T >> &  g,
U &&  dofs 
)
inline

Create a representation of a Dirichlet boundary condition where the space being constrained is the same as the function that defines the constraint values, i.e. share the same FunctionSpace.

Parameters
[in]gThe boundary condition value. The boundary condition can be applied to a function on the same space as g.
[in]dofsDegree-of-freedom block indices (std::vector<std::int32_t>) in the space of the boundary value function applied to V_dofs[i]. The dof block indices must be sorted.
Note
The indices in dofs are for blocks, e.g. a block index maps to 3 degrees-of-freedom if the dofmap associated with g has block size 3

◆ DirichletBC() [2/4]

template<typename T >
dolfinx::fem::DirichletBC< T >::DirichletBC ( const std::shared_ptr< const fem::Function< T >> &  g,
const std::array< std::vector< std::int32_t >, 2 > &  V_g_dofs,
std::shared_ptr< const fem::FunctionSpace V 
)
inline

Create a representation of a Dirichlet boundary condition where the space being constrained and the function that defines the constraint values do not share the same FunctionSpace. A typical examples is when applying a constraint on a subspace. The (sub)space and the constrain function must have the same finite element.

Parameters
[in]gThe boundary condition value
[in]V_g_dofsTwo arrays of degree-of-freedom indices. First array are indices in the space where boundary condition is applied (V), second array are indices in the space of the boundary condition value function g. The arrays must be sorted by the indices in the first array. The dof indices are unrolled, i.e. are not by dof block.
[in]VThe function (sub)space on which the boundary condition is applied
Note
The indices in dofs are unrolled and not for blocks

◆ DirichletBC() [3/4]

template<typename T >
dolfinx::fem::DirichletBC< T >::DirichletBC ( const DirichletBC< T > &  bc)
default

Copy constructor.

Parameters
[in]bcThe object to be copied

◆ DirichletBC() [4/4]

template<typename T >
dolfinx::fem::DirichletBC< T >::DirichletBC ( DirichletBC< T > &&  bc)
default

Move constructor.

Parameters
[in]bcThe object to be moved

Member Function Documentation

◆ dof_indices()

template<typename T >
std::pair<xtl::span<const std::int32_t>, std::int32_t> dolfinx::fem::DirichletBC< T >::dof_indices ( ) const
inline

Access dof indices (local indices, unrolled), including ghosts, to which a Dirichlet condition is applied, and the index to the first non-owned (ghost) index. The array of indices is sorted.

Returns
Sorted array of dof indices (unrolled) and index to the first entry in the dof index array that is not owned. Entries dofs[:pos] are owned and entries dofs[pos:] are ghosts.

◆ dof_values()

template<typename T >
void dolfinx::fem::DirichletBC< T >::dof_values ( xtl::span< T >  values) const
inline
Todo:
Review this function - it is almost identical to the 'DirichletBC::set' function

Set boundary condition value for entries with an applied boundary condition. Other entries are not modified.

Parameters
[in,out]valuesThe array in which to set the dof values. The array must be at least as long as the array associated with V1 (the space of the function that provides the dof values)

◆ function_space()

template<typename T >
std::shared_ptr<const fem::FunctionSpace> dolfinx::fem::DirichletBC< T >::function_space ( ) const
inline

The function space to which boundary conditions are applied.

Returns
The function space

◆ mark_dofs()

template<typename T >
void dolfinx::fem::DirichletBC< T >::mark_dofs ( std::vector< bool > &  markers) const
inline

Set markers[i] = true if dof i has a boundary condition applied. Value of markers[i] is not changed otherwise.

Parameters
[in,out]markersEntry makers[i] is set to true if dof i in V0 had a boundary condition applied, i.e. dofs which are fixed by a boundary condition. Other entries in markers are left unchanged.

◆ operator=()

template<typename T >
DirichletBC& dolfinx::fem::DirichletBC< T >::operator= ( const DirichletBC< T > &  bc)
default

Assignment operator.

Parameters
[in]bcAnother DirichletBC object

◆ set() [1/2]

template<typename T >
void dolfinx::fem::DirichletBC< T >::set ( xtl::span< T >  x,
const xtl::span< const T > &  x0,
double  scale = 1.0 
) const
inline

Set bc entries in x to scale * (x0 - x_bc)

Parameters
[in]xThe array in which to set scale * (x0 - x_bc)
[in]x0The array used in compute the value to set
[in]scaleThe scaling value to apply

◆ set() [2/2]

template<typename T >
void dolfinx::fem::DirichletBC< T >::set ( xtl::span< T >  x,
double  scale = 1.0 
) const
inline

Set bc entries in x to scale * x_bc

Parameters
[in]xThe array in which to set scale * x_bc[i], where x_bc[i] is the boundary value of x[i]. Entries in x that do not have a Dirichlet condition applied to them are unchanged. The length of x must be less than or equal to the index of the greatest boundary dof index. To set values only for degrees-of-freedom that are owned by the calling rank, the length of the array x should be equal to the number of dofs owned by this rank.
[in]scaleThe scaling value to apply

◆ value()

template<typename T >
std::shared_ptr<const fem::Function<T> > dolfinx::fem::DirichletBC< T >::value ( ) const
inline

Return boundary value function g.

Returns
The boundary values Function

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