AbstractFunctionalCalculator.hpp

00001 /*
00002 
00003 Copyright (C) University of Oxford, 2005-2011
00004 
00005 University of Oxford means the Chancellor, Masters and Scholars of the
00006 University of Oxford, having an administrative office at Wellington
00007 Square, Oxford OX1 2JD, UK.
00008 
00009 This file is part of Chaste.
00010 
00011 Chaste is free software: you can redistribute it and/or modify it
00012 under the terms of the GNU Lesser General Public License as published
00013 by the Free Software Foundation, either version 2.1 of the License, or
00014 (at your option) any later version.
00015 
00016 Chaste is distributed in the hope that it will be useful, but WITHOUT
00017 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018 FITNESS FOR A PARTICULAR PURPOSE.  See the GNU Lesser General Public
00019 License for more details. The offer of Chaste under the terms of the
00020 License is subject to the License being interpreted in accordance with
00021 English Law and subject to any action against the University of Oxford
00022 being under the jurisdiction of the English Courts.
00023 
00024 You should have received a copy of the GNU Lesser General Public License
00025 along with Chaste. If not, see <http://www.gnu.org/licenses/>.
00026 
00027 */
00028 
00029 #ifndef ABSTRACTFUNCTIONALCALCULATOR_HPP_
00030 #define ABSTRACTFUNCTIONALCALCULATOR_HPP_
00031 
00032 #include "LinearBasisFunction.hpp"
00033 #include "GaussianQuadratureRule.hpp"
00034 #include "AbstractTetrahedralMesh.hpp"
00035 #include "GaussianQuadratureRule.hpp"
00036 #include "ReplicatableVector.hpp"
00037 
00050 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00051 class AbstractFunctionalCalculator
00052 {
00053 private:
00054 
00056     ReplicatableVector mSolutionReplicated;
00057 
00065     virtual double GetIntegrand(ChastePoint<SPACE_DIM>& rX,
00066                                 c_vector<double,PROBLEM_DIM>& rU,
00067                                 c_matrix<double,PROBLEM_DIM,SPACE_DIM>& rGradU)=0;
00068 
00069 public:
00070 
00074     virtual ~AbstractFunctionalCalculator()
00075     {
00076     }
00077 
00088     double Calculate(AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>& rMesh, Vec solution);
00089 
00095     double CalculateOnElement(Element<ELEMENT_DIM,SPACE_DIM>& rElement);
00096 };
00097 
00099 // Implementation
00101 
00102 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00103 double AbstractFunctionalCalculator<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::CalculateOnElement(Element<ELEMENT_DIM, SPACE_DIM>& rElement)
00104 {
00105     double result_on_element = 0;
00106 
00107     GaussianQuadratureRule<ELEMENT_DIM> quad_rule(2);
00108 
00111     double jacobian_determinant;
00112     c_matrix<double, SPACE_DIM, ELEMENT_DIM> jacobian;
00113     c_matrix<double, ELEMENT_DIM, SPACE_DIM> inverse_jacobian;
00114     rElement.CalculateInverseJacobian(jacobian, jacobian_determinant, inverse_jacobian);
00115 
00116     const unsigned num_nodes = rElement.GetNumNodes();
00117 
00118     // Loop over Gauss points
00119     for (unsigned quad_index=0; quad_index < quad_rule.GetNumQuadPoints(); quad_index++)
00120     {
00121         const ChastePoint<ELEMENT_DIM>& quad_point = quad_rule.rGetQuadPoint(quad_index);
00122 
00123         c_vector<double, ELEMENT_DIM+1> phi;
00124         LinearBasisFunction<ELEMENT_DIM>::ComputeBasisFunctions(quad_point, phi);
00125         c_matrix<double, ELEMENT_DIM, ELEMENT_DIM+1> grad_phi;
00126         LinearBasisFunction<ELEMENT_DIM>::ComputeTransformedBasisFunctionDerivatives(quad_point, inverse_jacobian, grad_phi);
00127 
00128         // Location of the Gauss point in the original element will be stored in x
00129         ChastePoint<SPACE_DIM> x(0,0,0);
00130         c_vector<double,PROBLEM_DIM> u = zero_vector<double>(PROBLEM_DIM);
00131         c_matrix<double,PROBLEM_DIM,SPACE_DIM> grad_u = zero_matrix<double>(PROBLEM_DIM,SPACE_DIM);
00132 
00133         for (unsigned i=0; i<num_nodes; i++)
00134         {
00135             const c_vector<double, SPACE_DIM>& r_node_loc = rElement.GetNode(i)->rGetLocation();
00136 
00137             // Interpolate x
00138             x.rGetLocation() += phi(i)*r_node_loc;
00139 
00140             // Interpolate u and grad u
00141             unsigned node_global_index = rElement.GetNodeGlobalIndex(i);
00142             for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
00143             {
00144                 // NOTE - following assumes that, if say there are two unknowns u and v, they
00145                 // are stored in the current solution vector as
00146                 // [U1 V1 U2 V2 ... U_n V_n]
00147                 unsigned index_into_vec = PROBLEM_DIM*node_global_index + index_of_unknown;
00148 
00149                 double u_at_node = mSolutionReplicated[index_into_vec];
00150                 u(index_of_unknown) += phi(i)*u_at_node;
00151                 for (unsigned j=0; j<SPACE_DIM; j++)
00152                 {
00153                     grad_u(index_of_unknown,j) += grad_phi(j,i)*u_at_node;
00154                 }
00155             }
00156         }
00157 
00158         double wJ = jacobian_determinant * quad_rule.GetWeight(quad_index);
00159         result_on_element += GetIntegrand(x, u, grad_u) * wJ;
00160     }
00161 
00162     return result_on_element;
00163 }
00164 
00165 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00166 double AbstractFunctionalCalculator<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::Calculate(AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>& rMesh, Vec solution)
00167 {
00168     assert(solution);
00169     mSolutionReplicated.ReplicatePetscVector(solution);
00170     if (mSolutionReplicated.GetSize() != rMesh.GetNumNodes() * PROBLEM_DIM)
00171     {
00172         EXCEPTION("The solution size does not match the mesh");
00173     }
00174 
00175     double local_result = 0;
00176 
00177     try
00178     {
00179         for (typename AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::ElementIterator iter = rMesh.GetElementIteratorBegin();
00180              iter != rMesh.GetElementIteratorEnd();
00181              ++iter)
00182         {
00183             if (rMesh.CalculateDesignatedOwnershipOfElement((*iter).GetIndex()) == true)
00184             {
00185                 local_result += CalculateOnElement(*iter);
00186             }
00187         }
00188     }
00189     catch (Exception &exception_in_integral)
00190     {
00191         PetscTools::ReplicateException(true);
00192         throw exception_in_integral;
00193     }
00194     PetscTools::ReplicateException(false);
00195 
00196     double final_result;
00197     MPI_Allreduce(&local_result, &final_result, 1, MPI_DOUBLE, MPI_SUM, PETSC_COMM_WORLD);
00198     return final_result;
00199 }
00200 
00201 #endif /*ABSTRACTFUNCTIONALCALCULATOR_HPP_*/
Generated on Thu Dec 22 13:00:17 2011 for Chaste by  doxygen 1.6.3