LifeV
BCManageNormal< MatrixType > Class Template Reference

BCManageNormal - class for handling normal essential boundary conditions. More...

#include <BCManageNormal.hpp>

+ Collaboration diagram for BCManageNormal< MatrixType >:

Private Attributes

bool M_dataBuilt
 true when there are stored normals More...
 
matrixPtr_Type M_rotationMatrixPtr
 Shared pointer to the rotation matrix and its transpose. More...
 
matrixPtr_Type M_rotationMatrixTransposePtr
 
epetraMapPtr_Type M_localMapEpetraPtr
 Shared pointer to the local Map. More...
 
epetraVectorPtr_type M_firstTangentPtr
 Shared pointer to the vector of the first tangents to the domain boundary. More...
 
epetraVectorPtr_type M_secondTangentPtr
 Shared pointer to the vector of the second tangents to the domain boundary. More...
 
epetraVectorPtr_type M_normalPtr
 Shared pointer to the vector of the normals to the domain boundary. More...
 
epetraVectorPtr_type M_coordPtr
 Shared coordinates of the point. More...
 
UInt M_numDof
 Number of degree of freedom (of one component) More...
 
UInt M_numInvoledDof
 Number of degree of freedom (of one component) involved in boundary conditions. More...
 
flagsMap_Type M_flags
 flag of the boundary elements More...
 
versorsMap_Type M_givenVersors
 Versors that are given by the user. More...
 

Public Types

typedef std::map< ID, IDflagsMap_Type
 
typedef std::map< ID, VectorversorsMap_Type
 
typedef MatrixType matrix_Type
 
typedef std::shared_ptr< matrix_TypematrixPtr_Type
 
typedef std::shared_ptr< MapEpetraepetraMapPtr_Type
 
typedef std::shared_ptr< VectorEpetraepetraVectorPtr_type
 

Constructors & Destructor

 BCManageNormal ()
 Empty Constructor. More...
 
 BCManageNormal (BCManageNormal const &bcManageNormal)
 Copy constructor. More...
 
 ~BCManageNormal ()
 Destructor. More...
 

Operators

BCManageNormaloperator= (const BCManageNormal &bcManageNormal)
 Assignment operator. More...
 

Methods

void init (const BCBase &boundaryCondition, const Real &time)
 Store boundary DOFs id and coordinates related to the boundary conditions. More...
 
template<typename MeshType >
void build (const MeshType &mesh, const DOF &dof, CurrentFEManifold &currentBdFE, const MapEpetra &map, UInt offset)
 Build the rotation matrix. More...
 
template<typename VectorType >
void bcShiftToNormalTangentialCoordSystem (matrix_Type &systemMatrix, VectorType &rightHandSide) const
 This function modify the system matrix to apply a change of basis from the Cartesian coordinate system to the local coordinate system given by tangents and normals to the boundary. More...
 
template<typename VectorType >
void bcShiftToCartesianCoordSystem (matrix_Type &systemMatrix, VectorType &rightHandSide) const
 This function modify the system matrix to apply a change of basis from the local coordinate system given by tangents and normals to the boundary to the Cartesian coordinate system. More...
 
template<typename VectorType , typename MeshType >
void computeIntegratedNormals (const DOF &dof, CurrentFEManifold &currentBdFE, VectorType &normals, const MeshType &mesh)
 This method computes the normals to the domain boundary. More...
 
void exportToParaview (std::string fileName) const
 Export in the vtk format the tangential and normal vectors (t1, t2, n) for each DOF on the domain boundary. More...
 

Private Methods

template<typename MeshType >
void M_calculateCoordinates (const MeshType &mesh)
 stores the coordinates of the mesh vertices involved in the normal boundary conditions More...
 
void M_addBoundaryPoint (const ID &dofId, const ID &flag)
 Add point to the map of flags. More...
 
void M_addVersor (const ID &dofId, const Real &vx, const Real &vy, const Real &vz)
 Add point to the map of flags. More...
 
template<typename MeshType >
void M_calculateNormals (const MeshType &mesh, const DOF &dof, CurrentFEManifold &currentBdFE)
 Calculate the normal vectors. More...
 
void M_storeGivenVersors ()
 Store in *M_normalPtr the versors given by the user. More...
 
void M_calculateTangentVectors ()
 Compute the tangential vectors. More...
 
void M_buildRotationMatrix (const MapEpetra &systemMatrixMap, UInt offset=0)
 This function builds the rotation matrix for the change the basis. More...
 

Detailed Description

template<typename MatrixType>
class LifeV::BCManageNormal< MatrixType >

BCManageNormal - class for handling normal essential boundary conditions.

Author
Gwenol Grandperrin

The purpose of this class is to store the data and provide the methods needed to create the rotation matrix that help to impose normal essential boundary conditions.

Definition at line 56 of file BCManageNormal.hpp.

Member Typedef Documentation

◆ flagsMap_Type

typedef std::map<ID, ID> flagsMap_Type

Definition at line 62 of file BCManageNormal.hpp.

◆ versorsMap_Type

typedef std::map<ID, Vector > versorsMap_Type

Definition at line 63 of file BCManageNormal.hpp.

◆ matrix_Type

typedef MatrixType matrix_Type

Definition at line 64 of file BCManageNormal.hpp.

◆ matrixPtr_Type

typedef std::shared_ptr<matrix_Type> matrixPtr_Type

Definition at line 65 of file BCManageNormal.hpp.

◆ epetraMapPtr_Type

typedef std::shared_ptr<MapEpetra> epetraMapPtr_Type

Definition at line 66 of file BCManageNormal.hpp.

◆ epetraVectorPtr_type

typedef std::shared_ptr<VectorEpetra> epetraVectorPtr_type

Definition at line 67 of file BCManageNormal.hpp.

Constructor & Destructor Documentation

◆ BCManageNormal() [1/2]

Empty Constructor.

Definition at line 271 of file BCManageNormal.hpp.

◆ BCManageNormal() [2/2]

BCManageNormal ( BCManageNormal< MatrixType > const &  bcManageNormal)

Copy constructor.

All the stored data are copied so that the pointers of bcManageNormal and this do not share the same objects

Parameters
bcManageNormalBCManageNormal

Definition at line 281 of file BCManageNormal.hpp.

+ Here is the caller graph for this function:

◆ ~BCManageNormal()

Destructor.

Definition at line 299 of file BCManageNormal.hpp.

Member Function Documentation

◆ operator=()

BCManageNormal< MatrixType > & operator= ( const BCManageNormal< MatrixType > &  bcManageNormal)

Assignment operator.

Stored data are copied so that the pointers of bcManageNormal and of the returned class do not share the same objects

Parameters
bcManageNormalBCManageNormal
Returns
Reference to a new BCManageNormal with the same content of bcManageNormal

Definition at line 311 of file BCManageNormal.hpp.

◆ init()

void init ( const BCBase boundaryCondition,
const Real time 
)

Store boundary DOFs id and coordinates related to the boundary conditions.

This method store the coordinates and the flags of the DOFs on the boundary. In the directional case, it will also store the normals.

Parameters
boundaryConditionA BCBase object
timeThe actual time

Definition at line 335 of file BCManageNormal.hpp.

◆ build()

void build ( const MeshType mesh,
const DOF dof,
CurrentFEManifold currentBdFE,
const MapEpetra map,
UInt  offset 
)

Build the rotation matrix.

This method calculates the normal and tangential vectors and builds the rotation matrix

Parameters
meshThe mesh
dofThe DOF object
currentBdFEthe current boundary finite element
systemMatrixThe system matrix
offsetThe boundary condition offset
commPtrpointer to Epetra_Comm object

Definition at line 361 of file BCManageNormal.hpp.

◆ bcShiftToNormalTangentialCoordSystem()

void bcShiftToNormalTangentialCoordSystem ( matrix_Type systemMatrix,
VectorType &  rightHandSide 
) const

This function modify the system matrix to apply a change of basis from the Cartesian coordinate system to the local coordinate system given by tangents and normals to the boundary.

This method is used to change the basis of the DOF. Then it will be possible to prescribe an Essential boundary condition in the (t1, t2, n) direction by simply prescribing it along the x,y and z axis respectively. In fact it will perform the operations A := R*A*Rt and R*b, where A and be are the system matrix and right hand side, while R is the rotation matrix which transform the Cartesian system into the system of normals and tangents to the boundary.

Parameters
systemMatrixthe matrix of the problem
rightHandSidethe right hand side

Definition at line 414 of file BCManageNormal.hpp.

◆ bcShiftToCartesianCoordSystem()

void bcShiftToCartesianCoordSystem ( matrix_Type systemMatrix,
VectorType &  rightHandSide 
) const

This function modify the system matrix to apply a change of basis from the local coordinate system given by tangents and normals to the boundary to the Cartesian coordinate system.

This method is used to do the opposite operation of the method bcShiftToNormalTangentialCoordSystem. In fact it will perform the operations A := Rt*A*R and Rt*b, where A and be are the system matrix and right hand side, and R is the rotation matrix which transform the Cartesian system into the system of normals and tangents to the boundary.

Parameters
systemMatrixthe matrix of the problem
rightHandSidethe right hand side

Definition at line 437 of file BCManageNormal.hpp.

◆ computeIntegratedNormals()

void computeIntegratedNormals ( const DOF dof,
CurrentFEManifold currentBdFE,
VectorType &  normals,
const MeshType mesh 
)

This method computes the normals to the domain boundary.

Warning
the computed normals can be not compatible with the incompressibility condition
Parameters
dofThe DOF object containing the local-to-global map
currentBdFEThe current finite element on the boundary
normalsThe output vector
meshThe mesh

Definition at line 458 of file BCManageNormal.hpp.

◆ exportToParaview()

void exportToParaview ( std::string  fileName) const

Export in the vtk format the tangential and normal vectors (t1, t2, n) for each DOF on the domain boundary.

Parameters
fileNameThe name of the file where the data are exported

Definition at line 542 of file BCManageNormal.hpp.

◆ M_calculateCoordinates()

void M_calculateCoordinates ( const MeshType mesh)
private

stores the coordinates of the mesh vertices involved in the normal boundary conditions

Parameters
meshThe mesh

Definition at line 643 of file BCManageNormal.hpp.

◆ M_addBoundaryPoint()

void M_addBoundaryPoint ( const ID dofId,
const ID flag 
)
private

Add point to the map of flags.

Parameters
IDboundaryCondition A BCBase object
timeThe actual time

Definition at line 674 of file BCManageNormal.hpp.

+ Here is the caller graph for this function:

◆ M_addVersor()

void M_addVersor ( const ID dofId,
const Real vx,
const Real vy,
const Real vz 
)
private

Add point to the map of flags.

Parameters
IDdofId The point global id
nxThe x component of the unit vector
nyThe y component of the unit vector
nzThe z component of the unit vector

Definition at line 681 of file BCManageNormal.hpp.

+ Here is the caller graph for this function:

◆ M_calculateNormals()

void M_calculateNormals ( const MeshType mesh,
const DOF dof,
CurrentFEManifold currentBdFE 
)
private

Calculate the normal vectors.

Parameters
dofthe dof class
currentBdFEthe current boundary finite element

Definition at line 692 of file BCManageNormal.hpp.

◆ M_storeGivenVersors()

void M_storeGivenVersors ( )
private

Store in *M_normalPtr the versors given by the user.

Save the imposed normal vectors.

This function normalizes the versors in case they are not already normalized

Definition at line 711 of file BCManageNormal.hpp.

+ Here is the caller graph for this function:

◆ M_calculateTangentVectors()

void M_calculateTangentVectors ( )
private

Compute the tangential vectors.

Calculate the tangential vectors for each triad in the vector triad.

Definition at line 753 of file BCManageNormal.hpp.

+ Here is the caller graph for this function:

◆ M_buildRotationMatrix()

void M_buildRotationMatrix ( const MapEpetra systemMatrixMap,
UInt  offset = 0 
)
private

This function builds the rotation matrix for the change the basis.

Parameters
systemMatrixthe matrix of the problem
offsetthat will be used if there is more than one unknown to recover the global ID

Definition at line 841 of file BCManageNormal.hpp.

+ Here is the caller graph for this function:

Field Documentation

◆ M_dataBuilt

bool M_dataBuilt
private

true when there are stored normals

Definition at line 227 of file BCManageNormal.hpp.

◆ M_rotationMatrixPtr

matrixPtr_Type M_rotationMatrixPtr
private

Shared pointer to the rotation matrix and its transpose.

Definition at line 230 of file BCManageNormal.hpp.

◆ M_rotationMatrixTransposePtr

matrixPtr_Type M_rotationMatrixTransposePtr
private

Definition at line 231 of file BCManageNormal.hpp.

◆ M_localMapEpetraPtr

epetraMapPtr_Type M_localMapEpetraPtr
private

Shared pointer to the local Map.

Definition at line 234 of file BCManageNormal.hpp.

◆ M_firstTangentPtr

epetraVectorPtr_type M_firstTangentPtr
private

Shared pointer to the vector of the first tangents to the domain boundary.

Definition at line 237 of file BCManageNormal.hpp.

◆ M_secondTangentPtr

epetraVectorPtr_type M_secondTangentPtr
private

Shared pointer to the vector of the second tangents to the domain boundary.

Definition at line 240 of file BCManageNormal.hpp.

◆ M_normalPtr

epetraVectorPtr_type M_normalPtr
private

Shared pointer to the vector of the normals to the domain boundary.

Definition at line 243 of file BCManageNormal.hpp.

◆ M_coordPtr

epetraVectorPtr_type M_coordPtr
private

Shared coordinates of the point.

Definition at line 246 of file BCManageNormal.hpp.

◆ M_numDof

UInt M_numDof
private

Number of degree of freedom (of one component)

Definition at line 249 of file BCManageNormal.hpp.

◆ M_numInvoledDof

UInt M_numInvoledDof
private

Number of degree of freedom (of one component) involved in boundary conditions.

Definition at line 252 of file BCManageNormal.hpp.

◆ M_flags

flagsMap_Type M_flags
private

flag of the boundary elements

Definition at line 255 of file BCManageNormal.hpp.

◆ M_givenVersors

versorsMap_Type M_givenVersors
private

Versors that are given by the user.

Definition at line 258 of file BCManageNormal.hpp.


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