BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, PROBLEM_DIM > Class Template Reference

#include <BoundaryConditionsContainer.hpp>

Inheritance diagram for BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >:

Inheritance graph
[legend]
Collaboration diagram for BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >:

Collaboration graph
[legend]

List of all members.

Public Types

typedef std::map< const
BoundaryElement< ELEM_DIM-1,
SPACE_DIM > *, const
AbstractBoundaryCondition
< SPACE_DIM >
* >::const_iterator 
NeumannMapIterator
 Type of a read-only iterator over Neumann conditions.

Public Member Functions

 BoundaryConditionsContainer ()
 ~BoundaryConditionsContainer ()
void AddDirichletBoundaryCondition (const Node< SPACE_DIM > *pBoundaryNode, const AbstractBoundaryCondition< SPACE_DIM > *pBoundaryCondition, unsigned indexOfUnknown=0, bool checkIfBoundaryNode=true)
void AddNeumannBoundaryCondition (const BoundaryElement< ELEM_DIM-1, SPACE_DIM > *pBoundaryElement, const AbstractBoundaryCondition< SPACE_DIM > *pBoundaryCondition, unsigned indexOfUnknown=0)
void DefineZeroDirichletOnMeshBoundary (AbstractMesh< ELEM_DIM, SPACE_DIM > *pMesh, unsigned indexOfUnknown=0)
void DefineConstantDirichletOnMeshBoundary (AbstractMesh< ELEM_DIM, SPACE_DIM > *pMesh, double value, unsigned indexOfUnknown=0)
void DefineZeroNeumannOnMeshBoundary (AbstractMesh< ELEM_DIM, SPACE_DIM > *pMesh, unsigned indexOfUnknown=0)
void ApplyDirichletToLinearProblem (LinearSystem &rLinearSystem, bool applyToMatrix=true)
void ApplyDirichletToNonlinearResidual (const Vec currentSolution, Vec residual)
void ApplyDirichletToNonlinearJacobian (Mat jacobian)
bool Validate (AbstractMesh< ELEM_DIM, SPACE_DIM > *pMesh)
double GetNeumannBCValue (const BoundaryElement< ELEM_DIM-1, SPACE_DIM > *pSurfaceElement, const ChastePoint< SPACE_DIM > &x, unsigned indexOfUnknown=0)
bool HasNeumannBoundaryCondition (const BoundaryElement< ELEM_DIM-1, SPACE_DIM > *pSurfaceElement, unsigned indexOfUnknown=0)
bool AnyNonZeroNeumannConditions ()
NeumannMapIterator BeginNeumann ()
NeumannMapIterator EndNeumann ()

Private Attributes

std::map< const
BoundaryElement< ELEM_DIM-1,
SPACE_DIM > *, const
AbstractBoundaryCondition
< SPACE_DIM > * > * 
mpNeumannMap [PROBLEM_DIM]
NeumannMapIterator mLastNeumannCondition [PROBLEM_DIM]
bool mAnyNonZeroNeumannConditionsForUnknown [PROBLEM_DIM]
ConstBoundaryCondition
< SPACE_DIM > * 
mpZeroBoundaryCondition


Detailed Description

template<unsigned ELEM_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
class BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >

Boundary Conditions Container

This class contains a list of nodes on the dirichlet boundary and associated dirichlet boundary conditions, and a list of surface elements on the neumann boundary and associated neumann boundary conditions.

Todo:
Various operations are currently very inefficient - there is certainly scope for optimisation here!

Definition at line 58 of file BoundaryConditionsContainer.hpp.


Constructor & Destructor Documentation

template<unsigned ELEM_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::BoundaryConditionsContainer (  )  [inline]

Constructor calls base constuctor and allocates memory for the neumann boundary conditions lists.

Definition at line 40 of file BoundaryConditionsContainerImplementation.hpp.

References BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::mpNeumannMap, and BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::mpZeroBoundaryCondition.

template<unsigned ELEM_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::~BoundaryConditionsContainer (  )  [inline]

Note that the destructor will delete memory for each boundary condition object, as well as for the internal bookkeeping of this class.

Definition at line 56 of file BoundaryConditionsContainerImplementation.hpp.

References BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::mpNeumannMap, and BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::mpZeroBoundaryCondition.


Member Function Documentation

template<unsigned ELEM_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::AddDirichletBoundaryCondition ( const Node< SPACE_DIM > *  pBoundaryNode,
const AbstractBoundaryCondition< SPACE_DIM > *  pBoundaryCondition,
unsigned  indexOfUnknown = 0,
bool  checkIfBoundaryNode = true 
) [inline]

Add a dirichlet boundary condition specifying two parameters, a pointer to a node, and a pointer to a boundary condition object associated with that node.

The destructor for the BoundaryConditionsContainer will destroy the boundary conditions objects.

Parameters:
pBoundaryNode Pointer to a node on the boundary.
pBoundaryCondition Pointer to the dirichlet boundary condition at that node.

Definition at line 88 of file BoundaryConditionsContainerImplementation.hpp.

References Node< SPACE_DIM >::IsBoundaryNode().

Referenced by BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::DefineConstantDirichletOnMeshBoundary(), Electrodes< DIM >::Electrodes(), and BidomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::SetFixedExtracellularPotentialNodes().

template<unsigned ELEM_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::AddNeumannBoundaryCondition ( const BoundaryElement< ELEM_DIM-1, SPACE_DIM > *  pBoundaryElement,
const AbstractBoundaryCondition< SPACE_DIM > *  pBoundaryCondition,
unsigned  indexOfUnknown = 0 
) [inline]

Add a neumann boundary condition specifying two parameters, a pointer to a surface element, and a pointer to a boundary condition object associated with that element.

The destructor for the BoundaryConditionsContainer will destroy the boundary conditions objects.

Note that the value of a Neumann boundary condition should specify D * grad(u).n, not just grad(u).n.

Take care if using non-zero neumann boundary conditions in 1d. If applied at the left hand end you need to multiply the value by -1 to get the right answer.

Parameters:
pBoundaryElement Pointer to an element on the boundary.
pBoundaryCondition Pointer to the neumann boundary condition on that element.

Definition at line 103 of file BoundaryConditionsContainerImplementation.hpp.

References BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::mpNeumannMap, and BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::mpZeroBoundaryCondition.

Referenced by AbstractConvergenceTester< CELL, CARDIAC_PROBLEM, DIM, PROBLEM_DIM >::Converge(), BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::DefineZeroNeumannOnMeshBoundary(), and Electrodes< DIM >::Electrodes().

template<unsigned ELEM_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::DefineZeroDirichletOnMeshBoundary ( AbstractMesh< ELEM_DIM, SPACE_DIM > *  pMesh,
unsigned  indexOfUnknown = 0 
) [inline]

This function defines zero dirichlet boundary conditions on every boundary node of the mesh.

Parameters:
pMesh Pointer to a mesh object, from which we extract the boundary.

Definition at line 131 of file BoundaryConditionsContainerImplementation.hpp.

References BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::DefineConstantDirichletOnMeshBoundary().

template<unsigned ELEM_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::DefineConstantDirichletOnMeshBoundary ( AbstractMesh< ELEM_DIM, SPACE_DIM > *  pMesh,
double  value,
unsigned  indexOfUnknown = 0 
) [inline]

This function defines constant dirichlet boundary conditions on every boundary node of the mesh.

Parameters:
pMesh Pointer to a mesh object, from which we extract the boundary.
value the value of the constant Dirichlet boundary condition

Definition at line 138 of file BoundaryConditionsContainerImplementation.hpp.

References BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::AddDirichletBoundaryCondition().

Referenced by BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::DefineZeroDirichletOnMeshBoundary().

template<unsigned ELEM_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::DefineZeroNeumannOnMeshBoundary ( AbstractMesh< ELEM_DIM, SPACE_DIM > *  pMesh,
unsigned  indexOfUnknown = 0 
) [inline]

template<unsigned ELEM_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyDirichletToLinearProblem ( LinearSystem rLinearSystem,
bool  applyToMatrix = true 
) [inline]

Alter the given linear system to satisfy dirichlet boundary conditions

If the number of unknowns is greater than one, it is assumed the solution vector is of the form (in the case of two unknowns u and v, and N nodes): solnvec = (U_1, V_1, U_2, V_2, ...., U_N, V_N)

Parameters:
rLinearSystem Linear system on which to apply boundary conditions
applyToMatrix This optional parameter can be set as false to ensure that the matrix of the linear system is not updated. To be used when the matrix does not change between time steps.
Modifies a linear system to incorporate Dirichlet boundary conditions

The BCs are imposed in such a way as to ensure that a symmetric linear system remains symmetric. For each node with a boundary condition applied, both the corresponding row and column are zero'd and the RHS vector modified to take into account the zero'd column. See #577.

Suppose we have a matrix [a b c] [x] = [ b1 ] [d e f] [y] [ b2 ] [g h i] [z] [ b3 ] and we want to apply the boundary condition x=v without losing symmetry if the matrix is symmetric. We apply the boundary condition [1 0 0] [x] = [ v ] [d e f] [y] [ b2 ] [g h i] [z] [ b3 ] and then zero the column as well, adding a term to the RHS to take account for the zero-matrix components [1 0 0] [x] = [ v ] - v[ 0 ] [0 e f] [y] [ b2 ] [ d ] [0 h i] [z] [ b3 ] [ g ] Note the last term is the first column of the matrix, with one component zeroed, and multiplied by the boundary condition. This last term is the stored in rLinearSystem.rGetDirichletBoundaryConditionsVector(), and in general form is the SUM_{d=1..D} v_d a'_d where v_d is the boundary value of boundary condition d (d an index into the matrix), and a'_d is the dth-column of the matrix but with the d-th component zeroed, and where there are D boundary conditions

Definition at line 209 of file BoundaryConditionsContainerImplementation.hpp.

References LinearSystem::AssembleFinalLinearSystem(), LinearSystem::rGetDirichletBoundaryConditionsVector(), LinearSystem::rGetLhsMatrix(), LinearSystem::rGetRhsVector(), LinearSystem::SetMatrixElement(), LinearSystem::SetRhsVectorElement(), LinearSystem::ZeroMatrixColumn(), and LinearSystem::ZeroMatrixRow().

Referenced by AbstractLinearAssembler< DIM, DIM, 1, false, MonodomainRhsMatrixAssembler< DIM > >::ApplyDirichletConditions().

template<unsigned ELEM_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyDirichletToNonlinearResidual ( const Vec  currentSolution,
Vec  residual 
) [inline]

Alter the residual vector for a nonlinear system to satisfy dirichlet boundary conditions.

If the number of unknowns is greater than one, it is assumed the solution vector is of the form (in the case of two unknowns u and v, and N nodes): solnvec = (U_1, V_1, U_2, V_2, ...., U_N, V_N)

Definition at line 296 of file BoundaryConditionsContainerImplementation.hpp.

References DistributedVector::IsGlobalIndexLocal(), and DistributedVector::Restore().

Referenced by AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, 1, SimpleNonlinearEllipticAssembler< ELEMENT_DIM, SPACE_DIM > >::ApplyDirichletConditions().

template<unsigned ELEM_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyDirichletToNonlinearJacobian ( Mat  jacobian  )  [inline]

Alter the jacobian matrix vector for a nonlinear system to satisfy dirichlet boundary conditions.

If the number of unknowns is greater than one, it is assumed the solution vector is of the form (in the case of two unknowns u and v, and N nodes): solnvec = (U_1, V_1, U_2, V_2, ...., U_N, V_N)

Definition at line 327 of file BoundaryConditionsContainerImplementation.hpp.

Referenced by AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, 1, SimpleNonlinearEllipticAssembler< ELEMENT_DIM, SPACE_DIM > >::ApplyDirichletConditions().

template<unsigned ELEM_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
bool BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::Validate ( AbstractMesh< ELEM_DIM, SPACE_DIM > *  pMesh  )  [inline]

Check that we have boundary conditions defined everywhere on mesh boundary.

We iterate over all surface elements, and check either that they have an associated Neumann condition, or that each node in the element has an associated Dirichlet condition.

Todo:
Might we want to throw an exception specifying which node failed? What about checking for multiple conditions at a point (might be intentional)?
Parameters:
pMesh Pointer to the mesh to check for validity.
Returns:
true iff all boundaries have boundary conditions defined.

Definition at line 355 of file BoundaryConditionsContainerImplementation.hpp.

References BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::HasNeumannBoundaryCondition().

template<unsigned ELEM_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
double BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::GetNeumannBCValue ( const BoundaryElement< ELEM_DIM-1, SPACE_DIM > *  pSurfaceElement,
const ChastePoint< SPACE_DIM > &  x,
unsigned  indexOfUnknown = 0 
) [inline]

Obtain value of neumann boundary condition at a specified point in a given surface element

It is up to the user to ensure that the point x is contained in the surface element.

Definition at line 384 of file BoundaryConditionsContainerImplementation.hpp.

References BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::mpNeumannMap.

Referenced by SimpleNonlinearEllipticAssembler< ELEMENT_DIM, SPACE_DIM >::ComputeVectorSurfaceTerm(), and BidomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::ComputeVectorSurfaceTerm().

template<unsigned ELEM_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
bool BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::HasNeumannBoundaryCondition ( const BoundaryElement< ELEM_DIM-1, SPACE_DIM > *  pSurfaceElement,
unsigned  indexOfUnknown = 0 
) [inline]

Test if there is a Neumann boundary condition defined on the given element. Used by SimpleLinearEllipticAssembler.

Todo:
This is a horrendously inefficient fix. Perhaps have flag in element object?

Definition at line 408 of file BoundaryConditionsContainerImplementation.hpp.

References BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::mpNeumannMap.

Referenced by BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::Validate().


Member Data Documentation

template<unsigned ELEM_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
std::map< const BoundaryElement<ELEM_DIM-1, SPACE_DIM> *, const AbstractBoundaryCondition<SPACE_DIM>* >* BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::mpNeumannMap[PROBLEM_DIM] [private]

template<unsigned ELEM_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
ConstBoundaryCondition<SPACE_DIM>* BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::mpZeroBoundaryCondition [private]


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

Generated on Wed Mar 18 12:52:23 2009 for Chaste by  doxygen 1.5.5