AbstractFeCableIntegralAssembler.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 ABSTRACTFECABLEINTEGRALASSEMBLER_HPP_
00030 #define ABSTRACTFECABLEINTEGRALASSEMBLER_HPP_
00031 
00032 #include "AbstractFeAssemblerCommon.hpp"
00033 #include "GaussianQuadratureRule.hpp"
00034 #include "PetscVecTools.hpp"
00035 #include "PetscMatTools.hpp"
00036 
00043 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
00044 class AbstractFeCableIntegralAssembler : public AbstractFeAssemblerCommon<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM,CAN_ASSEMBLE_VECTOR,CAN_ASSEMBLE_MATRIX,INTERPOLATION_LEVEL>
00045 {
00046 protected:
00048     static const unsigned CABLE_ELEMENT_DIM = 1;
00049 
00051     static const unsigned NUM_CABLE_ELEMENT_NODES = 2;
00052 
00054     MixedDimensionMesh<ELEMENT_DIM, SPACE_DIM>* mpMesh;
00055 
00057     GaussianQuadratureRule<1>* mpCableQuadRule;
00058 
00060     typedef LinearBasisFunction<1> CableBasisFunction;
00061 
00076     void ComputeTransformedBasisFunctionDerivatives(const ChastePoint<CABLE_ELEMENT_DIM>& rPoint,
00077                                                     const c_matrix<double, CABLE_ELEMENT_DIM, SPACE_DIM>& rInverseJacobian,
00078                                                     c_matrix<double, SPACE_DIM, NUM_CABLE_ELEMENT_NODES>& rReturnValue);
00079 
00085     void DoAssemble();
00086 
00109     virtual c_matrix<double,PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES,PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES> ComputeCableMatrixTerm(
00110         c_vector<double, NUM_CABLE_ELEMENT_NODES>& rPhi,
00111         c_matrix<double, SPACE_DIM, NUM_CABLE_ELEMENT_NODES>& rGradPhi,
00112         ChastePoint<SPACE_DIM>& rX,
00113         c_vector<double,PROBLEM_DIM>& rU,
00114         c_matrix<double, PROBLEM_DIM, SPACE_DIM>& rGradU,
00115         Element<CABLE_ELEMENT_DIM,SPACE_DIM>* pElement)
00116     {
00117         // If this line is reached this means this method probably hasn't been over-ridden correctly in
00118         // the concrete class
00119         NEVER_REACHED;
00120         return zero_matrix<double>(PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES,PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES);
00121     }
00122 
00143     virtual c_vector<double,PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES> ComputeCableVectorTerm(
00144         c_vector<double, NUM_CABLE_ELEMENT_NODES>& rPhi,
00145         c_matrix<double, SPACE_DIM, NUM_CABLE_ELEMENT_NODES>& rGradPhi,
00146         ChastePoint<SPACE_DIM>& rX,
00147         c_vector<double,PROBLEM_DIM>& rU,
00148         c_matrix<double, PROBLEM_DIM, SPACE_DIM>& rGradU,
00149         Element<CABLE_ELEMENT_DIM,SPACE_DIM>* pElement)
00150     {
00151         // If this line is reached this means this method probably hasn't been over-ridden correctly in
00152         // the concrete class
00153         NEVER_REACHED;
00154         return zero_vector<double>(PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES);
00155     }
00156 
00171     virtual void AssembleOnCableElement(Element<CABLE_ELEMENT_DIM,SPACE_DIM>& rElement,
00172                                         c_matrix<double, PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES, PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES >& rAElem,
00173                                         c_vector<double, PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES>& rBElem);
00174 
00175 
00183     virtual bool ElementAssemblyCriterion(Element<CABLE_ELEMENT_DIM,SPACE_DIM>& rElement)
00184     {
00185         return true;
00186     }
00187 
00188 public:
00189 
00197     AbstractFeCableIntegralAssembler(MixedDimensionMesh<ELEMENT_DIM,SPACE_DIM>* pMesh, unsigned numQuadPoints=2);
00198 
00202     virtual ~AbstractFeCableIntegralAssembler()
00203     {
00204         delete mpCableQuadRule;
00205     }
00206 };
00207 
00208 
00209 
00210 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
00211 AbstractFeCableIntegralAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL>::AbstractFeCableIntegralAssembler(
00212             MixedDimensionMesh<ELEMENT_DIM,SPACE_DIM>* pMesh, unsigned numQuadPoints)
00213     : AbstractFeAssemblerCommon<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL>(),
00214       mpMesh(pMesh)
00215 {
00216     assert(pMesh);
00217     assert(numQuadPoints > 0);
00218     assert(CAN_ASSEMBLE_VECTOR || CAN_ASSEMBLE_MATRIX);
00219 
00220     mpCableQuadRule = new GaussianQuadratureRule<CABLE_ELEMENT_DIM>(numQuadPoints);
00221 
00222     // Not supporting this yet - if a nonlinear assembler on cable elements is required, uncomment code
00223     // in AssembleOnCableElement below (search for NONLINEAR)
00224     assert(INTERPOLATION_LEVEL!=NONLINEAR);
00225 }
00226 
00227 
00228 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
00229 void AbstractFeCableIntegralAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL>::DoAssemble()
00230 {
00231     HeartEventHandler::EventType assemble_event;
00232     if (this->mAssembleMatrix)
00233     {
00234         assemble_event = HeartEventHandler::ASSEMBLE_SYSTEM;
00235     }
00236     else
00237     {
00238         assemble_event = HeartEventHandler::ASSEMBLE_RHS;
00239     }
00240 
00241     if (this->mAssembleMatrix && this->mMatrixToAssemble==NULL)
00242     {
00243         EXCEPTION("Matrix to be assembled has not been set");
00244     }
00245     if (this->mAssembleVector && this->mVectorToAssemble==NULL)
00246     {
00247         EXCEPTION("Vector to be assembled has not been set");
00248     }
00249 
00250     // This has to be below PrepareForAssembleSystem as in that method the ODEs are solved in cardiac problems
00251     HeartEventHandler::BeginEvent(assemble_event);
00252 
00253     // Zero the matrix/vector if it is to be assembled
00254     if (this->mAssembleVector && this->mZeroVectorBeforeAssembly)
00255     {
00256         PetscVecTools::Zero(this->mVectorToAssemble);
00257     }
00258     if (this->mAssembleMatrix && this->mZeroMatrixBeforeAssembly)
00259     {
00260         PetscMatTools::Zero(this->mMatrixToAssemble);
00261     }
00262 
00263     const size_t STENCIL_SIZE=PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES;
00264     c_matrix<double, STENCIL_SIZE, STENCIL_SIZE> a_elem;
00265     c_vector<double, STENCIL_SIZE> b_elem;
00266 
00267     // Loop over elements
00268     if (this->mAssembleMatrix || this->mAssembleVector)
00269     {
00270         for (typename MixedDimensionMesh<CABLE_ELEMENT_DIM, SPACE_DIM>::CableElementIterator iter = mpMesh->GetCableElementIteratorBegin();
00271              iter != mpMesh->GetCableElementIteratorEnd();
00272              ++iter)
00273         {
00274             Element<CABLE_ELEMENT_DIM, SPACE_DIM>& r_element = *(*iter);
00275 
00276             // Test for ownership first, since it's pointless to test the criterion on something which we might know nothing about.
00277             if ( r_element.GetOwnership() == true && ElementAssemblyCriterion(r_element)==true )
00278             {
00279                 AssembleOnCableElement(r_element, a_elem, b_elem);
00280 
00281                 unsigned p_indices[STENCIL_SIZE];
00282                 r_element.GetStiffnessMatrixGlobalIndices(PROBLEM_DIM, p_indices);
00283 
00284                 if (this->mAssembleMatrix)
00285                 {
00286                     PetscMatTools::AddMultipleValues<STENCIL_SIZE>(this->mMatrixToAssemble, p_indices, a_elem);
00287                 }
00288 
00289                 if (this->mAssembleVector)
00290                 {
00291                     PetscVecTools::AddMultipleValues<STENCIL_SIZE>(this->mVectorToAssemble, p_indices, b_elem);
00292                 }
00293             }
00294         }
00295     }
00296 
00297     HeartEventHandler::EndEvent(assemble_event);
00298 }
00299 
00301 // Implementation - AssembleOnCableElement and smaller
00303 
00304 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
00305 void AbstractFeCableIntegralAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL>::ComputeTransformedBasisFunctionDerivatives(
00306         const ChastePoint<CABLE_ELEMENT_DIM>& rPoint,
00307         const c_matrix<double, CABLE_ELEMENT_DIM, SPACE_DIM>& rInverseJacobian,
00308         c_matrix<double, SPACE_DIM, NUM_CABLE_ELEMENT_NODES>& rReturnValue)
00309 {
00310     static c_matrix<double, CABLE_ELEMENT_DIM, NUM_CABLE_ELEMENT_NODES> grad_phi;
00311 
00312     LinearBasisFunction<CABLE_ELEMENT_DIM>::ComputeBasisFunctionDerivatives(rPoint, grad_phi);
00313     rReturnValue = prod(trans(rInverseJacobian), grad_phi);
00314 }
00315 
00316 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
00317 void AbstractFeCableIntegralAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL>::AssembleOnCableElement(
00318     Element<CABLE_ELEMENT_DIM,SPACE_DIM>& rElement,
00319     c_matrix<double, PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES, PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES >& rAElem,
00320     c_vector<double, PROBLEM_DIM*NUM_CABLE_ELEMENT_NODES>& rBElem)
00321 {
00327     c_matrix<double, SPACE_DIM, CABLE_ELEMENT_DIM> jacobian;
00328     c_matrix<double, CABLE_ELEMENT_DIM, SPACE_DIM> inverse_jacobian;
00329     double jacobian_determinant;
00330 
00333     // mpMesh->GetInverseJacobianForElement(rElement.GetIndex(), jacobian, jacobian_determinant, inverse_jacobian);
00335     rElement.CalculateInverseJacobian(jacobian, jacobian_determinant, inverse_jacobian);
00336 
00337     if (this->mAssembleMatrix)
00338     {
00339         rAElem.clear();
00340     }
00341 
00342     if (this->mAssembleVector)
00343     {
00344         rBElem.clear();
00345     }
00346 
00347     const unsigned num_nodes = rElement.GetNumNodes();
00348 
00349     // Allocate memory for the basis functions values and derivative values
00350     c_vector<double, NUM_CABLE_ELEMENT_NODES> phi;
00351     c_matrix<double, SPACE_DIM, NUM_CABLE_ELEMENT_NODES> grad_phi;
00352 
00353     // Loop over Gauss points
00354     for (unsigned quad_index=0; quad_index < mpCableQuadRule->GetNumQuadPoints(); quad_index++)
00355     {
00356         const ChastePoint<CABLE_ELEMENT_DIM>& quad_point = mpCableQuadRule->rGetQuadPoint(quad_index);
00357 
00358         CableBasisFunction::ComputeBasisFunctions(quad_point, phi);
00359 
00360         if ( this->mAssembleMatrix || INTERPOLATION_LEVEL==NONLINEAR )
00361         {
00362             ComputeTransformedBasisFunctionDerivatives(quad_point, inverse_jacobian, grad_phi);
00363         }
00364 
00365         // Location of the Gauss point in the original element will be stored in x
00366         // Where applicable, u will be set to the value of the current solution at x
00367         ChastePoint<SPACE_DIM> x(0,0,0);
00368 
00369         c_vector<double,PROBLEM_DIM> u = zero_vector<double>(PROBLEM_DIM);
00370         c_matrix<double,PROBLEM_DIM,SPACE_DIM> grad_u = zero_matrix<double>(PROBLEM_DIM,SPACE_DIM);
00371 
00372         // Allow the concrete version of the assembler to interpolate any desired quantities
00373         this->ResetInterpolatedQuantities();
00374 
00375         // Interpolation
00376         for (unsigned i=0; i<num_nodes; i++)
00377         {
00378             const Node<SPACE_DIM>* p_node = rElement.GetNode(i);
00379 
00380             if (INTERPOLATION_LEVEL != CARDIAC) // don't interpolate X if cardiac problem
00381             {
00382                 const c_vector<double, SPACE_DIM>& r_node_loc = p_node->rGetLocation();
00383                 // interpolate x
00384                 x.rGetLocation() += phi(i)*r_node_loc;
00385             }
00386 
00387             // Interpolate u and grad u if a current solution or guess exists
00388             unsigned node_global_index = rElement.GetNodeGlobalIndex(i);
00389             if (this->mCurrentSolutionOrGuessReplicated.GetSize() > 0)
00390             {
00391                 for (unsigned index_of_unknown=0; index_of_unknown<(INTERPOLATION_LEVEL!=CARDIAC ? PROBLEM_DIM : 1); index_of_unknown++)
00392                 {
00393                     /*
00394                      * If we have a solution (e.g. this is a dynamic problem) then
00395                      * interpolate the value at this quadrature point.
00396                      *
00397                      * NOTE: the following assumes that if, say, there are two unknowns
00398                      * u and v, they are stored in the current solution vector as
00399                      * [U1 V1 U2 V2 ... U_n V_n].
00400                      */
00401                     double u_at_node = this->GetCurrentSolutionOrGuessValue(node_global_index, index_of_unknown);
00402                     u(index_of_unknown) += phi(i)*u_at_node;
00403 
00405 //                    if (INTERPOLATION_LEVEL==NONLINEAR) // don't need to construct grad_phi or grad_u in other cases
00406 //                    {
00407 //                        for (unsigned j=0; j<SPACE_DIM; j++)
00408 //                        {
00409 //                            grad_u(index_of_unknown,j) += grad_phi(j,i)*u_at_node;
00410 //                        }
00411 //                    }
00412                 }
00413             }
00414 
00415             // Allow the concrete version of the assembler to interpolate any desired quantities
00416             IncrementInterpolatedQuantities(phi(i), p_node);
00417         }
00418 
00419         double wJ = jacobian_determinant * mpCableQuadRule->GetWeight(quad_index);
00420 
00421         // Create rAElem and rBElem
00422         if (this->mAssembleMatrix)
00423         {
00424             noalias(rAElem) += ComputeCableMatrixTerm(phi, grad_phi, x, u, grad_u, &rElement) * wJ;
00425         }
00426 
00427         if (this->mAssembleVector)
00428         {
00429             noalias(rBElem) += ComputeCableVectorTerm(phi, grad_phi, x, u, grad_u, &rElement) * wJ;
00430         }
00431     }
00432 }
00433 
00434 #endif /*ABSTRACTFECABLEINTEGRALASSEMBLER_HPP_*/
Generated on Thu Dec 22 13:00:17 2011 for Chaste by  doxygen 1.6.3