Chaste Commit::f2ff7ee04e70ac9d06c57344df8d017dbb12b97b
BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM > Class Template Reference

#include <BoundaryConditionsContainer.hpp>

+ Inheritance diagram for BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >:
+ Collaboration diagram for BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >:

Public Types

typedef std::map< constBoundaryElement< ELEMENT_DIM-1, SPACE_DIM > *, constAbstractBoundaryCondition< SPACE_DIM > * >::const_iterator NeumannMapIterator
 
typedef AbstractBoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM > BaseClassType
 

Public Member Functions

 BoundaryConditionsContainer (bool deleteConditions=true)
 
 ~BoundaryConditionsContainer ()
 
void AddDirichletBoundaryCondition (const Node< SPACE_DIM > *pBoundaryNode, const AbstractBoundaryCondition< SPACE_DIM > *pBoundaryCondition, unsigned indexOfUnknown=0, bool checkIfBoundaryNode=true)
 
void AddNeumannBoundaryCondition (const BoundaryElement< ELEMENT_DIM-1, SPACE_DIM > *pBoundaryElement, const AbstractBoundaryCondition< SPACE_DIM > *pBoundaryCondition, unsigned indexOfUnknown=0)
 
void AddPeriodicBoundaryCondition (const Node< SPACE_DIM > *pNode1, const Node< SPACE_DIM > *pNode2)
 
void DefineZeroDirichletOnMeshBoundary (AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh, unsigned indexOfUnknown=0)
 
void DefineConstantDirichletOnMeshBoundary (AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh, double value, unsigned indexOfUnknown=0)
 
void DefineZeroNeumannOnMeshBoundary (AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh, unsigned indexOfUnknown=0)
 
void ApplyDirichletToLinearProblem (LinearSystem &rLinearSystem, bool applyToMatrix=true, bool applyToRhsVector=true)
 
void ApplyPeriodicBcsToLinearProblem (LinearSystem &rLinearSystem, bool applyToMatrix=true, bool applyToRhsVector=true)
 
void ApplyDirichletToNonlinearResidual (const Vec currentSolution, Vec residual, DistributedVectorFactory &rFactory)
 
void ApplyDirichletToNonlinearJacobian (Mat jacobian)
 
bool Validate (AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh)
 
double GetNeumannBCValue (const BoundaryElement< ELEMENT_DIM-1, SPACE_DIM > *pSurfaceElement, const ChastePoint< SPACE_DIM > &rX, unsigned indexOfUnknown=0)
 
bool HasNeumannBoundaryCondition (const BoundaryElement< ELEMENT_DIM-1, SPACE_DIM > *pSurfaceElement, unsigned indexOfUnknown=0)
 
bool AnyNonZeroNeumannConditions ()
 
NeumannMapIterator BeginNeumann ()
 
NeumannMapIterator EndNeumann ()
 
template<class Archive >
void LoadFromArchive (Archive &archive, AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh)
 
template<class Archive >
void MergeFromArchive (Archive &archive, AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh)
 
- Public Member Functions inherited from AbstractBoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >
 AbstractBoundaryConditionsContainer (bool deleteConditions=true)
 
 ~AbstractBoundaryConditionsContainer ()
 
bool HasDirichletBoundaryConditions ()
 
double GetDirichletBCValue (const Node< SPACE_DIM > *pBoundaryNode, unsigned indexOfUnknown=0)
 
bool HasDirichletBoundaryCondition (const Node< SPACE_DIM > *pNode, unsigned indexOfUnknown=0)
 
void ResetDirichletCommunication ()
 

Private Member Functions

template<class Archive >
void save (Archive &archive, const unsigned int version) const
 
template<class Archive >
void load (Archive &archive, const unsigned int version)
 

Private Attributes

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

Friends

class boost::serialization::access
 

Additional Inherited Members

- Protected Types inherited from AbstractBoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >
typedef std::map< const Node< SPACE_DIM > *, const AbstractBoundaryCondition< SPACE_DIM > *, LessThanNode< SPACE_DIM > > DirichletMapType
 
typedef std::map< constNode< SPACE_DIM > *, constAbstractBoundaryCondition< SPACE_DIM > *, LessThanNode< SPACE_DIM > >::const_iterator DirichletIteratorType
 
- Protected Member Functions inherited from AbstractBoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >
void DeleteDirichletBoundaryConditions (std::set< const AbstractBoundaryCondition< SPACE_DIM > * > alreadyDeletedConditions=std::set< const AbstractBoundaryCondition< SPACE_DIM > * >())
 
- Protected Attributes inherited from AbstractBoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >
DirichletMapTypempDirichletMap [PROBLEM_DIM]
 
DirichletIteratorType mDirichIterator
 
bool mHasDirichletBCs
 
bool mCheckedAndCommunicatedIfDirichletBcs
 
bool mDeleteConditions
 

Detailed Description

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
class BoundaryConditionsContainer< ELEMENT_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:
#1321 Various operations are currently very inefficient - there is certainly scope for optimisation here!

Definition at line 67 of file BoundaryConditionsContainer.hpp.

Member Typedef Documentation

◆ BaseClassType

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
typedef AbstractBoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM> BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::BaseClassType

Base class type.

Definition at line 76 of file BoundaryConditionsContainer.hpp.

◆ NeumannMapIterator

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
typedef std::map<constBoundaryElement<ELEMENT_DIM-1,SPACE_DIM>*,constAbstractBoundaryCondition<SPACE_DIM>*>::const_iterator BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::NeumannMapIterator

Type of a read-only iterator over Neumann boundary conditions.

Definition at line 73 of file BoundaryConditionsContainer.hpp.

Constructor & Destructor Documentation

◆ BoundaryConditionsContainer()

◆ ~BoundaryConditionsContainer()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::~BoundaryConditionsContainer ( )

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 69 of file BoundaryConditionsContainerImplementation.hpp.

Member Function Documentation

◆ AddDirichletBoundaryCondition()

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

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
pBoundaryNodePointer to a node on the boundary.
pBoundaryConditionPointer to the Dirichlet boundary condition at that node.
indexOfUnknowndefaults to 0
checkIfBoundaryNodedefaults to true

Definition at line 106 of file BoundaryConditionsContainerImplementation.hpp.

◆ AddNeumannBoundaryCondition()

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

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
pBoundaryElementPointer to an element on the boundary
pBoundaryConditionPointer to the Neumann boundary condition on that element
indexOfUnknowndefaults to 0

Definition at line 137 of file BoundaryConditionsContainerImplementation.hpp.

References AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNode(), and ConstBoundaryCondition< SPACE_DIM >::GetValue().

◆ AddPeriodicBoundaryCondition()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::AddPeriodicBoundaryCondition ( const Node< SPACE_DIM > *  pNode1,
const Node< SPACE_DIM > *  pNode2 
)

Add a periodic boundary condition: provide two nodes to be identified when solving

Parameters
pNode1node 1
pNode2node 2

This method identifies the nodes for all unknowns, so doesn't have to be called for each unknown.

Definition at line 121 of file BoundaryConditionsContainerImplementation.hpp.

References Node< SPACE_DIM >::IsBoundaryNode().

◆ AnyNonZeroNeumannConditions()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
bool BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::AnyNonZeroNeumannConditions ( )
Returns
whether there are any non-zero Neumann boundary conditions

Definition at line 613 of file BoundaryConditionsContainerImplementation.hpp.

◆ ApplyDirichletToLinearProblem()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyDirichletToLinearProblem ( LinearSystem rLinearSystem,
bool  applyToMatrix = true,
bool  applyToRhsVector = true 
)

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
rLinearSystemLinear system on which to apply boundary conditions
applyToMatrixThis 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.
applyToRhsVectorSimilarly, whether to apply the changes to the RHS vector (b in Ax=b).

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 value. This last term is then 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 257 of file BoundaryConditionsContainerImplementation.hpp.

References GenericEventHandler< 16, HeartEventHandler >::BeginEvent().

◆ ApplyDirichletToNonlinearJacobian()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyDirichletToNonlinearJacobian ( Mat  jacobian)

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)

Parameters
jacobian

Definition at line 522 of file BoundaryConditionsContainerImplementation.hpp.

References PetscMatTools::Finalise(), and PetscMatTools::ZeroRowsWithValueOnDiagonal().

◆ ApplyDirichletToNonlinearResidual()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyDirichletToNonlinearResidual ( const Vec  currentSolution,
Vec  residual,
DistributedVectorFactory rFactory 
)

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)

Parameters
currentSolution
residual
rFactorythe factory to use to create DistributedVector objects

Definition at line 488 of file BoundaryConditionsContainerImplementation.hpp.

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

◆ ApplyPeriodicBcsToLinearProblem()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyPeriodicBcsToLinearProblem ( LinearSystem rLinearSystem,
bool  applyToMatrix = true,
bool  applyToRhsVector = true 
)

Alter the given linear system to satisfy periodic boundary conditions.

For one of the two nodes that have been identified, the row corresponding to the unknown which has periodic BCs, for one of the nodes, ie replaced with [0 0 0 ... -1 0 0 .. 0 1 0 .. 0] where the 1 is the diagonal entry and the -1 on the column corresponding to that unknown and the other node.

The entry in the RHS vector is zeroed.

Parameters
rLinearSystemLinear system on which to apply boundary conditions
applyToMatrixThis 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.
applyToRhsVectorSimilarly, whether to apply the changes to the RHS vector (b in Ax=b).

Definition at line 422 of file BoundaryConditionsContainerImplementation.hpp.

References EXCEPT_IF_NOT, LinearSystem::rGetLhsMatrix(), PetscMatTools::SetElement(), LinearSystem::SetRhsVectorElement(), and LinearSystem::ZeroMatrixRowsWithValueOnDiagonal().

◆ BeginNeumann()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::NeumannMapIterator BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::BeginNeumann ( )
Returns
iterator pointing to the first Neumann boundary condition

Definition at line 627 of file BoundaryConditionsContainerImplementation.hpp.

◆ DefineConstantDirichletOnMeshBoundary()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::DefineConstantDirichletOnMeshBoundary ( AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *  pMesh,
double  value,
unsigned  indexOfUnknown = 0 
)

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

Parameters
pMeshPointer to a mesh object, from which we extract the boundary
valuethe value of the constant Dirichlet boundary condition
indexOfUnknowndefaults to 0

Definition at line 186 of file BoundaryConditionsContainerImplementation.hpp.

◆ DefineZeroDirichletOnMeshBoundary()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::DefineZeroDirichletOnMeshBoundary ( AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *  pMesh,
unsigned  indexOfUnknown = 0 
)

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

Parameters
pMeshPointer to a mesh object, from which we extract the boundary
indexOfUnknowndefaults to 0

Definition at line 179 of file BoundaryConditionsContainerImplementation.hpp.

◆ DefineZeroNeumannOnMeshBoundary()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::DefineZeroNeumannOnMeshBoundary ( AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *  pMesh,
unsigned  indexOfUnknown = 0 
)

This function defines zero Neumann boundary conditions on every boundary element of the mesh.

Parameters
pMeshPointer to a mesh object, from which we extract the boundary
indexOfUnknowndefaults to 0

Definition at line 206 of file BoundaryConditionsContainerImplementation.hpp.

References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNumBoundaryElements().

◆ EndNeumann()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::NeumannMapIterator BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::EndNeumann ( )
Returns
iterator pointing to one past the last Neumann boundary condition

Definition at line 634 of file BoundaryConditionsContainerImplementation.hpp.

Referenced by AbstractFeSurfaceIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::DoAssemble().

◆ GetNeumannBCValue()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
double BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::GetNeumannBCValue ( const BoundaryElement< ELEMENT_DIM-1, SPACE_DIM > *  pSurfaceElement,
const ChastePoint< SPACE_DIM > &  rX,
unsigned  indexOfUnknown = 0 
)
Returns
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.

Parameters
pSurfaceElementpointer to a boundary element
rXa point
indexOfUnknowndefaults to 0

Definition at line 578 of file BoundaryConditionsContainerImplementation.hpp.

◆ HasNeumannBoundaryCondition()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
bool BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::HasNeumannBoundaryCondition ( const BoundaryElement< ELEMENT_DIM-1, SPACE_DIM > *  pSurfaceElement,
unsigned  indexOfUnknown = 0 
)
Returns
true if there is a Neumann boundary condition defined on the given element.
Todo:
#1321 This is a horrendously inefficient fix. Perhaps have flag in element object?
Parameters
pSurfaceElementpointer to a boundary element
indexOfUnknowndefaults to 0

Definition at line 602 of file BoundaryConditionsContainerImplementation.hpp.

◆ load()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
template<class Archive >
void BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::load ( Archive &  archive,
const unsigned int  version 
)
inlineprivate

Load this container, but not its content.

Objects loading a boundary conditions container should call LoadFromArchive on the new object immediately after loading it from the archive.

Note that boundary conditions should be saved to the ProcessSpecificArchive, since if a DistributedTetrahedralMesh is used each process will only know a portion of the mesh, and hence a portion of the boundary conditions.

Extra care needs to be taken when migrating to ensure that boundary conditions are loaded appropriately. See BidomainProblem::LoadExtraArchiveForBidomain and AbstractCardiacProblem::LoadExtraArchive for examples.

Parameters
archive
version

Definition at line 389 of file BoundaryConditionsContainer.hpp.

◆ LoadFromArchive()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
template<class Archive >
void BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::LoadFromArchive ( Archive &  archive,
AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *  pMesh 
)
inline

Load a collection of boundary conditions from an archive.

Note
We assume this collection is empty prior to being called. If it is not, any boundary conditions already present may get replaced by conditions loaded from the archive, which may lead to a memory leak.

This method only loads data if mLoadedFromArchive is false, to allow for multiple pointers to the same container to be handled correctly. It sets mLoadedFromArchive when done.

Parameters
archivethe archive to load from
pMeshthe mesh to use to resolve Node and BoundaryElement indices

Definition at line 335 of file BoundaryConditionsContainer.hpp.

References BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::MergeFromArchive(), and BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mLoadedFromArchive.

◆ MergeFromArchive()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
template<class Archive >
void BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::MergeFromArchive ( Archive &  archive,
AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *  pMesh 
)

Load extra boundary conditions from an archive to add to this collection.

Multiple pointers to the same container need to be handled by the caller - we assume there will be conditions to load. Sets mLoadedFromArchive when done.

Parameters
archivethe archive to load from
pMeshthe mesh to use to resolve Node and BoundaryElement indices

Definition at line 437 of file BoundaryConditionsContainer.hpp.

References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetBoundaryElement(), and AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNodeFromPrePermutationIndex().

Referenced by BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::LoadFromArchive().

◆ save()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
template<class Archive >
void BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::save ( Archive &  archive,
const unsigned int  version 
) const
private

Save this container and its contents.

Parameters
archive
version

Definition at line 398 of file BoundaryConditionsContainer.hpp.

◆ Validate()

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
bool BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::Validate ( AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *  pMesh)

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.

Parameters
pMeshPointer to the mesh to check for validity.
Returns
true iff all boundaries have boundary conditions defined.

Definition at line 549 of file BoundaryConditionsContainerImplementation.hpp.

References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetBoundaryElementIteratorBegin(), and AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetBoundaryElementIteratorEnd().

Friends And Related Symbol Documentation

◆ boost::serialization::access

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
friend class boost::serialization::access
friend

Needed for serialization.

Definition at line 360 of file BoundaryConditionsContainer.hpp.

Member Data Documentation

◆ mAnyNonZeroNeumannConditionsForUnknown

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
bool BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mAnyNonZeroNeumannConditionsForUnknown[PROBLEM_DIM]
private

Array storing whether there are any Neumann boundary conditions for each unknown.

Definition at line 94 of file BoundaryConditionsContainer.hpp.

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

◆ mLastNeumannCondition

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
NeumannMapIterator BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mLastNeumannCondition[PROBLEM_DIM]
private

◆ mLoadedFromArchive

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
bool BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mLoadedFromArchive
private

◆ mpNeumannMap

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

List (map) of Neumann boundary conditions.

Definition at line 81 of file BoundaryConditionsContainer.hpp.

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

◆ mpPeriodicBcMap

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
std::map< const Node<SPACE_DIM> *, const Node<SPACE_DIM> * >* BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpPeriodicBcMap[PROBLEM_DIM]
private

Nodes to be identified when periodic boundary conditions are applied

Definition at line 84 of file BoundaryConditionsContainer.hpp.

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

◆ mpZeroBoundaryCondition

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
ConstBoundaryCondition<SPACE_DIM>* BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpZeroBoundaryCondition
private

A zero boundary condition, used for other unknowns in ApplyNeumannBoundaryCondition

Definition at line 97 of file BoundaryConditionsContainer.hpp.

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


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