DirichletBC¶
-
class
dolfin.cpp.fem.
DirichletBC
(*args)¶ Bases:
dolfin.cpp.fem.HierarchicalDirichletBC
,dolfin.cpp.common.Variable
This class specifies the interface for setting (strong) Dirichlet boundary conditions for partial differential equations,
\[u = g \hbox{ 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.
Overloaded versions
DirichletBC(V, g, sub_domain, method=”topological”, check_midpoint=true)
Create boundary condition for subdomain
- Arguments
- V (
FunctionSpace
) The function space.
- g (
GenericFunction
) The value.
- sub_domain (
SubDomain
) The subdomain.
- method (str)
Optional argument: A string specifying the method to identify dofs.
- V (
DirichletBC(V, g, sub_domain, method=”topological”, check_midpoint=true)
Create boundary condition for subdomain
- Arguments
- V (
FunctionSpace
) The function space
- g (
GenericFunction
) The value
- sub_domain (
SubDomain
) The subdomain
- method (str)
Optional argument: A string specifying the method to identify dofs
- V (
DirichletBC(V, g, sub_domains, sub_domain, method=”topological”)
Create boundary condition for subdomain specified by index
- Arguments
- V (
FunctionSpace
) The function space.
- g (
GenericFunction
) The value.
- sub_domains (
MeshFunction
) Subdomain markers
- sub_domain (int)
The subdomain index (number)
- method (str)
Optional argument: A string specifying the method to identify dofs.
- V (
DirichletBC(V, g, sub_domains, sub_domain, method=”topological”)
Create boundary condition for subdomain specified by index
- Arguments
- V (
FunctionSpace
) The function space.
- g (
GenericFunction
) The value.
- sub_domains (
MeshFunction
) Subdomain markers
- sub_domain (int)
The subdomain index (number)
- method (str)
Optional argument: A string specifying the method to identify dofs.
- V (
DirichletBC(V, g, sub_domain, method=”topological”)
Create boundary condition for boundary data included in the mesh
- Arguments
- V (
FunctionSpace
) The function space.
- g (
GenericFunction
) The value.
- sub_domain (int)
The subdomain index (number)
- method (str)
Optional argument: A string specifying the method to identify dofs.
- V (
DirichletBC(V, g, sub_domain, method=”topological”)
Create boundary condition for boundary data included in the mesh
- Arguments
- V (
FunctionSpace
) The function space.
- g (
GenericFunction
) The value.
- sub_domain (int)
The subdomain index (number)
- method (str)
Optional argument: A string specifying the method to identify dofs.
- V (
DirichletBC(V, g, markers, method=”topological”)
Create boundary condition for subdomain by boundary markers (cells, local facet numbers)
- Arguments
- V (
FunctionSpace
) The function space.
- g (
GenericFunction
) The value.
- markers (numpy.array(int))
Subdomain markers (facet index local to process)
- method (str)
Optional argument: A string specifying the method to identify dofs.
- V (
DirichletBC(bc)
Copy constructor. Either cached DOF data are copied.
- Arguments
- bc (
DirichletBC
) The object to be copied.
- bc (
-
apply
()¶ Overloaded versions
apply(A)
Apply boundary condition to a matrix
- Arguments
- A (
GenericMatrix
) The matrix to apply boundary condition to.
- A (
apply(b)
Apply boundary condition to a vector
- Arguments
- b (
GenericVector
) The vector to apply boundary condition to.
- b (
apply(A, b)
Apply boundary condition to a linear system
- Arguments
- A (
GenericMatrix
) The matrix to apply boundary condition to.
- b (
GenericVector
) The vector to apply boundary condition to.
- A (
apply(b, x)
Apply boundary condition to vectors for a nonlinear problem
- Arguments
- b (
GenericVector
) The vector to apply boundary conditions to.
- x (
GenericVector
) Another vector (nonlinear problem).
- b (
apply(A, b, x)
Apply boundary condition to a linear system for a nonlinear problem
- Arguments
- A (
GenericMatrix
) The matrix to apply boundary conditions to.
- b (
GenericVector
) The vector to apply boundary conditions to.
- x (
GenericVector
) Another vector (nonlinear problem).
- A (
-
static
default_parameters
()¶ Default parameter values
-
function_space
()¶ Return the FunctionSpace
-
get_boundary_values
()¶ 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.
- Arguments
- boundary_values (std::unordered_map<std::size_t, double>)
- Map from dof to boundary value.
- method (str)
- Optional argument: A string specifying which method to use.
-
homogenize
()¶ Set value to 0.0
-
is_compatible
()¶ Check if given function is compatible with boundary condition (checking only vertex values)
- Arguments
- v (
GenericFunction
) - The function to check for compatibility with boundary condition.
- v (
- Returns
- bool
- True if compatible.
-
markers
()¶ Return boundary markers
- Returns
- numpy.array((int, int))
- Boundary markers (facets stored as pairs of cells and local facet numbers).
-
method
()¶ Return method used for computing Dirichlet dofs
- Returns
- str
- Method used for computing Dirichlet dofs (“topological”, “geometric” or “pointwise”).
-
set_value
()¶ Overloaded versions
set_value(g)
Set value g for boundary condition, domain remains unchanged
- Arguments
- g (
GenericFunction
) The value.
- g (
set_value(g)
Set value g for boundary condition, domain remains unchanged
- Arguments
- g (
GenericFunction
) The value.
- g (
-
thisown
¶ The membership flag
-
user_sub_domain
()¶ Return shared pointer to subdomain
- Returns
SubDomain
- Shared pointer to subdomain.
-
value
()¶ Return boundary value g
- Returns
GenericFunction
- The boundary values.
-
zero
()¶ Make rows of matrix associated with boundary condition zero, useful for non-diagonal matrices in a block matrix.
- Arguments
- A (
GenericMatrix
) - The matrix
- A (
-
zero_columns
()¶ Make columns of matrix associated with boundary condition zero, and update a (right-hand side) vector to reflect the changes. Useful for non-diagonals.
- Arguments
- A (
GenericMatrix
) - The matrix
- b (
GenericVector
) - The vector
- diag_val (float)
- This parameter would normally be -1, 0 or 1.
- A (