Chaste Commit::f2ff7ee04e70ac9d06c57344df8d017dbb12b97b
StokesFlowAssembler< DIM > Class Template Reference

#include <StokesFlowAssembler.hpp>

+ Inheritance diagram for StokesFlowAssembler< DIM >:
+ Collaboration diagram for StokesFlowAssembler< DIM >:

Public Member Functions

 StokesFlowAssembler (AbstractTetrahedralMesh< DIM, DIM > *pMesh, StokesFlowProblemDefinition< DIM > *pProblemDefinition)
 
- Public Member Functions inherited from AbstractContinuumMechanicsAssembler< DIM, true, true >
 AbstractContinuumMechanicsAssembler (AbstractTetrahedralMesh< DIM, DIM > *pMesh)
 
virtual ~AbstractContinuumMechanicsAssembler ()
 
- Public Member Functions inherited from AbstractFeAssemblerInterface< CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX >
 AbstractFeAssemblerInterface ()
 
void SetMatrixToAssemble (Mat &rMatToAssemble, bool zeroMatrixBeforeAssembly=true)
 
void SetVectorToAssemble (Vec &rVecToAssemble, bool zeroVectorBeforeAssembly)
 
void Assemble ()
 
void AssembleMatrix ()
 
void AssembleVector ()
 
virtual ~AbstractFeAssemblerInterface ()
 

Private Member Functions

c_matrix< double, SPATIAL_BLOCK_SIZE_ELEMENTAL, SPATIAL_BLOCK_SIZE_ELEMENTALComputeSpatialSpatialMatrixTerm (c_vector< double, NUM_NODES_PER_ELEMENT > &rQuadPhi, c_matrix< double, DIM, NUM_NODES_PER_ELEMENT > &rGradQuadPhi, c_vector< double, DIM > &rX, Element< DIM, DIM > *pElement)
 
c_matrix< double, SPATIAL_BLOCK_SIZE_ELEMENTAL, PRESSURE_BLOCK_SIZE_ELEMENTALComputeSpatialPressureMatrixTerm (c_vector< double, NUM_NODES_PER_ELEMENT > &rQuadPhi, c_matrix< double, DIM, NUM_NODES_PER_ELEMENT > &rGradQuadPhi, c_vector< double, NUM_VERTICES_PER_ELEMENT > &rLinearPhi, c_matrix< double, DIM, NUM_VERTICES_PER_ELEMENT > &rGradLinearPhi, c_vector< double, DIM > &rX, Element< DIM, DIM > *pElement)
 
c_vector< double, SPATIAL_BLOCK_SIZE_ELEMENTALComputeSpatialVectorTerm (c_vector< double, NUM_NODES_PER_ELEMENT > &rQuadPhi, c_matrix< double, DIM, NUM_NODES_PER_ELEMENT > &rGradQuadPhi, c_vector< double, DIM > &rX, Element< DIM, DIM > *pElement)
 

Private Attributes

StokesFlowProblemDefinition< DIM > * mpProblemDefinition
 
double mScaleFactor
 

Static Private Attributes

static const unsigned NUM_VERTICES_PER_ELEMENT = DIM+1
 
static const unsigned NUM_NODES_PER_ELEMENT = (DIM+1)*(DIM+2)/2
 
static const unsigned SPATIAL_BLOCK_SIZE_ELEMENTAL = DIM*NUM_NODES_PER_ELEMENT
 
static const unsigned PRESSURE_BLOCK_SIZE_ELEMENTAL = NUM_VERTICES_PER_ELEMENT
 

Friends

class TestStokesFlowAssembler
 

Additional Inherited Members

- Protected Member Functions inherited from AbstractContinuumMechanicsAssembler< DIM, true, true >
void DoAssemble ()
 
virtual c_matrix< double, PRESSURE_BLOCK_SIZE_ELEMENTAL, PRESSURE_BLOCK_SIZE_ELEMENTALComputePressurePressureMatrixTerm (c_vector< double, NUM_VERTICES_PER_ELEMENT > &rLinearPhi, c_matrix< double, DIM, NUM_VERTICES_PER_ELEMENT > &rGradLinearPhi, c_vector< double, DIM > &rX, Element< DIM, DIM > *pElement)
 
virtual c_vector< double, PRESSURE_BLOCK_SIZE_ELEMENTALComputePressureVectorTerm (c_vector< double, NUM_VERTICES_PER_ELEMENT > &rLinearPhi, c_matrix< double, DIM, NUM_VERTICES_PER_ELEMENT > &rGradLinearPhi, c_vector< double, DIM > &rX, Element< DIM, DIM > *pElement)
 
void AssembleOnElement (Element< DIM, DIM > &rElement, c_matrix< double, STENCIL_SIZE, STENCIL_SIZE > &rAElem, c_vector< double, STENCIL_SIZE > &rBElem)
 
- Protected Attributes inherited from AbstractContinuumMechanicsAssembler< DIM, true, true >
AbstractTetrahedralMesh< DIM, DIM > * mpMesh
 
GaussianQuadratureRule< DIM > * mpQuadRule
 
- Protected Attributes inherited from AbstractFeAssemblerInterface< CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX >
Vec mVectorToAssemble
 
Mat mMatrixToAssemble
 
bool mAssembleMatrix
 
bool mAssembleVector
 
bool mZeroMatrixBeforeAssembly
 
bool mZeroVectorBeforeAssembly
 
PetscInt mOwnershipRangeLo
 
PetscInt mOwnershipRangeHi
 
- Static Protected Attributes inherited from AbstractContinuumMechanicsAssembler< DIM, true, true >
static const bool BLOCK_SYMMETRIC_MATRIX
 
static const unsigned NUM_VERTICES_PER_ELEMENT
 
static const unsigned NUM_NODES_PER_ELEMENT
 
static const unsigned SPATIAL_BLOCK_SIZE_ELEMENTAL
 
static const unsigned PRESSURE_BLOCK_SIZE_ELEMENTAL
 
static const unsigned STENCIL_SIZE
 

Detailed Description

template<unsigned DIM>
class StokesFlowAssembler< DIM >

Assembler for setting up (volume-integral parts of) the matrix and vector used in the FEM discretisation of Stokes' Flow.

The matrix has the block-form (except see comment below) [A B] [B^T 0] and the vector has the block form (except see comment below) [b] [0]

NOTE: The elemental matrix and vector is as above. The full matrix and vector uses a completely different ordering: for parallelisation reasons the pressure variables are interleaved with the spatial variables and dummy pressure variables are used for internal nodes. For example, in 2d, the ordering is [U1 V1 P1 , .. , Un Vn, Pn] where n is the total number of nodes.

Definition at line 62 of file StokesFlowAssembler.hpp.

Constructor & Destructor Documentation

◆ StokesFlowAssembler()

template<unsigned DIM>
StokesFlowAssembler< DIM >::StokesFlowAssembler ( AbstractTetrahedralMesh< DIM, DIM > *  pMesh,
StokesFlowProblemDefinition< DIM > *  pProblemDefinition 
)
inline

Constructor

Parameters
pMesh
pProblemDefinition

Definition at line 251 of file StokesFlowAssembler.hpp.

Member Function Documentation

◆ ComputeSpatialPressureMatrixTerm()

template<unsigned DIM>
c_matrix< double, SPATIAL_BLOCK_SIZE_ELEMENTAL, PRESSURE_BLOCK_SIZE_ELEMENTAL > StokesFlowAssembler< DIM >::ComputeSpatialPressureMatrixTerm ( c_vector< double, NUM_NODES_PER_ELEMENT > &  rQuadPhi,
c_matrix< double, DIM, NUM_NODES_PER_ELEMENT > &  rGradQuadPhi,
c_vector< double, NUM_VERTICES_PER_ELEMENT > &  rLinearPhi,
c_matrix< double, DIM, NUM_VERTICES_PER_ELEMENT > &  rGradLinearPhi,
c_vector< double, DIM > &  rX,
Element< DIM, DIM > *  pElement 
)
inlineprivatevirtual

The matrix has the form (except see comments about ordering above) [A B ] [B^T 0 ] The function is related to the spatial-pressure block, ie matrix B.

For the (volume-integral) contribution to B from a given element, this method returns the INTEGRAND in the definition of B.

Parameters
rQuadPhiAll the quadratic basis functions on this element, evaluated at the current quad point
rGradQuadPhiGradients of all the quadratic basis functions on this element, evaluated at the current quad point
rLinearPhiAll the linear basis functions on this element, evaluated at the current quad point
rGradLinearPhiGradients of all the linear basis functions on this element, evaluated at the current quad point
rXCurrent location (physical position corresponding to quad point)
pElementCurrent element
Returns
stencil matrix

Reimplemented from AbstractContinuumMechanicsAssembler< DIM, true, true >.

Definition at line 168 of file StokesFlowAssembler.hpp.

References StokesFlowAssembler< DIM >::NUM_NODES_PER_ELEMENT, StokesFlowAssembler< DIM >::NUM_VERTICES_PER_ELEMENT, StokesFlowAssembler< DIM >::PRESSURE_BLOCK_SIZE_ELEMENTAL, and StokesFlowAssembler< DIM >::SPATIAL_BLOCK_SIZE_ELEMENTAL.

◆ ComputeSpatialSpatialMatrixTerm()

template<unsigned DIM>
c_matrix< double, SPATIAL_BLOCK_SIZE_ELEMENTAL, SPATIAL_BLOCK_SIZE_ELEMENTAL > StokesFlowAssembler< DIM >::ComputeSpatialSpatialMatrixTerm ( c_vector< double, NUM_NODES_PER_ELEMENT > &  rQuadPhi,
c_matrix< double, DIM, NUM_NODES_PER_ELEMENT > &  rGradQuadPhi,
c_vector< double, DIM > &  rX,
Element< DIM, DIM > *  pElement 
)
inlineprivatevirtual

The matrix has the form (except see comments about ordering above) [A B ] [B^T 0 ] The function is related to the spatial-spatial block, ie matrix A.

For the (volume-integral) contribution to A from a given element, this method returns the INTEGRAND in the definition of A.

Parameters
rQuadPhiAll the quadratic basis functions on this element, evaluated at the current quad point
rGradQuadPhiGradients of all the quadratic basis functions on this element, evaluated at the current quad point
rXCurrent location (physical position corresponding to quad point)
pElementCurrent element
Returns
stencil matrix

Reimplemented from AbstractContinuumMechanicsAssembler< DIM, true, true >.

Definition at line 113 of file StokesFlowAssembler.hpp.

References StokesFlowAssembler< DIM >::mpProblemDefinition, StokesFlowAssembler< DIM >::mScaleFactor, StokesFlowAssembler< DIM >::NUM_NODES_PER_ELEMENT, and StokesFlowAssembler< DIM >::SPATIAL_BLOCK_SIZE_ELEMENTAL.

◆ ComputeSpatialVectorTerm()

template<unsigned DIM>
c_vector< double, SPATIAL_BLOCK_SIZE_ELEMENTAL > StokesFlowAssembler< DIM >::ComputeSpatialVectorTerm ( c_vector< double, NUM_NODES_PER_ELEMENT > &  rQuadPhi,
c_matrix< double, DIM, NUM_NODES_PER_ELEMENT > &  rGradQuadPhi,
c_vector< double, DIM > &  rX,
Element< DIM, DIM > *  pElement 
)
inlineprivatevirtual

The matrix has the form (except see comments about ordering above) [A B ] [B^T 0 ] and the vector has the form [b1] [b2] The function is related to the spatial-block in the vector, ie b1.

For the contribution to b1 from a given element, this method should return the INTEGRAND in the definition of b1.

Parameters
rQuadPhiAll the quadratic basis functions on this element, evaluated at the current quad point
rGradQuadPhiGradients of all the quadratic basis functions on this element, evaluated at the current quad point
rXCurrent location (physical position)
pElementCurrent element
Returns
stencil vector

Implements AbstractContinuumMechanicsAssembler< DIM, true, true >.

Definition at line 217 of file StokesFlowAssembler.hpp.

References StokesFlowAssembler< DIM >::mpProblemDefinition, StokesFlowAssembler< DIM >::NUM_NODES_PER_ELEMENT, and StokesFlowAssembler< DIM >::SPATIAL_BLOCK_SIZE_ELEMENTAL.

Friends And Related Symbol Documentation

◆ TestStokesFlowAssembler

template<unsigned DIM>
friend class TestStokesFlowAssembler
friend

Definition at line 64 of file StokesFlowAssembler.hpp.

Member Data Documentation

◆ mpProblemDefinition

template<unsigned DIM>
StokesFlowProblemDefinition<DIM>* StokesFlowAssembler< DIM >::mpProblemDefinition
private

◆ mScaleFactor

template<unsigned DIM>
double StokesFlowAssembler< DIM >::mScaleFactor
private

This variable is initialised to 1.0 and almost never changed, and is used in the spatial-spatial matrix term. One test (see TestStokesFlowAssembler) sets it to 0.0 before assembling the matrix. Basically, a different weak form USED to be implemented here (corresponding to one kind of Neumann boundary condition), and the matrix for that weak form was compared against exact solutions in this test. mScaleFactor = 0.0 corresponds to old weak form, mScaleFactor = 1 corresponds to new weak form as documented in fem implementations pdf.

Definition at line 95 of file StokesFlowAssembler.hpp.

Referenced by StokesFlowAssembler< DIM >::ComputeSpatialSpatialMatrixTerm().

◆ NUM_NODES_PER_ELEMENT

template<unsigned DIM>
const unsigned StokesFlowAssembler< DIM >::NUM_NODES_PER_ELEMENT = (DIM+1)*(DIM+2)/2
staticprivate

◆ NUM_VERTICES_PER_ELEMENT

template<unsigned DIM>
const unsigned StokesFlowAssembler< DIM >::NUM_VERTICES_PER_ELEMENT = DIM+1
staticprivate

Number of vertices per element

Definition at line 68 of file StokesFlowAssembler.hpp.

Referenced by StokesFlowAssembler< DIM >::ComputeSpatialPressureMatrixTerm().

◆ PRESSURE_BLOCK_SIZE_ELEMENTAL

template<unsigned DIM>
const unsigned StokesFlowAssembler< DIM >::PRESSURE_BLOCK_SIZE_ELEMENTAL = NUM_VERTICES_PER_ELEMENT
staticprivate

Size of the pressure block, per element, equal num_vertices times 1 (as p solved for at each vertex.

Definition at line 83 of file StokesFlowAssembler.hpp.

Referenced by StokesFlowAssembler< DIM >::ComputeSpatialPressureMatrixTerm().

◆ SPATIAL_BLOCK_SIZE_ELEMENTAL

template<unsigned DIM>
const unsigned StokesFlowAssembler< DIM >::SPATIAL_BLOCK_SIZE_ELEMENTAL = DIM*NUM_NODES_PER_ELEMENT
staticprivate

Size of the spatial block, per element, equal num_nodes times DIM (as say in 2d (u,v) solved for at each node

Definition at line 77 of file StokesFlowAssembler.hpp.

Referenced by StokesFlowAssembler< DIM >::ComputeSpatialPressureMatrixTerm(), StokesFlowAssembler< DIM >::ComputeSpatialSpatialMatrixTerm(), and StokesFlowAssembler< DIM >::ComputeSpatialVectorTerm().


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