Chaste Release::3.1
|
#include <AbstractFeCableIntegralAssembler.hpp>
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 |
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.
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.
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.
pMesh | The mesh |
numQuadPoints | The number of quadratures points (in each dimension) to use per element. Defaults to 2. |
Definition at line 218 of file AbstractFeCableIntegralAssembler.hpp.
virtual AbstractFeCableIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::~AbstractFeCableIntegralAssembler | ( | ) | [inline, virtual] |
Destructor.
Definition at line 209 of file AbstractFeCableIntegralAssembler.hpp.
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.
rElement | The element to assemble on. |
rAElem | The 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. |
rBElem | The 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.
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().
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.
rPhi | The basis functions, rPhi(i) = phi_i, i=1..numBases. |
rGradPhi | Basis gradients, rGradPhi(i,j) = d(phi_j)/d(X_i). |
rX | The point in space. |
rU | The unknown as a vector, u(i) = u_i. |
rGradU | The gradient of the unknown as a matrix, rGradU(i,j) = d(u_i)/d(X_j). |
pElement | Pointer 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.
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.
rPhi | The basis functions, rPhi(i) = phi_i, i=1..numBases |
rGradPhi | Basis gradients, rGradPhi(i,j) = d(phi_j)/d(X_i) |
rX | The point in space |
rU | The unknown as a vector, u(i) = u_i |
rGradU | The gradient of the unknown as a matrix, rGradU(i,j) = d(u_i)/d(X_j) |
pElement | Pointer 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.
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.
rPoint | The point at which to compute the basis functions. The results are undefined if this is not within the canonical element. |
rInverseJacobian | The inverse of the Jacobian matrix mapping the real element into the canonical element. |
rReturnValue | A reference to a vector, to be filled in |
Definition at line 312 of file AbstractFeCableIntegralAssembler.hpp.
References LinearBasisFunction< ELEMENT_DIM >::ComputeBasisFunctionDerivatives().
void AbstractFeCableIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::DoAssemble | ( | ) | [protected, virtual] |
The main assembly method. Protected, should only be called through Assemble(), AssembleMatrix() or AssembleVector() which set mAssembleMatrix, mAssembleVector accordingly.
Implements AbstractFeAssemblerInterface< CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX >.
Definition at line 236 of file AbstractFeCableIntegralAssembler.hpp.
References GenericEventHandler< 16, HeartEventHandler >::BeginEvent(), GenericEventHandler< 16, HeartEventHandler >::EndEvent(), EXCEPTION, AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetOwnership(), AbstractTetrahedralElement< ELEMENT_DIM, SPACE_DIM >::GetStiffnessMatrixGlobalIndices(), PetscMatTools::Zero(), and PetscVecTools::Zero().
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.
rElement | the element |
Definition at line 190 of file AbstractFeCableIntegralAssembler.hpp.
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.
GaussianQuadratureRule<1>* AbstractFeCableIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::mpCableQuadRule [protected] |
Quadrature rule for use on cable elements.
Definition at line 64 of file AbstractFeCableIntegralAssembler.hpp.
Referenced by AbstractFeCableIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::AbstractFeCableIntegralAssembler(), and AbstractFeCableIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::~AbstractFeCableIntegralAssembler().
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.
const unsigned AbstractFeCableIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::NUM_CABLE_ELEMENT_NODES = 2 [static, protected] |
Number of nodes in a cable element.
Definition at line 58 of file AbstractFeCableIntegralAssembler.hpp.
Referenced by AbstractFeCableIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::ComputeCableMatrixTerm(), and AbstractFeCableIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL >::ComputeCableVectorTerm().