Chaste Release::3.1
AbstractFeCableIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL > Class Template Reference

#include <AbstractFeCableIntegralAssembler.hpp>

Inheritance diagram for AbstractFeCableIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >:
Collaboration diagram for AbstractFeCableIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >:

List of all members.

Public Member Functions

 AbstractFeCableIntegralAssembler (MixedDimensionMesh< ELEMENT_DIM, SPACE_DIM > *pMesh, unsigned numQuadPoints=2)
virtual ~AbstractFeCableIntegralAssembler ()

Protected Types

typedef LinearBasisFunction< 1 > CableBasisFunction

Protected Member Functions

void ComputeTransformedBasisFunctionDerivatives (const ChastePoint< CABLE_ELEMENT_DIM > &rPoint, const c_matrix< double, CABLE_ELEMENT_DIM, SPACE_DIM > &rInverseJacobian, c_matrix< double, SPACE_DIM, NUM_CABLE_ELEMENT_NODES > &rReturnValue)
void DoAssemble ()
virtual c_matrix< double,
PROBLEM_DIM
*NUM_CABLE_ELEMENT_NODES,
PROBLEM_DIM
*NUM_CABLE_ELEMENT_NODES
ComputeCableMatrixTerm (c_vector< double, NUM_CABLE_ELEMENT_NODES > &rPhi, c_matrix< double, SPACE_DIM, NUM_CABLE_ELEMENT_NODES > &rGradPhi, ChastePoint< SPACE_DIM > &rX, c_vector< double, PROBLEM_DIM > &rU, c_matrix< double, PROBLEM_DIM, SPACE_DIM > &rGradU, Element< CABLE_ELEMENT_DIM, SPACE_DIM > *pElement)
virtual c_vector< double,
PROBLEM_DIM
*NUM_CABLE_ELEMENT_NODES
ComputeCableVectorTerm (c_vector< double, NUM_CABLE_ELEMENT_NODES > &rPhi, c_matrix< double, SPACE_DIM, NUM_CABLE_ELEMENT_NODES > &rGradPhi, ChastePoint< SPACE_DIM > &rX, c_vector< double, PROBLEM_DIM > &rU, c_matrix< double, PROBLEM_DIM, SPACE_DIM > &rGradU, Element< CABLE_ELEMENT_DIM, SPACE_DIM > *pElement)
virtual void AssembleOnCableElement (Element< CABLE_ELEMENT_DIM, SPACE_DIM > &rElement, c_matrix< double, PROBLEM_DIM *NUM_CABLE_ELEMENT_NODES, PROBLEM_DIM *NUM_CABLE_ELEMENT_NODES > &rAElem, c_vector< double, PROBLEM_DIM *NUM_CABLE_ELEMENT_NODES > &rBElem)
virtual bool ElementAssemblyCriterion (Element< CABLE_ELEMENT_DIM, SPACE_DIM > &rElement)

Protected Attributes

MixedDimensionMesh
< ELEMENT_DIM, SPACE_DIM > * 
mpMesh
GaussianQuadratureRule< 1 > * mpCableQuadRule

Static Protected Attributes

static const unsigned CABLE_ELEMENT_DIM = 1
static const unsigned NUM_CABLE_ELEMENT_NODES = 2

Detailed Description

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
class AbstractFeCableIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >

The class in similar to AbstractFeVolumeIntegralAssembler (see documentation for this), but is for creating a finite element matrices or vectors that involve integrals over CABLES, ie 1d regions in a 2d/3d mesh. Required for cardiac simulations with a Purkinje network. Uses a MixedDimensionMesh, which is composed of the normal mesh plus cables.

Definition at line 51 of file AbstractFeCableIntegralAssembler.hpp.


Member Typedef Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
typedef LinearBasisFunction<1> AbstractFeCableIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::CableBasisFunction [protected]

Basis function for use with normal elements.

Definition at line 67 of file AbstractFeCableIntegralAssembler.hpp.


Constructor & Destructor Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
AbstractFeCableIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::AbstractFeCableIntegralAssembler ( MixedDimensionMesh< ELEMENT_DIM, SPACE_DIM > *  pMesh,
unsigned  numQuadPoints = 2 
)

Constructor.

Parameters:
pMeshThe mesh
numQuadPointsThe number of quadratures points (in each dimension) to use per element. Defaults to 2.

Definition at line 218 of file AbstractFeCableIntegralAssembler.hpp.

References AbstractFeCableIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mpCableQuadRule.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
virtual AbstractFeCableIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::~AbstractFeCableIntegralAssembler ( ) [inline, virtual]

Member Function Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
void AbstractFeCableIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::AssembleOnCableElement ( Element< CABLE_ELEMENT_DIM, SPACE_DIM > &  rElement,
c_matrix< double, PROBLEM_DIM *NUM_CABLE_ELEMENT_NODES, PROBLEM_DIM *NUM_CABLE_ELEMENT_NODES > &  rAElem,
c_vector< double, PROBLEM_DIM *NUM_CABLE_ELEMENT_NODES > &  rBElem 
) [protected, virtual]

Calculate the contribution of a single cable element to the linear system.

Parameters:
rElementThe element to assemble on.
rAElemThe element's contribution to the LHS matrix is returned in this n by n matrix, where n is the no. of nodes in this element. There is no need to zero this matrix before calling.
rBElemThe element's contribution to the RHS vector is returned in this vector of length n, the no. of nodes in this element. There is no need to zero this vector before calling.

Called by AssembleSystem(). Calls ComputeCableMatrixTerm() etc.

Todo:
#1320 This assumes that the Jacobian is constant on an element. This is true for linear basis functions, but not for any other type of basis function.

Definition at line 324 of file AbstractFeCableIntegralAssembler.hpp.

References AbstractTetrahedralElement< ELEMENT_DIM, SPACE_DIM >::CalculateInverseJacobian(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNode(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNodeGlobalIndex(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNumNodes(), ChastePoint< DIM >::rGetLocation(), and Node< SPACE_DIM >::rGetLocation().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
virtual c_matrix<double,PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES,PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES> AbstractFeCableIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::ComputeCableMatrixTerm ( c_vector< double, NUM_CABLE_ELEMENT_NODES > &  rPhi,
c_matrix< double, SPACE_DIM, NUM_CABLE_ELEMENT_NODES > &  rGradPhi,
ChastePoint< SPACE_DIM > &  rX,
c_vector< double, PROBLEM_DIM > &  rU,
c_matrix< double, PROBLEM_DIM, SPACE_DIM > &  rGradU,
Element< CABLE_ELEMENT_DIM, SPACE_DIM > *  pElement 
) [inline, protected, virtual]

This method returns the matrix to be added to element stiffness matrix for a given Gauss point, ie, essentially the INTEGRAND in the integral definition of the matrix (integral over cable region).

The arguments are the bases, bases gradients, x and current solution computed at the Gauss point. The returned matrix will be multiplied by the Gauss weight and Jacobian determinant and added to the element stiffness matrix (see AssembleOnElement()).

** This method has to be implemented in the concrete class if CAN_ASSEMBLE_MATRIX is true. **

NOTE: for linear problems rGradU is NOT set up correctly because it should not be needed.

Parameters:
rPhiThe basis functions, rPhi(i) = phi_i, i=1..numBases.
rGradPhiBasis gradients, rGradPhi(i,j) = d(phi_j)/d(X_i).
rXThe point in space.
rUThe unknown as a vector, u(i) = u_i.
rGradUThe gradient of the unknown as a matrix, rGradU(i,j) = d(u_i)/d(X_j).
pElementPointer to the element.

Definition at line 116 of file AbstractFeCableIntegralAssembler.hpp.

References NEVER_REACHED, and AbstractFeCableIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::NUM_CABLE_ELEMENT_NODES.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
virtual c_vector<double,PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES> AbstractFeCableIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::ComputeCableVectorTerm ( c_vector< double, NUM_CABLE_ELEMENT_NODES > &  rPhi,
c_matrix< double, SPACE_DIM, NUM_CABLE_ELEMENT_NODES > &  rGradPhi,
ChastePoint< SPACE_DIM > &  rX,
c_vector< double, PROBLEM_DIM > &  rU,
c_matrix< double, PROBLEM_DIM, SPACE_DIM > &  rGradU,
Element< CABLE_ELEMENT_DIM, SPACE_DIM > *  pElement 
) [inline, protected, virtual]

This method returns the vector to be added to element stiffness vector for a given Gauss point, ie, essentially the INTEGRAND in the integral definition of the vector. The arguments are the bases, x and current solution computed at the Gauss point. The returned vector will be multiplied by the Gauss weight and Jacobian determinant and added to the element stiffness matrix (see AssembleOnElement()).

** This method has to be implemented in the concrete class if CAN_ASSEMBLE_VECTOR is true. **

NOTE: for linear problems rGradPhi and rGradU are NOT set up correctly because they should not be needed.

Parameters:
rPhiThe basis functions, rPhi(i) = phi_i, i=1..numBases
rGradPhiBasis gradients, rGradPhi(i,j) = d(phi_j)/d(X_i)
rXThe point in space
rUThe unknown as a vector, u(i) = u_i
rGradUThe gradient of the unknown as a matrix, rGradU(i,j) = d(u_i)/d(X_j)
pElementPointer to the element

Definition at line 150 of file AbstractFeCableIntegralAssembler.hpp.

References NEVER_REACHED, and AbstractFeCableIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::NUM_CABLE_ELEMENT_NODES.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
void AbstractFeCableIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::ComputeTransformedBasisFunctionDerivatives ( const ChastePoint< CABLE_ELEMENT_DIM > &  rPoint,
const c_matrix< double, CABLE_ELEMENT_DIM, SPACE_DIM > &  rInverseJacobian,
c_matrix< double, SPACE_DIM, NUM_CABLE_ELEMENT_NODES > &  rReturnValue 
) [protected]

Compute the derivatives of all basis functions at a point within a cable element. This method will transform the results, for use within Gaussian quadrature for example.

Parameters:
rPointThe point at which to compute the basis functions. The results are undefined if this is not within the canonical element.
rInverseJacobianThe inverse of the Jacobian matrix mapping the real element into the canonical element.
rReturnValueA reference to a vector, to be filled in
Returns:
The derivatives of the basis functions, in local index order. Each entry is a vector (c_vector<double, SPACE_DIM> instance) giving the derivative along each axis.

Definition at line 312 of file AbstractFeCableIntegralAssembler.hpp.

References LinearBasisFunction< ELEMENT_DIM >::ComputeBasisFunctionDerivatives().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
void AbstractFeCableIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::DoAssemble ( ) [protected, virtual]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
virtual bool AbstractFeCableIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::ElementAssemblyCriterion ( Element< CABLE_ELEMENT_DIM, SPACE_DIM > &  rElement) [inline, protected, virtual]

Whether to include this (cable) element when assembling. Returns true here but can be overridden by the concrete assembler if not all elements should be included.

Parameters:
rElementthe element

Definition at line 190 of file AbstractFeCableIntegralAssembler.hpp.


Member Data Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
const unsigned AbstractFeCableIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::CABLE_ELEMENT_DIM = 1 [static, protected]

Cable element dimension.

Definition at line 55 of file AbstractFeCableIntegralAssembler.hpp.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
GaussianQuadratureRule<1>* AbstractFeCableIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mpCableQuadRule [protected]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
MixedDimensionMesh<ELEMENT_DIM, SPACE_DIM>* AbstractFeCableIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mpMesh [protected]

Mesh to be solved on.

Definition at line 61 of file AbstractFeCableIntegralAssembler.hpp.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
const unsigned AbstractFeCableIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::NUM_CABLE_ELEMENT_NODES = 2 [static, protected]

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