Chaste Release::3.1
AbstractFeVolumeIntegralAssembler.hpp
00001 /*
00002 
00003 Copyright (c) 2005-2012, 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 ABSTRACTFEVOLUMEINTEGRALASSEMBLER_HPP_
00037 #define ABSTRACTFEVOLUMEINTEGRALASSEMBLER_HPP_
00038 
00039 #include "AbstractFeAssemblerCommon.hpp"
00040 #include "GaussianQuadratureRule.hpp"
00041 #include "BoundaryConditionsContainer.hpp"
00042 #include "PetscVecTools.hpp"
00043 #include "PetscMatTools.hpp"
00044 
00077 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
00078 class AbstractFeVolumeIntegralAssembler :
00079      public AbstractFeAssemblerCommon<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM,CAN_ASSEMBLE_VECTOR,CAN_ASSEMBLE_MATRIX,INTERPOLATION_LEVEL>
00080 {
00081 protected:
00083     AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>* mpMesh;
00084 
00086     GaussianQuadratureRule<ELEMENT_DIM>* mpQuadRule;
00087 
00089     typedef LinearBasisFunction<ELEMENT_DIM> BasisFunction;
00090 
00110     void ComputeTransformedBasisFunctionDerivatives(const ChastePoint<ELEMENT_DIM>& rPoint,
00111                                                     const c_matrix<double, ELEMENT_DIM, SPACE_DIM>& rInverseJacobian,
00112                                                     c_matrix<double, SPACE_DIM, ELEMENT_DIM+1>& rReturnValue);
00113 
00119     void DoAssemble();
00120 
00121 protected:
00122 
00142     virtual c_matrix<double,PROBLEM_DIM*(ELEMENT_DIM+1),PROBLEM_DIM*(ELEMENT_DIM+1)> ComputeMatrixTerm(
00143         c_vector<double, ELEMENT_DIM+1>& rPhi,
00144         c_matrix<double, SPACE_DIM, ELEMENT_DIM+1>& rGradPhi,
00145         ChastePoint<SPACE_DIM>& rX,
00146         c_vector<double,PROBLEM_DIM>& rU,
00147         c_matrix<double, PROBLEM_DIM, SPACE_DIM>& rGradU,
00148         Element<ELEMENT_DIM,SPACE_DIM>* pElement)
00149     {
00150         // If this line is reached this means this method probably hasn't been over-ridden correctly in
00151         // the concrete class
00152         NEVER_REACHED;
00153         return zero_matrix<double>(PROBLEM_DIM*(ELEMENT_DIM+1),PROBLEM_DIM*(ELEMENT_DIM+1));
00154     }
00155 
00175     virtual c_vector<double,PROBLEM_DIM*(ELEMENT_DIM+1)> ComputeVectorTerm(
00176         c_vector<double, ELEMENT_DIM+1>& rPhi,
00177         c_matrix<double, SPACE_DIM, ELEMENT_DIM+1>& rGradPhi,
00178         ChastePoint<SPACE_DIM>& rX,
00179         c_vector<double,PROBLEM_DIM>& rU,
00180         c_matrix<double, PROBLEM_DIM, SPACE_DIM>& rGradU,
00181         Element<ELEMENT_DIM,SPACE_DIM>* pElement)
00182     {
00183         // If this line is reached this means this method probably hasn't been over-ridden correctly in
00184         // the concrete class
00185         NEVER_REACHED;
00186         return zero_vector<double>(PROBLEM_DIM*(ELEMENT_DIM+1));
00187     }
00188 
00189 
00204     virtual void AssembleOnElement(Element<ELEMENT_DIM,SPACE_DIM>& rElement,
00205                                    c_matrix<double, PROBLEM_DIM*(ELEMENT_DIM+1), PROBLEM_DIM*(ELEMENT_DIM+1) >& rAElem,
00206                                    c_vector<double, PROBLEM_DIM*(ELEMENT_DIM+1)>& rBElem);
00207 
00215     virtual bool ElementAssemblyCriterion(Element<ELEMENT_DIM,SPACE_DIM>& rElement)
00216     {
00217         return true;
00218     }
00219 
00220 
00221 public:
00222 
00230     AbstractFeVolumeIntegralAssembler(AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh, unsigned numQuadPoints=2);
00231 
00235     virtual ~AbstractFeVolumeIntegralAssembler()
00236     {
00237         delete mpQuadRule;
00238     }
00239 };
00240 
00241 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
00242 AbstractFeVolumeIntegralAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL>::AbstractFeVolumeIntegralAssembler(
00243             AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh, unsigned numQuadPoints)
00244     : AbstractFeAssemblerCommon<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL>(),
00245       mpMesh(pMesh)
00246 {
00247     assert(pMesh);
00248     assert(numQuadPoints > 0);
00249 
00250     mpQuadRule = new GaussianQuadratureRule<ELEMENT_DIM>(numQuadPoints);
00251 }
00252 
00253 
00254 
00255 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
00256 void AbstractFeVolumeIntegralAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL>::DoAssemble()
00257 {
00258     assert(this->mAssembleMatrix || this->mAssembleVector);
00259 
00260     HeartEventHandler::EventType assemble_event;
00261     if (this->mAssembleMatrix)
00262     {
00263         assemble_event = HeartEventHandler::ASSEMBLE_SYSTEM;
00264     }
00265     else
00266     {
00267         assemble_event = HeartEventHandler::ASSEMBLE_RHS;
00268     }
00269 
00270     if (this->mAssembleMatrix && this->mMatrixToAssemble==NULL)
00271     {
00272         EXCEPTION("Matrix to be assembled has not been set");
00273     }
00274     if (this->mAssembleVector && this->mVectorToAssemble==NULL)
00275     {
00276         EXCEPTION("Vector to be assembled has not been set");
00277     }
00278 
00279     HeartEventHandler::BeginEvent(assemble_event);
00280 
00281     // Zero the matrix/vector if it is to be assembled
00282     if (this->mAssembleVector && this->mZeroVectorBeforeAssembly)
00283     {
00284         PetscVecTools::Zero(this->mVectorToAssemble);
00285     }
00286     if (this->mAssembleMatrix && this->mZeroMatrixBeforeAssembly)
00287     {
00288         PetscMatTools::Zero(this->mMatrixToAssemble);
00289     }
00290 
00291     const size_t STENCIL_SIZE=PROBLEM_DIM*(ELEMENT_DIM+1);
00292     c_matrix<double, STENCIL_SIZE, STENCIL_SIZE> a_elem;
00293     c_vector<double, STENCIL_SIZE> b_elem;
00294 
00295     // Loop over elements
00296     for (typename AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::ElementIterator iter = mpMesh->GetElementIteratorBegin();
00297          iter != mpMesh->GetElementIteratorEnd();
00298          ++iter)
00299     {
00300         Element<ELEMENT_DIM, SPACE_DIM>& r_element = *iter;
00301 
00302         // Test for ownership first, since it's pointless to test the criterion on something which we might know nothing about.
00303         if ( r_element.GetOwnership() == true && ElementAssemblyCriterion(r_element)==true )
00304         {
00305             AssembleOnElement(r_element, a_elem, b_elem);
00306 
00307             unsigned p_indices[STENCIL_SIZE];
00308             r_element.GetStiffnessMatrixGlobalIndices(PROBLEM_DIM, p_indices);
00309 
00310             if (this->mAssembleMatrix)
00311             {
00312                 PetscMatTools::AddMultipleValues<STENCIL_SIZE>(this->mMatrixToAssemble, p_indices, a_elem);
00313             }
00314 
00315             if (this->mAssembleVector)
00316             {
00317                 PetscVecTools::AddMultipleValues<STENCIL_SIZE>(this->mVectorToAssemble, p_indices, b_elem);
00318             }
00319         }
00320     }
00321 
00322     HeartEventHandler::EndEvent(assemble_event);
00323 }
00324 
00325 
00327 // Implementation - AssembleOnElement and smaller
00329 
00330 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
00331 void AbstractFeVolumeIntegralAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL>::ComputeTransformedBasisFunctionDerivatives(
00332         const ChastePoint<ELEMENT_DIM>& rPoint,
00333         const c_matrix<double, ELEMENT_DIM, SPACE_DIM>& rInverseJacobian,
00334         c_matrix<double, SPACE_DIM, ELEMENT_DIM+1>& rReturnValue)
00335 {
00336     assert(ELEMENT_DIM < 4 && ELEMENT_DIM > 0);
00337     static c_matrix<double, ELEMENT_DIM, ELEMENT_DIM+1> grad_phi;
00338 
00339     LinearBasisFunction<ELEMENT_DIM>::ComputeBasisFunctionDerivatives(rPoint, grad_phi);
00340     rReturnValue = prod(trans(rInverseJacobian), grad_phi);
00341 }
00342 
00343 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
00344 void AbstractFeVolumeIntegralAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL>::AssembleOnElement(
00345     Element<ELEMENT_DIM,SPACE_DIM>& rElement,
00346     c_matrix<double, PROBLEM_DIM*(ELEMENT_DIM+1), PROBLEM_DIM*(ELEMENT_DIM+1) >& rAElem,
00347     c_vector<double, PROBLEM_DIM*(ELEMENT_DIM+1)>& rBElem)
00348 {
00354     c_matrix<double, SPACE_DIM, ELEMENT_DIM> jacobian;
00355     c_matrix<double, ELEMENT_DIM, SPACE_DIM> inverse_jacobian;
00356     double jacobian_determinant;
00357 
00358     mpMesh->GetInverseJacobianForElement(rElement.GetIndex(), jacobian, jacobian_determinant, inverse_jacobian);
00359 
00360     if (this->mAssembleMatrix)
00361     {
00362         rAElem.clear();
00363     }
00364 
00365     if (this->mAssembleVector)
00366     {
00367         rBElem.clear();
00368     }
00369 
00370     const unsigned num_nodes = rElement.GetNumNodes();
00371 
00372     // Allocate memory for the basis functions values and derivative values
00373     c_vector<double, ELEMENT_DIM+1> phi;
00374     c_matrix<double, SPACE_DIM, ELEMENT_DIM+1> grad_phi;
00375 
00376     // Loop over Gauss points
00377     for (unsigned quad_index=0; quad_index < mpQuadRule->GetNumQuadPoints(); quad_index++)
00378     {
00379         const ChastePoint<ELEMENT_DIM>& quad_point = mpQuadRule->rGetQuadPoint(quad_index);
00380 
00381         BasisFunction::ComputeBasisFunctions(quad_point, phi);
00382 
00383         if ( this->mAssembleMatrix || INTERPOLATION_LEVEL==NONLINEAR )
00384         {
00385             ComputeTransformedBasisFunctionDerivatives(quad_point, inverse_jacobian, grad_phi);
00386         }
00387 
00388         // Location of the Gauss point in the original element will be stored in x
00389         // Where applicable, u will be set to the value of the current solution at x
00390         ChastePoint<SPACE_DIM> x(0,0,0);
00391 
00392         c_vector<double,PROBLEM_DIM> u = zero_vector<double>(PROBLEM_DIM);
00393         c_matrix<double,PROBLEM_DIM,SPACE_DIM> grad_u = zero_matrix<double>(PROBLEM_DIM,SPACE_DIM);
00394 
00395         // Allow the concrete version of the assembler to interpolate any desired quantities
00396         this->ResetInterpolatedQuantities();
00397 
00398         // Interpolation
00399         for (unsigned i=0; i<num_nodes; i++)
00400         {
00401             const Node<SPACE_DIM>* p_node = rElement.GetNode(i);
00402 
00403             if (INTERPOLATION_LEVEL != CARDIAC) // don't even interpolate X if cardiac problem
00404             {
00405                 const c_vector<double, SPACE_DIM>& r_node_loc = p_node->rGetLocation();
00406                 // interpolate x
00407                 x.rGetLocation() += phi(i)*r_node_loc;
00408             }
00409 
00410             // Interpolate u and grad u if a current solution or guess exists
00411             unsigned node_global_index = rElement.GetNodeGlobalIndex(i);
00412             if (this->mCurrentSolutionOrGuessReplicated.GetSize() > 0)
00413             {
00414                 for (unsigned index_of_unknown=0; index_of_unknown<(INTERPOLATION_LEVEL!=CARDIAC ? PROBLEM_DIM : 1); index_of_unknown++)
00415                 {
00416                     /*
00417                      * If we have a solution (e.g. this is a dynamic problem) then
00418                      * interpolate the value at this quadrature point.
00419                      *
00420                      * NOTE: the following assumes that if, say, there are two unknowns
00421                      * u and v, they are stored in the current solution vector as
00422                      * [U1 V1 U2 V2 ... U_n V_n].
00423                      */
00424                     double u_at_node = this->GetCurrentSolutionOrGuessValue(node_global_index, index_of_unknown);
00425                     u(index_of_unknown) += phi(i)*u_at_node;
00426 
00427                     if (INTERPOLATION_LEVEL==NONLINEAR) // don't need to construct grad_phi or grad_u in other cases
00428                     {
00429                         for (unsigned j=0; j<SPACE_DIM; j++)
00430                         {
00431                             grad_u(index_of_unknown,j) += grad_phi(j,i)*u_at_node;
00432                         }
00433                     }
00434                 }
00435             }
00436 
00437             // Allow the concrete version of the assembler to interpolate any desired quantities
00438             this->IncrementInterpolatedQuantities(phi(i), p_node);
00439             if ( this->mAssembleMatrix || INTERPOLATION_LEVEL==NONLINEAR )
00440             {
00441                 this->IncrementInterpolatedGradientQuantities(grad_phi, i, p_node);
00442             }
00443         }
00444 
00445         double wJ = jacobian_determinant * mpQuadRule->GetWeight(quad_index);
00446 
00447         // Create rAElem and rBElem
00448         if (this->mAssembleMatrix)
00449         {
00450             noalias(rAElem) += ComputeMatrixTerm(phi, grad_phi, x, u, grad_u, &rElement) * wJ;
00451         }
00452 
00453         if (this->mAssembleVector)
00454         {
00455             noalias(rBElem) += ComputeVectorTerm(phi, grad_phi, x, u, grad_u, &rElement) * wJ;
00456         }
00457     }
00458 }
00459 
00460 
00461 #endif /*ABSTRACTFEVOLUMEINTEGRALASSEMBLER_HPP_*/