36 #ifndef ABSTRACTFUNCTIONALCALCULATOR_HPP_ 37 #define ABSTRACTFUNCTIONALCALCULATOR_HPP_ 39 #include "LinearBasisFunction.hpp" 40 #include "GaussianQuadratureRule.hpp" 41 #include "AbstractTetrahedralMesh.hpp" 42 #include "ReplicatableVector.hpp" 57 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM,
unsigned PROBLEM_DIM>
73 c_vector<double,PROBLEM_DIM>& rU,
74 c_matrix<double,PROBLEM_DIM,SPACE_DIM>& rGradU)=0;
117 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM,
unsigned PROBLEM_DIM>
120 double result_on_element = 0;
127 double jacobian_determinant;
128 c_matrix<double, SPACE_DIM, ELEMENT_DIM> jacobian;
129 c_matrix<double, ELEMENT_DIM, SPACE_DIM> inverse_jacobian;
135 for (
unsigned quad_index=0; quad_index < quad_rule.
GetNumQuadPoints(); quad_index++)
139 c_vector<double, ELEMENT_DIM+1> phi;
141 c_matrix<double, ELEMENT_DIM, ELEMENT_DIM+1> grad_phi;
146 c_vector<double,PROBLEM_DIM> u = zero_vector<double>(PROBLEM_DIM);
147 c_matrix<double,PROBLEM_DIM,SPACE_DIM> grad_u = zero_matrix<double>(PROBLEM_DIM,SPACE_DIM);
149 for (
unsigned i=0; i<num_nodes; i++)
151 const c_vector<double, SPACE_DIM>& r_node_loc = rElement.
GetNode(i)->rGetLocation();
158 for (
unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
163 unsigned index_into_vec = PROBLEM_DIM*node_global_index + index_of_unknown;
166 u(index_of_unknown) += phi(i)*u_at_node;
167 for (
unsigned j=0; j<SPACE_DIM; j++)
169 grad_u(index_of_unknown,j) += grad_phi(j,i)*u_at_node;
174 double wJ = jacobian_determinant * quad_rule.
GetWeight(quad_index);
178 return result_on_element;
181 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM,
unsigned PROBLEM_DIM>
188 EXCEPTION(
"The solution size does not match the mesh");
191 double local_result = 0;
208 throw exception_in_integral;
213 MPI_Allreduce(&local_result, &final_result, 1, MPI_DOUBLE, MPI_SUM, PETSC_COMM_WORLD);
217 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM,
unsigned PROBLEM_DIM>
ElementIterator GetElementIteratorBegin(bool skipDeletedElements=true)
virtual bool ShouldSkipThisElement(Element< ELEMENT_DIM, SPACE_DIM > &rElement)
c_vector< double, DIM > & rGetLocation()
unsigned GetNodeGlobalIndex(unsigned localIndex) const
ElementIterator GetElementIteratorEnd()
#define EXCEPTION(message)
double CalculateOnElement(Element< ELEMENT_DIM, SPACE_DIM > &rElement)
virtual unsigned GetNumNodes() const
virtual ~AbstractFunctionalCalculator()
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
static void ComputeTransformedBasisFunctionDerivatives(const ChastePoint< ELEMENT_DIM > &rPoint, const c_matrix< double, ELEMENT_DIM, ELEMENT_DIM > &rInverseJacobian, c_matrix< double, ELEMENT_DIM, ELEMENT_DIM+1 > &rReturnValue)
double Calculate(AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > &rMesh, Vec solution)
unsigned GetNumNodes() const
static void ComputeBasisFunctions(const ChastePoint< ELEMENT_DIM > &rPoint, c_vector< double, ELEMENT_DIM+1 > &rReturnValue)
ReplicatableVector mSolutionReplicated
void ReplicatePetscVector(Vec vec)
double GetWeight(unsigned index) const
virtual bool CalculateDesignatedOwnershipOfElement(unsigned elementIndex)
virtual double GetIntegrand(ChastePoint< SPACE_DIM > &rX, c_vector< double, PROBLEM_DIM > &rU, c_matrix< double, PROBLEM_DIM, SPACE_DIM > &rGradU)=0
const ChastePoint< ELEMENT_DIM > & rGetQuadPoint(unsigned index) const