36 #ifndef ABSTRACTFECABLEINTEGRALASSEMBLER_HPP_ 37 #define ABSTRACTFECABLEINTEGRALASSEMBLER_HPP_ 39 #include "AbstractFeAssemblerCommon.hpp" 40 #include "GaussianQuadratureRule.hpp" 41 #include "PetscVecTools.hpp" 42 #include "PetscMatTools.hpp" 50 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM,
unsigned PROBLEM_DIM,
bool CAN_ASSEMBLE_VECTOR,
bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
84 const c_matrix<double, CABLE_ELEMENT_DIM, SPACE_DIM>& rInverseJacobian,
85 c_matrix<double, SPACE_DIM, NUM_CABLE_ELEMENT_NODES>& rReturnValue);
118 c_vector<double, NUM_CABLE_ELEMENT_NODES>& rPhi,
119 c_matrix<double, SPACE_DIM, NUM_CABLE_ELEMENT_NODES>& rGradPhi,
121 c_vector<double,PROBLEM_DIM>& rU,
122 c_matrix<double, PROBLEM_DIM, SPACE_DIM>& rGradU,
154 c_vector<double, NUM_CABLE_ELEMENT_NODES>& rPhi,
155 c_matrix<double, SPACE_DIM, NUM_CABLE_ELEMENT_NODES>& rGradPhi,
157 c_vector<double,PROBLEM_DIM>& rU,
158 c_matrix<double, PROBLEM_DIM, SPACE_DIM>& rGradU,
183 c_matrix<double, PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES, PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES >& rAElem,
184 c_vector<double, PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES>& rBElem);
217 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM,
unsigned PROBLEM_DIM,
bool CAN_ASSEMBLE_VECTOR,
bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
220 :
AbstractFeAssemblerCommon<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL>(),
224 assert(CAN_ASSEMBLE_VECTOR || CAN_ASSEMBLE_MATRIX);
232 assert(INTERPOLATION_LEVEL!=NONLINEAR);
236 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM,
unsigned PROBLEM_DIM,
bool CAN_ASSEMBLE_VECTOR,
bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
242 assemble_event = HeartEventHandler::ASSEMBLE_SYSTEM;
246 assemble_event = HeartEventHandler::ASSEMBLE_RHS;
251 EXCEPTION(
"Matrix to be assembled has not been set");
255 EXCEPTION(
"Vector to be assembled has not been set");
272 c_matrix<double, STENCIL_SIZE, STENCIL_SIZE> a_elem;
273 c_vector<double, STENCIL_SIZE> b_elem;
279 iter !=
mpMesh->GetCableElementIteratorEnd();
289 unsigned p_indices[STENCIL_SIZE];
294 PetscMatTools::AddMultipleValues<STENCIL_SIZE>(this->
mMatrixToAssemble, p_indices, a_elem);
299 PetscVecTools::AddMultipleValues<STENCIL_SIZE>(this->
mVectorToAssemble, p_indices, b_elem);
312 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM,
unsigned PROBLEM_DIM,
bool CAN_ASSEMBLE_VECTOR,
bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
315 const c_matrix<double, CABLE_ELEMENT_DIM, SPACE_DIM>& rInverseJacobian,
316 c_matrix<double, SPACE_DIM, NUM_CABLE_ELEMENT_NODES>& rReturnValue)
318 static c_matrix<double, CABLE_ELEMENT_DIM, NUM_CABLE_ELEMENT_NODES> grad_phi;
321 rReturnValue = prod(trans(rInverseJacobian), grad_phi);
324 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM,
unsigned PROBLEM_DIM,
bool CAN_ASSEMBLE_VECTOR,
bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
327 c_matrix<double, PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES, PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES >& rAElem,
328 c_vector<double, PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES>& rBElem)
335 c_matrix<double, SPACE_DIM, CABLE_ELEMENT_DIM> jacobian;
336 c_matrix<double, CABLE_ELEMENT_DIM, SPACE_DIM> inverse_jacobian;
337 double jacobian_determinant;
358 c_vector<double, NUM_CABLE_ELEMENT_NODES> phi;
359 c_matrix<double, SPACE_DIM, NUM_CABLE_ELEMENT_NODES> grad_phi;
377 c_vector<double,PROBLEM_DIM> u = zero_vector<double>(PROBLEM_DIM);
378 c_matrix<double,PROBLEM_DIM,SPACE_DIM> grad_u = zero_matrix<double>(PROBLEM_DIM,SPACE_DIM);
384 for (
unsigned i=0; i<num_nodes; i++)
388 if (INTERPOLATION_LEVEL != CARDIAC)
390 const c_vector<double, SPACE_DIM>& r_node_loc = p_node->
rGetLocation();
399 for (
unsigned index_of_unknown=0; index_of_unknown<(INTERPOLATION_LEVEL!=CARDIAC ? PROBLEM_DIM : 1); index_of_unknown++)
410 u(index_of_unknown) += phi(i)*u_at_node;
bool mZeroVectorBeforeAssembly
c_vector< double, DIM > & rGetLocation()
bool mZeroMatrixBeforeAssembly
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)
unsigned GetNodeGlobalIndex(unsigned localIndex) const
virtual double GetCurrentSolutionOrGuessValue(unsigned nodeIndex, unsigned indexOfUnknown)
#define EXCEPTION(message)
static void BeginEvent(unsigned event)
static const unsigned NUM_CABLE_ELEMENT_NODES
virtual void ResetInterpolatedQuantities()
static void ComputeBasisFunctionDerivatives(const ChastePoint< ELEMENT_DIM > &rPoint, c_matrix< double, ELEMENT_DIM, ELEMENT_DIM+1 > &rReturnValue)
unsigned GetNumQuadPoints() const
void CalculateInverseJacobian(c_matrix< double, SPACE_DIM, ELEMENT_DIM > &rJacobian, double &rJacobianDeterminant, c_matrix< double, ELEMENT_DIM, SPACE_DIM > &rInverseJacobian)
Node< SPACE_DIM > * GetNode(unsigned localIndex) const
virtual ~AbstractFeCableIntegralAssembler()
bool GetOwnership() const
static const unsigned CABLE_ELEMENT_DIM
virtual bool ElementAssemblyCriterion(Element< CABLE_ELEMENT_DIM, SPACE_DIM > &rElement)
unsigned GetNumNodes() const
void GetStiffnessMatrixGlobalIndices(unsigned problemDim, unsigned *pIndices) const
LinearBasisFunction< 1 > CableBasisFunction
static void ComputeBasisFunctions(const ChastePoint< ELEMENT_DIM > &rPoint, c_vector< double, ELEMENT_DIM+1 > &rReturnValue)
MixedDimensionMesh< ELEMENT_DIM, SPACE_DIM > * mpMesh
const c_vector< double, SPACE_DIM > & rGetLocation() const
virtual void IncrementInterpolatedQuantities(double phiI, const Node< SPACE_DIM > *pNode)
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)
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)
double GetWeight(unsigned index) const
AbstractFeCableIntegralAssembler(MixedDimensionMesh< ELEMENT_DIM, SPACE_DIM > *pMesh)
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)
static void EndEvent(unsigned event)
const ChastePoint< ELEMENT_DIM > & rGetQuadPoint(unsigned index) const
ReplicatableVector mCurrentSolutionOrGuessReplicated
GaussianQuadratureRule< 1 > * mpCableQuadRule
std::vector< Element< 1, SPACE_DIM > * >::const_iterator CableElementIterator