Chaste Commit::ca8ccdedf819b6e02855bc0e8e6f50bdecbc5208
|
#include <ExtendedBidomainMassMatrixAssembler.hpp>
Protected Member Functions | |
virtual c_matrix< double, 3 *(DIM+1), 3 *(DIM+1)> | ComputeMatrixTerm (c_vector< double, DIM+1 > &rPhi, c_matrix< double, DIM, DIM+1 > &rGradPhi, ChastePoint< DIM > &rX, c_vector< double, 3 > &rU, c_matrix< double, 3, DIM > &rGradU, Element< DIM, DIM > *pElement) |
Protected Member Functions inherited from AbstractFeVolumeIntegralAssembler< DIM, DIM, 3, false, true, NORMAL > | |
void | ComputeTransformedBasisFunctionDerivatives (const ChastePoint< ELEMENT_DIM > &rPoint, const c_matrix< double, ELEMENT_DIM, SPACE_DIM > &rInverseJacobian, c_matrix< double, SPACE_DIM, ELEMENT_DIM+1 > &rReturnValue) |
void | DoAssemble () |
virtual c_matrix< double, PROBLEM_DIM *(ELEMENT_DIM+1), PROBLEM_DIM *(ELEMENT_DIM+1)> | ComputeMatrixTerm (c_vector< double, ELEMENT_DIM+1 > &rPhi, c_matrix< double, SPACE_DIM, ELEMENT_DIM+1 > &rGradPhi, ChastePoint< SPACE_DIM > &rX, c_vector< double, PROBLEM_DIM > &rU, c_matrix< double, PROBLEM_DIM, SPACE_DIM > &rGradU, Element< ELEMENT_DIM, SPACE_DIM > *pElement) |
virtual c_vector< double, PROBLEM_DIM *(ELEMENT_DIM+1)> | ComputeVectorTerm (c_vector< double, ELEMENT_DIM+1 > &rPhi, c_matrix< double, SPACE_DIM, ELEMENT_DIM+1 > &rGradPhi, ChastePoint< SPACE_DIM > &rX, c_vector< double, PROBLEM_DIM > &rU, c_matrix< double, PROBLEM_DIM, SPACE_DIM > &rGradU, Element< ELEMENT_DIM, SPACE_DIM > *pElement) |
virtual void | AssembleOnElement (Element< ELEMENT_DIM, SPACE_DIM > &rElement, c_matrix< double, PROBLEM_DIM *(ELEMENT_DIM+1), PROBLEM_DIM *(ELEMENT_DIM+1) > &rAElem, c_vector< double, PROBLEM_DIM *(ELEMENT_DIM+1)> &rBElem) |
virtual bool | ElementAssemblyCriterion (Element< ELEMENT_DIM, SPACE_DIM > &rElement) |
Protected Member Functions inherited from AbstractFeAssemblerCommon< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL > | |
virtual double | GetCurrentSolutionOrGuessValue (unsigned nodeIndex, unsigned indexOfUnknown) |
virtual void | ResetInterpolatedQuantities () |
virtual void | IncrementInterpolatedQuantities (double phiI, const Node< SPACE_DIM > *pNode) |
virtual void | IncrementInterpolatedGradientQuantities (const c_matrix< double, SPACE_DIM, ELEMENT_DIM+1 > &rGradPhi, unsigned phiIndex, const Node< SPACE_DIM > *pNode) |
Additional Inherited Members | |
Protected Types inherited from AbstractFeVolumeIntegralAssembler< DIM, DIM, 3, false, true, NORMAL > | |
typedef LinearBasisFunction< ELEMENT_DIM > | BasisFunction |
Protected Attributes inherited from AbstractFeVolumeIntegralAssembler< DIM, DIM, 3, false, true, NORMAL > | |
AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > * | mpMesh |
GaussianQuadratureRule< ELEMENT_DIM > * | mpQuadRule |
Protected Attributes inherited from AbstractFeAssemblerCommon< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL > | |
ReplicatableVector | mCurrentSolutionOrGuessReplicated |
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 |
Constructs a matrix with the mass matrix in the voltage-voltage block.
Ie. IF the extended bidomain unknowns were ordered [phi1_1,..,phi1_n, phi2_1, ..., phi2_n, phie_1,..,phie_n], the matrix would be, in block form
[ M 0 0] [ 0 M 0] [ 0 0 M]
where M is the standard nxn mass matrix.
Since the bidomain ordering is not [phi1_1,..,phi1_n,phi2_1,..,phi2_n, phie_1,...phie_n] but [phi1_1,phi2_1,phie_1,..,phi1_n,phi2_n,phie_n], the matrix has a different form.
Definition at line 60 of file ExtendedBidomainMassMatrixAssembler.hpp.
|
inline |
Constructor
pMesh | pointer to the mesh |
Definition at line 93 of file ExtendedBidomainMassMatrixAssembler.hpp.
|
inlinevirtual |
Destructor.
Definition at line 101 of file ExtendedBidomainMassMatrixAssembler.hpp.
|
protectedvirtual |
This method is called by AssembleOnElement() and tells the assembler the contribution to add to the element stiffness matrix.
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 41 of file ExtendedBidomainMassMatrixAssembler.cpp.
References AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetUnsignedAttribute(), and HeartRegionCode::IsRegionBath().