AbstractFunctionalCalculator.hpp

00001 /*
00002 
00003 Copyright (c) 2005-2015, University of Oxford.
00004 All rights reserved.
00005 
00006 University of Oxford means the Chancellor, Masters and Scholars of the
00007 University of Oxford, having an administrative office at Wellington
00008 Square, Oxford OX1 2JD, UK.
00009 
00010 This file is part of Chaste.
00011 
00012 Redistribution and use in source and binary forms, with or without
00013 modification, are permitted provided that the following conditions are met:
00014  * Redistributions of source code must retain the above copyright notice,
00015    this list of conditions and the following disclaimer.
00016  * Redistributions in binary form must reproduce the above copyright notice,
00017    this list of conditions and the following disclaimer in the documentation
00018    and/or other materials provided with the distribution.
00019  * Neither the name of the University of Oxford nor the names of its
00020    contributors may be used to endorse or promote products derived from this
00021    software without specific prior written permission.
00022 
00023 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
00024 AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00025 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
00026 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
00027 LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
00028 CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
00029 GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
00030 HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
00031 LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
00032 OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
00033 
00034 */
00035 
00036 #ifndef ABSTRACTFUNCTIONALCALCULATOR_HPP_
00037 #define ABSTRACTFUNCTIONALCALCULATOR_HPP_
00038 
00039 #include "LinearBasisFunction.hpp"
00040 #include "GaussianQuadratureRule.hpp"
00041 #include "AbstractTetrahedralMesh.hpp"
00042 #include "ReplicatableVector.hpp"
00043 
00057 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00058 class AbstractFunctionalCalculator
00059 {
00060 private:
00061 
00063     ReplicatableVector mSolutionReplicated;
00064 
00072     virtual double GetIntegrand(ChastePoint<SPACE_DIM>& rX,
00073                                 c_vector<double,PROBLEM_DIM>& rU,
00074                                 c_matrix<double,PROBLEM_DIM,SPACE_DIM>& rGradU)=0;
00075 
00082     virtual bool ShouldSkipThisElement(Element<ELEMENT_DIM,SPACE_DIM>& rElement);
00083 
00084 public:
00085 
00089     virtual ~AbstractFunctionalCalculator()
00090     {
00091     }
00092 
00103     double Calculate(AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>& rMesh, Vec solution);
00104 
00110     double CalculateOnElement(Element<ELEMENT_DIM,SPACE_DIM>& rElement);
00111 };
00112 
00114 // Implementation
00116 
00117 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00118 double AbstractFunctionalCalculator<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::CalculateOnElement(Element<ELEMENT_DIM, SPACE_DIM>& rElement)
00119 {
00120     double result_on_element = 0;
00121 
00122     // Third order quadrature.  Note that the functional may be non-polynomial (see documentation of class).
00123     GaussianQuadratureRule<ELEMENT_DIM> quad_rule(3);
00124 
00127     double jacobian_determinant;
00128     c_matrix<double, SPACE_DIM, ELEMENT_DIM> jacobian;
00129     c_matrix<double, ELEMENT_DIM, SPACE_DIM> inverse_jacobian;
00130     rElement.CalculateInverseJacobian(jacobian, jacobian_determinant, inverse_jacobian);
00131 
00132     const unsigned num_nodes = rElement.GetNumNodes();
00133 
00134     // Loop over Gauss points
00135     for (unsigned quad_index=0; quad_index < quad_rule.GetNumQuadPoints(); quad_index++)
00136     {
00137         const ChastePoint<ELEMENT_DIM>& quad_point = quad_rule.rGetQuadPoint(quad_index);
00138 
00139         c_vector<double, ELEMENT_DIM+1> phi;
00140         LinearBasisFunction<ELEMENT_DIM>::ComputeBasisFunctions(quad_point, phi);
00141         c_matrix<double, ELEMENT_DIM, ELEMENT_DIM+1> grad_phi;
00142         LinearBasisFunction<ELEMENT_DIM>::ComputeTransformedBasisFunctionDerivatives(quad_point, inverse_jacobian, grad_phi);
00143 
00144         // Location of the Gauss point in the original element will be stored in x
00145         ChastePoint<SPACE_DIM> x(0,0,0);
00146         c_vector<double,PROBLEM_DIM> u = zero_vector<double>(PROBLEM_DIM);
00147         c_matrix<double,PROBLEM_DIM,SPACE_DIM> grad_u = zero_matrix<double>(PROBLEM_DIM,SPACE_DIM);
00148 
00149         for (unsigned i=0; i<num_nodes; i++)
00150         {
00151             const c_vector<double, SPACE_DIM>& r_node_loc = rElement.GetNode(i)->rGetLocation();
00152 
00153             // Interpolate x
00154             x.rGetLocation() += phi(i)*r_node_loc;
00155 
00156             // Interpolate u and grad u
00157             unsigned node_global_index = rElement.GetNodeGlobalIndex(i);
00158             for (unsigned index_of_unknown=0; index_of_unknown<PROBLEM_DIM; index_of_unknown++)
00159             {
00160                 // NOTE - following assumes that, if say there are two unknowns u and v, they
00161                 // are stored in the current solution vector as
00162                 // [U1 V1 U2 V2 ... U_n V_n]
00163                 unsigned index_into_vec = PROBLEM_DIM*node_global_index + index_of_unknown;
00164 
00165                 double u_at_node = mSolutionReplicated[index_into_vec];
00166                 u(index_of_unknown) += phi(i)*u_at_node;
00167                 for (unsigned j=0; j<SPACE_DIM; j++)
00168                 {
00169                     grad_u(index_of_unknown,j) += grad_phi(j,i)*u_at_node;
00170                 }
00171             }
00172         }
00173 
00174         double wJ = jacobian_determinant * quad_rule.GetWeight(quad_index);
00175         result_on_element += GetIntegrand(x, u, grad_u) * wJ;
00176     }
00177 
00178     return result_on_element;
00179 }
00180 
00181 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00182 double AbstractFunctionalCalculator<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::Calculate(AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>& rMesh, Vec solution)
00183 {
00184     assert(solution);
00185     mSolutionReplicated.ReplicatePetscVector(solution);
00186     if (mSolutionReplicated.GetSize() != rMesh.GetNumNodes() * PROBLEM_DIM)
00187     {
00188         EXCEPTION("The solution size does not match the mesh");
00189     }
00190 
00191     double local_result = 0;
00192 
00193     try
00194     {
00195         for (typename AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::ElementIterator iter = rMesh.GetElementIteratorBegin();
00196              iter != rMesh.GetElementIteratorEnd();
00197              ++iter)
00198         {
00199             if (rMesh.CalculateDesignatedOwnershipOfElement((*iter).GetIndex()) == true && !ShouldSkipThisElement(*iter))
00200             {
00201                 local_result += CalculateOnElement(*iter);
00202             }
00203         }
00204     }
00205     catch (Exception &exception_in_integral)
00206     {
00207         PetscTools::ReplicateException(true);
00208         throw exception_in_integral;
00209     }
00210     PetscTools::ReplicateException(false);
00211 
00212     double final_result;
00213     MPI_Allreduce(&local_result, &final_result, 1, MPI_DOUBLE, MPI_SUM, PETSC_COMM_WORLD);
00214     return final_result;
00215 }
00216 
00217 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00218 bool AbstractFunctionalCalculator<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::ShouldSkipThisElement(Element<ELEMENT_DIM,SPACE_DIM>& rElement)
00219 {
00220     return false;
00221 }
00222 
00223 #endif /*ABSTRACTFUNCTIONALCALCULATOR_HPP_*/

Generated by  doxygen 1.6.2