36 #ifndef ABSTRACTFEVOLUMEINTEGRALASSEMBLER_HPP_ 37 #define ABSTRACTFEVOLUMEINTEGRALASSEMBLER_HPP_ 39 #include "AbstractFeAssemblerCommon.hpp" 40 #include "GaussianQuadratureRule.hpp" 41 #include "BoundaryConditionsContainer.hpp" 42 #include "PetscVecTools.hpp" 43 #include "PetscMatTools.hpp" 77 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM,
unsigned PROBLEM_DIM,
bool CAN_ASSEMBLE_VECTOR,
bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
79 public AbstractFeAssemblerCommon<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM,CAN_ASSEMBLE_VECTOR,CAN_ASSEMBLE_MATRIX,INTERPOLATION_LEVEL>
111 const c_matrix<double, ELEMENT_DIM, SPACE_DIM>& rInverseJacobian,
112 c_matrix<double, SPACE_DIM, ELEMENT_DIM+1>& rReturnValue);
143 virtual c_matrix<double,PROBLEM_DIM*(ELEMENT_DIM+1),PROBLEM_DIM*(ELEMENT_DIM+1)>
ComputeMatrixTerm(
144 c_vector<double, ELEMENT_DIM+1>& rPhi,
145 c_matrix<double, SPACE_DIM, ELEMENT_DIM+1>& rGradPhi,
147 c_vector<double,PROBLEM_DIM>& rU,
148 c_matrix<double, PROBLEM_DIM, SPACE_DIM>& rGradU,
154 return zero_matrix<double>(PROBLEM_DIM*(ELEMENT_DIM+1),PROBLEM_DIM*(ELEMENT_DIM+1));
179 c_vector<double, ELEMENT_DIM+1>& rPhi,
180 c_matrix<double, SPACE_DIM, ELEMENT_DIM+1>& rGradPhi,
182 c_vector<double,PROBLEM_DIM>& rU,
183 c_matrix<double, PROBLEM_DIM, SPACE_DIM>& rGradU,
189 return zero_vector<double>(PROBLEM_DIM*(ELEMENT_DIM+1));
209 c_matrix<
double, PROBLEM_DIM*(ELEMENT_DIM+1), PROBLEM_DIM*(ELEMENT_DIM+1) >& rAElem,
210 c_vector<
double, PROBLEM_DIM*(ELEMENT_DIM+1)>& rBElem);
243 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM,
unsigned PROBLEM_DIM,
bool CAN_ASSEMBLE_VECTOR,
bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
246 :
AbstractFeAssemblerCommon<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL>(),
256 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM,
unsigned PROBLEM_DIM,
bool CAN_ASSEMBLE_VECTOR,
bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
264 assemble_event = HeartEventHandler::ASSEMBLE_SYSTEM;
268 assemble_event = HeartEventHandler::ASSEMBLE_RHS;
273 EXCEPTION(
"Matrix to be assembled has not been set");
277 EXCEPTION(
"Vector to be assembled has not been set");
292 const size_t STENCIL_SIZE=PROBLEM_DIM*(ELEMENT_DIM+1);
293 c_matrix<double, STENCIL_SIZE, STENCIL_SIZE> a_elem;
294 c_vector<double, STENCIL_SIZE> b_elem;
298 iter !=
mpMesh->GetElementIteratorEnd();
308 unsigned p_indices[STENCIL_SIZE];
313 PetscMatTools::AddMultipleValues<STENCIL_SIZE>(this->
mMatrixToAssemble, p_indices, a_elem);
318 PetscVecTools::AddMultipleValues<STENCIL_SIZE>(this->
mVectorToAssemble, p_indices, b_elem);
331 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM,
unsigned PROBLEM_DIM,
bool CAN_ASSEMBLE_VECTOR,
bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
334 const c_matrix<double, ELEMENT_DIM, SPACE_DIM>& rInverseJacobian,
335 c_matrix<double, SPACE_DIM, ELEMENT_DIM+1>& rReturnValue)
337 assert(ELEMENT_DIM < 4 && ELEMENT_DIM > 0);
338 static c_matrix<double, ELEMENT_DIM, ELEMENT_DIM+1> grad_phi;
341 rReturnValue = prod(trans(rInverseJacobian), grad_phi);
344 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM,
unsigned PROBLEM_DIM,
bool CAN_ASSEMBLE_VECTOR,
bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
347 c_matrix<
double, PROBLEM_DIM*(ELEMENT_DIM+1), PROBLEM_DIM*(ELEMENT_DIM+1) >& rAElem,
348 c_vector<
double, PROBLEM_DIM*(ELEMENT_DIM+1)>& rBElem)
355 c_matrix<double, SPACE_DIM, ELEMENT_DIM> jacobian;
356 c_matrix<double, ELEMENT_DIM, SPACE_DIM> inverse_jacobian;
357 double jacobian_determinant;
359 mpMesh->GetInverseJacobianForElement(rElement.
GetIndex(), jacobian, jacobian_determinant, inverse_jacobian);
374 c_vector<double, ELEMENT_DIM+1> phi;
375 c_matrix<double, SPACE_DIM, ELEMENT_DIM+1> grad_phi;
378 for (
unsigned quad_index=0; quad_index <
mpQuadRule->GetNumQuadPoints(); quad_index++)
393 c_vector<double,PROBLEM_DIM> u = zero_vector<double>(PROBLEM_DIM);
394 c_matrix<double,PROBLEM_DIM,SPACE_DIM> grad_u = zero_matrix<double>(PROBLEM_DIM,SPACE_DIM);
400 for (
unsigned i=0; i<num_nodes; i++)
404 if (INTERPOLATION_LEVEL != CARDIAC)
406 const c_vector<double, SPACE_DIM>& r_node_loc = p_node->
rGetLocation();
415 for (
unsigned index_of_unknown=0; index_of_unknown<(INTERPOLATION_LEVEL!=CARDIAC ? PROBLEM_DIM : 1); index_of_unknown++)
426 u(index_of_unknown) += phi(i)*u_at_node;
428 if (INTERPOLATION_LEVEL==NONLINEAR)
430 for (
unsigned j=0; j<SPACE_DIM; j++)
432 grad_u(index_of_unknown,j) += grad_phi(j,i)*u_at_node;
446 double wJ = jacobian_determinant *
mpQuadRule->GetWeight(quad_index);
451 noalias(rAElem) +=
ComputeMatrixTerm(phi, grad_phi, x, u, grad_u, &rElement) * wJ;
456 noalias(rBElem) +=
ComputeVectorTerm(phi, grad_phi, x, u, grad_u, &rElement) * wJ;
bool mZeroVectorBeforeAssembly
virtual ~AbstractFeVolumeIntegralAssembler()
c_vector< double, DIM > & rGetLocation()
bool mZeroMatrixBeforeAssembly
unsigned GetNodeGlobalIndex(unsigned localIndex) const
virtual double GetCurrentSolutionOrGuessValue(unsigned nodeIndex, unsigned indexOfUnknown)
#define EXCEPTION(message)
static void BeginEvent(unsigned event)
virtual void ResetInterpolatedQuantities()
virtual void IncrementInterpolatedGradientQuantities(const c_matrix< double, SPACE_DIM, ELEMENT_DIM+1 > &rGradPhi, unsigned phiIndex, const Node< SPACE_DIM > *pNode)
LinearBasisFunction< ELEMENT_DIM > BasisFunction
static void ComputeBasisFunctionDerivatives(const ChastePoint< ELEMENT_DIM > &rPoint, c_matrix< double, ELEMENT_DIM, ELEMENT_DIM+1 > &rReturnValue)
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)
Node< SPACE_DIM > * GetNode(unsigned localIndex) const
AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > * mpMesh
bool GetOwnership() const
GaussianQuadratureRule< ELEMENT_DIM > * mpQuadRule
unsigned GetNumNodes() const
void GetStiffnessMatrixGlobalIndices(unsigned problemDim, unsigned *pIndices) const
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)
static void ComputeBasisFunctions(const ChastePoint< ELEMENT_DIM > &rPoint, c_vector< double, ELEMENT_DIM+1 > &rReturnValue)
virtual bool ElementAssemblyCriterion(Element< ELEMENT_DIM, SPACE_DIM > &rElement)
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)
const c_vector< double, SPACE_DIM > & rGetLocation() const
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 void IncrementInterpolatedQuantities(double phiI, const Node< SPACE_DIM > *pNode)
unsigned GetIndex() const
static void EndEvent(unsigned event)
AbstractFeVolumeIntegralAssembler(AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh)
ReplicatableVector mCurrentSolutionOrGuessReplicated