AbstractFeVolumeIntegralAssembler.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 ABSTRACTFEVOLUMEINTEGRALASSEMBLER_HPP_
00030 #define ABSTRACTFEVOLUMEINTEGRALASSEMBLER_HPP_
00031 
00032 #include "AbstractFeAssemblerCommon.hpp"
00033 #include "GaussianQuadratureRule.hpp"
00034 #include "BoundaryConditionsContainer.hpp"
00035 #include "PetscVecTools.hpp"
00036 #include "PetscMatTools.hpp"
00037 
00070 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
00071 class AbstractFeVolumeIntegralAssembler :
00072      public AbstractFeAssemblerCommon<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM,CAN_ASSEMBLE_VECTOR,CAN_ASSEMBLE_MATRIX,INTERPOLATION_LEVEL>
00073 {
00074 protected:
00076     AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>* mpMesh;
00077 
00079     GaussianQuadratureRule<ELEMENT_DIM>* mpQuadRule;
00080 
00082     typedef LinearBasisFunction<ELEMENT_DIM> BasisFunction;
00083 
00103     void ComputeTransformedBasisFunctionDerivatives(const ChastePoint<ELEMENT_DIM>& rPoint,
00104                                                     const c_matrix<double, ELEMENT_DIM, SPACE_DIM>& rInverseJacobian,
00105                                                     c_matrix<double, SPACE_DIM, ELEMENT_DIM+1>& rReturnValue);
00106 
00112     void DoAssemble();
00113 
00114 protected:
00115 
00135     virtual c_matrix<double,PROBLEM_DIM*(ELEMENT_DIM+1),PROBLEM_DIM*(ELEMENT_DIM+1)> ComputeMatrixTerm(
00136         c_vector<double, ELEMENT_DIM+1>& rPhi,
00137         c_matrix<double, SPACE_DIM, ELEMENT_DIM+1>& rGradPhi,
00138         ChastePoint<SPACE_DIM>& rX,
00139         c_vector<double,PROBLEM_DIM>& rU,
00140         c_matrix<double, PROBLEM_DIM, SPACE_DIM>& rGradU,
00141         Element<ELEMENT_DIM,SPACE_DIM>* pElement)
00142     {
00143         // If this line is reached this means this method probably hasn't been over-ridden correctly in
00144         // the concrete class
00145         NEVER_REACHED;
00146         return zero_matrix<double>(PROBLEM_DIM*(ELEMENT_DIM+1),PROBLEM_DIM*(ELEMENT_DIM+1));
00147     }
00148 
00168     virtual c_vector<double,PROBLEM_DIM*(ELEMENT_DIM+1)> ComputeVectorTerm(
00169         c_vector<double, ELEMENT_DIM+1>& rPhi,
00170         c_matrix<double, SPACE_DIM, ELEMENT_DIM+1>& rGradPhi,
00171         ChastePoint<SPACE_DIM>& rX,
00172         c_vector<double,PROBLEM_DIM>& rU,
00173         c_matrix<double, PROBLEM_DIM, SPACE_DIM>& rGradU,
00174         Element<ELEMENT_DIM,SPACE_DIM>* pElement)
00175     {
00176         // If this line is reached this means this method probably hasn't been over-ridden correctly in
00177         // the concrete class
00178         NEVER_REACHED;
00179         return zero_vector<double>(PROBLEM_DIM*(ELEMENT_DIM+1));
00180     }
00181 
00182 
00197     virtual void AssembleOnElement(Element<ELEMENT_DIM,SPACE_DIM>& rElement,
00198                                    c_matrix<double, PROBLEM_DIM*(ELEMENT_DIM+1), PROBLEM_DIM*(ELEMENT_DIM+1) >& rAElem,
00199                                    c_vector<double, PROBLEM_DIM*(ELEMENT_DIM+1)>& rBElem);
00200 
00208     virtual bool ElementAssemblyCriterion(Element<ELEMENT_DIM,SPACE_DIM>& rElement)
00209     {
00210         return true;
00211     }
00212 
00213 
00214 public:
00215 
00223     AbstractFeVolumeIntegralAssembler(AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh, unsigned numQuadPoints=2);
00224 
00228     virtual ~AbstractFeVolumeIntegralAssembler()
00229     {
00230         delete mpQuadRule;
00231     }
00232 };
00233 
00234 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
00235 AbstractFeVolumeIntegralAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL>::AbstractFeVolumeIntegralAssembler(
00236             AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh, unsigned numQuadPoints)
00237     : AbstractFeAssemblerCommon<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL>(),
00238       mpMesh(pMesh)
00239 {
00240     assert(pMesh);
00241     assert(numQuadPoints > 0);
00242 
00243     mpQuadRule = new GaussianQuadratureRule<ELEMENT_DIM>(numQuadPoints);
00244 }
00245 
00246 
00247 
00248 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
00249 void AbstractFeVolumeIntegralAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL>::DoAssemble()
00250 {
00251     assert(this->mAssembleMatrix || this->mAssembleVector);
00252 
00253     HeartEventHandler::EventType assemble_event;
00254     if (this->mAssembleMatrix)
00255     {
00256         assemble_event = HeartEventHandler::ASSEMBLE_SYSTEM;
00257     }
00258     else
00259     {
00260         assemble_event = HeartEventHandler::ASSEMBLE_RHS;
00261     }
00262 
00263     if (this->mAssembleMatrix && this->mMatrixToAssemble==NULL)
00264     {
00265         EXCEPTION("Matrix to be assembled has not been set");
00266     }
00267     if (this->mAssembleVector && this->mVectorToAssemble==NULL)
00268     {
00269         EXCEPTION("Vector to be assembled has not been set");
00270     }
00271 
00272     HeartEventHandler::BeginEvent(assemble_event);
00273 
00274     // Zero the matrix/vector if it is to be assembled
00275     if (this->mAssembleVector && this->mZeroVectorBeforeAssembly)
00276     {
00277         PetscVecTools::Zero(this->mVectorToAssemble);
00278     }
00279     if (this->mAssembleMatrix && this->mZeroMatrixBeforeAssembly)
00280     {
00281         PetscMatTools::Zero(this->mMatrixToAssemble);
00282     }
00283 
00284     const size_t STENCIL_SIZE=PROBLEM_DIM*(ELEMENT_DIM+1);
00285     c_matrix<double, STENCIL_SIZE, STENCIL_SIZE> a_elem;
00286     c_vector<double, STENCIL_SIZE> b_elem;
00287 
00288     // Loop over elements
00289     for (typename AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::ElementIterator iter = mpMesh->GetElementIteratorBegin();
00290          iter != mpMesh->GetElementIteratorEnd();
00291          ++iter)
00292     {
00293         Element<ELEMENT_DIM, SPACE_DIM>& r_element = *iter;
00294 
00295         // Test for ownership first, since it's pointless to test the criterion on something which we might know nothing about.
00296         if ( r_element.GetOwnership() == true && ElementAssemblyCriterion(r_element)==true )
00297         {
00298             AssembleOnElement(r_element, a_elem, b_elem);
00299 
00300             unsigned p_indices[STENCIL_SIZE];
00301             r_element.GetStiffnessMatrixGlobalIndices(PROBLEM_DIM, p_indices);
00302 
00303             if (this->mAssembleMatrix)
00304             {
00305                 PetscMatTools::AddMultipleValues<STENCIL_SIZE>(this->mMatrixToAssemble, p_indices, a_elem);
00306             }
00307 
00308             if (this->mAssembleVector)
00309             {
00310                 PetscVecTools::AddMultipleValues<STENCIL_SIZE>(this->mVectorToAssemble, p_indices, b_elem);
00311             }
00312         }
00313     }
00314 
00315     HeartEventHandler::EndEvent(assemble_event);
00316 }
00317 
00318 
00320 // Implementation - AssembleOnElement and smaller
00322 
00323 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
00324 void AbstractFeVolumeIntegralAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL>::ComputeTransformedBasisFunctionDerivatives(
00325         const ChastePoint<ELEMENT_DIM>& rPoint,
00326         const c_matrix<double, ELEMENT_DIM, SPACE_DIM>& rInverseJacobian,
00327         c_matrix<double, SPACE_DIM, ELEMENT_DIM+1>& rReturnValue)
00328 {
00329     assert(ELEMENT_DIM < 4 && ELEMENT_DIM > 0);
00330     static c_matrix<double, ELEMENT_DIM, ELEMENT_DIM+1> grad_phi;
00331 
00332     LinearBasisFunction<ELEMENT_DIM>::ComputeBasisFunctionDerivatives(rPoint, grad_phi);
00333     rReturnValue = prod(trans(rInverseJacobian), grad_phi);
00334 }
00335 
00336 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
00337 void AbstractFeVolumeIntegralAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL>::AssembleOnElement(
00338     Element<ELEMENT_DIM,SPACE_DIM>& rElement,
00339     c_matrix<double, PROBLEM_DIM*(ELEMENT_DIM+1), PROBLEM_DIM*(ELEMENT_DIM+1) >& rAElem,
00340     c_vector<double, PROBLEM_DIM*(ELEMENT_DIM+1)>& rBElem)
00341 {
00347     c_matrix<double, SPACE_DIM, ELEMENT_DIM> jacobian;
00348     c_matrix<double, ELEMENT_DIM, SPACE_DIM> inverse_jacobian;
00349     double jacobian_determinant;
00350 
00351     mpMesh->GetInverseJacobianForElement(rElement.GetIndex(), jacobian, jacobian_determinant, inverse_jacobian);
00352 
00353     if (this->mAssembleMatrix)
00354     {
00355         rAElem.clear();
00356     }
00357 
00358     if (this->mAssembleVector)
00359     {
00360         rBElem.clear();
00361     }
00362 
00363     const unsigned num_nodes = rElement.GetNumNodes();
00364 
00365     // Allocate memory for the basis functions values and derivative values
00366     c_vector<double, ELEMENT_DIM+1> phi;
00367     c_matrix<double, SPACE_DIM, ELEMENT_DIM+1> grad_phi;
00368 
00369     // Loop over Gauss points
00370     for (unsigned quad_index=0; quad_index < mpQuadRule->GetNumQuadPoints(); quad_index++)
00371     {
00372         const ChastePoint<ELEMENT_DIM>& quad_point = mpQuadRule->rGetQuadPoint(quad_index);
00373 
00374         BasisFunction::ComputeBasisFunctions(quad_point, phi);
00375 
00376         if ( this->mAssembleMatrix || INTERPOLATION_LEVEL==NONLINEAR )
00377         {
00378             ComputeTransformedBasisFunctionDerivatives(quad_point, inverse_jacobian, grad_phi);
00379         }
00380 
00381         // Location of the Gauss point in the original element will be stored in x
00382         // Where applicable, u will be set to the value of the current solution at x
00383         ChastePoint<SPACE_DIM> x(0,0,0);
00384 
00385         c_vector<double,PROBLEM_DIM> u = zero_vector<double>(PROBLEM_DIM);
00386         c_matrix<double,PROBLEM_DIM,SPACE_DIM> grad_u = zero_matrix<double>(PROBLEM_DIM,SPACE_DIM);
00387 
00388         // Allow the concrete version of the assembler to interpolate any desired quantities
00389         this->ResetInterpolatedQuantities();
00390 
00391         // Interpolation
00392         for (unsigned i=0; i<num_nodes; i++)
00393         {
00394             const Node<SPACE_DIM>* p_node = rElement.GetNode(i);
00395 
00396             if (INTERPOLATION_LEVEL != CARDIAC) // don't even interpolate X if cardiac problem
00397             {
00398                 const c_vector<double, SPACE_DIM>& r_node_loc = p_node->rGetLocation();
00399                 // interpolate x
00400                 x.rGetLocation() += phi(i)*r_node_loc;
00401             }
00402 
00403             // Interpolate u and grad u if a current solution or guess exists
00404             unsigned node_global_index = rElement.GetNodeGlobalIndex(i);
00405             if (this->mCurrentSolutionOrGuessReplicated.GetSize() > 0)
00406             {
00407                 for (unsigned index_of_unknown=0; index_of_unknown<(INTERPOLATION_LEVEL!=CARDIAC ? PROBLEM_DIM : 1); index_of_unknown++)
00408                 {
00409                     /*
00410                      * If we have a solution (e.g. this is a dynamic problem) then
00411                      * interpolate the value at this quadrature point.
00412                      *
00413                      * NOTE: the following assumes that if, say, there are two unknowns
00414                      * u and v, they are stored in the current solution vector as
00415                      * [U1 V1 U2 V2 ... U_n V_n].
00416                      */
00417                     double u_at_node = this->GetCurrentSolutionOrGuessValue(node_global_index, index_of_unknown);
00418                     u(index_of_unknown) += phi(i)*u_at_node;
00419 
00420                     if (INTERPOLATION_LEVEL==NONLINEAR) // don't need to construct grad_phi or grad_u in other cases
00421                     {
00422                         for (unsigned j=0; j<SPACE_DIM; j++)
00423                         {
00424                             grad_u(index_of_unknown,j) += grad_phi(j,i)*u_at_node;
00425                         }
00426                     }
00427                 }
00428             }
00429 
00430             // Allow the concrete version of the assembler to interpolate any desired quantities
00431             this->IncrementInterpolatedQuantities(phi(i), p_node);
00432         }
00433 
00434         double wJ = jacobian_determinant * mpQuadRule->GetWeight(quad_index);
00435 
00436         // Create rAElem and rBElem
00437         if (this->mAssembleMatrix)
00438         {
00439             noalias(rAElem) += ComputeMatrixTerm(phi, grad_phi, x, u, grad_u, &rElement) * wJ;
00440         }
00441 
00442         if (this->mAssembleVector)
00443         {
00444             noalias(rBElem) += ComputeVectorTerm(phi, grad_phi, x, u, grad_u, &rElement) * wJ;
00445         }
00446     }
00447 }
00448 
00449 
00450 #endif /*ABSTRACTFEVOLUMEINTEGRALASSEMBLER_HPP_*/
Generated on Thu Dec 22 13:00:17 2011 for Chaste by  doxygen 1.6.3