AbstractStaticAssembler.hpp

00001 /*
00002 
00003 Copyright (C) University of Oxford, 2005-2010
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 #ifndef _ABSTRACTSTATICASSEMBLER_HPP_
00029 #define _ABSTRACTSTATICASSEMBLER_HPP_
00030 
00031 #include "AbstractAssembler.hpp"
00032 #include "LinearBasisFunction.hpp"
00033 #include "GaussianQuadratureRule.hpp"
00034 #include "AbstractTetrahedralMesh.hpp"
00035 #include "BoundaryConditionsContainer.hpp"
00036 #include "LinearSystem.hpp"
00037 #include "GaussianQuadratureRule.hpp"
00038 #include "ReplicatableVector.hpp"
00039 #include "DistributedVector.hpp"
00040 #include "PetscTools.hpp"
00041 #include "HeartEventHandler.hpp"
00042 #include <iostream>
00043 
00044 #include <boost/mpl/void.hpp>
00045 
00062 template<class T>
00063 struct AssemblerTraits
00064 {
00066     typedef T CVT_CLASS;
00068     typedef T CMT_CLASS;
00070     typedef AbstractAssembler<T::E_DIM, T::S_DIM, T::P_DIM> INTERPOLATE_CLASS;
00071 };
00072 
00074 template<>
00075 struct AssemblerTraits<boost::mpl::void_>
00076 {
00078     typedef boost::mpl::void_ CVT_CLASS;
00080     typedef boost::mpl::void_ CMT_CLASS;
00082     typedef boost::mpl::void_ INTERPOLATE_CLASS;
00083 };
00084 
00103 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool NON_HEART, class CONCRETE>
00104 class AbstractStaticAssembler : virtual public AbstractAssembler<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>
00105 {
00106 protected:
00107 
00109     AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>* mpMesh;
00110 
00112     GaussianQuadratureRule<ELEMENT_DIM>* mpQuadRule;
00113 
00115     GaussianQuadratureRule<ELEMENT_DIM-1>* mpSurfaceQuadRule;
00116 
00118     typedef LinearBasisFunction<ELEMENT_DIM> BasisFunction;
00120     typedef LinearBasisFunction<ELEMENT_DIM-1> SurfaceBasisFunction;
00121 
00126     ReplicatableVector mCurrentSolutionOrGuessReplicated;
00127 
00132     LinearSystem* mpLinearSystem;
00133 
00153     void ComputeTransformedBasisFunctionDerivatives(const ChastePoint<ELEMENT_DIM>& rPoint,
00154                                                     const c_matrix<double, ELEMENT_DIM, SPACE_DIM>& rInverseJacobian,
00155                                                     c_matrix<double, SPACE_DIM, ELEMENT_DIM+1>& rReturnValue);
00174     virtual void AssembleOnElement(Element<ELEMENT_DIM,SPACE_DIM>& rElement,
00175                                    c_matrix<double, PROBLEM_DIM*(ELEMENT_DIM+1), PROBLEM_DIM*(ELEMENT_DIM+1) >& rAElem,
00176                                    c_vector<double, PROBLEM_DIM*(ELEMENT_DIM+1)>& rBElem,
00177                                    bool assembleVector,
00178                                    bool assembleMatrix);
00179 
00189     virtual void AssembleOnSurfaceElement(const BoundaryElement<ELEMENT_DIM-1,SPACE_DIM>& rSurfaceElement,
00190                                           c_vector<double, PROBLEM_DIM*ELEMENT_DIM>& rBSurfElem);
00191 
00219     virtual void AssembleSystem(bool assembleVector, bool assembleMatrix,
00220                                 Vec currentSolutionOrGuess=NULL, double currentTime=0.0);
00221 
00226     virtual void PrepareForSolve();
00227 
00228 public:
00233     LinearSystem** GetLinearSystem();
00234 
00235 protected:
00240     ReplicatableVector& rGetCurrentSolutionOrGuess();
00241 
00248     virtual double GetCurrentSolutionOrGuessValue(unsigned nodeIndex, unsigned indexOfUnknown);
00249 
00250 public:
00251 
00257     AbstractStaticAssembler(unsigned numQuadPoints=2);
00258 
00267     void SetNumberOfQuadraturePointsPerDimension(unsigned numQuadPoints);
00268 
00274     void SetMesh(AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh);
00275 
00279     virtual ~AbstractStaticAssembler();
00280 
00281 };
00282 
00283 
00285 // Implementation
00287 
00288 
00289 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool NON_HEART, class CONCRETE>
00290 void AbstractStaticAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE>::ComputeTransformedBasisFunctionDerivatives(
00291         const ChastePoint<ELEMENT_DIM>& rPoint,
00292         const c_matrix<double, ELEMENT_DIM, SPACE_DIM>& rInverseJacobian,
00293         c_matrix<double, SPACE_DIM, ELEMENT_DIM+1>& rReturnValue)
00294 {
00295     assert(ELEMENT_DIM < 4 && ELEMENT_DIM > 0);
00296     static c_matrix<double, ELEMENT_DIM, ELEMENT_DIM+1> grad_phi;
00297 
00298     LinearBasisFunction<ELEMENT_DIM>::ComputeBasisFunctionDerivatives(rPoint, grad_phi);
00299     rReturnValue = prod(trans(rInverseJacobian), grad_phi);
00300 }
00301 
00302 
00303 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool NON_HEART, class CONCRETE>
00304 void AbstractStaticAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE>::AssembleOnElement(Element<ELEMENT_DIM,SPACE_DIM>& rElement,
00305                                 c_matrix<double, PROBLEM_DIM*(ELEMENT_DIM+1), PROBLEM_DIM*(ELEMENT_DIM+1) >& rAElem,
00306                                 c_vector<double, PROBLEM_DIM*(ELEMENT_DIM+1)>& rBElem,
00307                                 bool assembleVector,
00308                                 bool assembleMatrix)
00309 {
00310     GaussianQuadratureRule<ELEMENT_DIM>& quad_rule =
00311         *(AbstractStaticAssembler<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM, NON_HEART, CONCRETE>::mpQuadRule);
00312 
00318     c_matrix<double, SPACE_DIM, ELEMENT_DIM> jacobian;
00319     c_matrix<double, ELEMENT_DIM, SPACE_DIM> inverse_jacobian;
00320     double jacobian_determinant;
00321 
00322     mpMesh->GetInverseJacobianForElement(rElement.GetIndex(), jacobian, jacobian_determinant, inverse_jacobian);
00323 
00324 // With the new signature of GetInverseJacobianForElement, inverse and jacobian are returned at the same time
00325 //        // Initialise element contributions to zero
00326 //        if ( assembleMatrix || this->ProblemIsNonlinear() ) // don't need to construct grad_phi or grad_u in that case
00327 //        {
00328 //            this->mpMesh->GetInverseJacobianForElement(rElement.GetIndex(), inverse_jacobian);
00329 //        }
00330 
00331     if (assembleMatrix)
00332     {
00333         rAElem.clear();
00334     }
00335 
00336     if (assembleVector)
00337     {
00338         rBElem.clear();
00339     }
00340 
00341     const unsigned num_nodes = rElement.GetNumNodes();
00342 
00343     // allocate memory for the basis functions values and derivative values
00344     c_vector<double, ELEMENT_DIM+1> phi;
00345     c_matrix<double, SPACE_DIM, ELEMENT_DIM+1> grad_phi;
00346 
00347     // loop over Gauss points
00348     for (unsigned quad_index=0; quad_index < quad_rule.GetNumQuadPoints(); quad_index++)
00349     {
00350         const ChastePoint<ELEMENT_DIM>& quad_point = quad_rule.rGetQuadPoint(quad_index);
00351 
00352         BasisFunction::ComputeBasisFunctions(quad_point, phi);
00353 
00354         if ( assembleMatrix || this->ProblemIsNonlinear() )
00355         {
00356             ComputeTransformedBasisFunctionDerivatives(quad_point, inverse_jacobian, grad_phi);
00357         }
00358 
00359         // Location of the gauss point in the original element will be stored in x
00360         // Where applicable, u will be set to the value of the current solution at x
00361         ChastePoint<SPACE_DIM> x(0,0,0);
00362 
00363         c_vector<double,PROBLEM_DIM> u = zero_vector<double>(PROBLEM_DIM);
00364         c_matrix<double,PROBLEM_DIM,SPACE_DIM> grad_u = zero_matrix<double>(PROBLEM_DIM,SPACE_DIM);
00365 
00366         // allow the concrete version of the assembler to interpolate any
00367         // desired quantities
00368         static_cast<typename AssemblerTraits<CONCRETE>::INTERPOLATE_CLASS *>(this)->ResetInterpolatedQuantities();
00369 
00371         // interpolation
00373         for (unsigned i=0; i<num_nodes; i++)
00374         {
00375             const Node<SPACE_DIM>* p_node = rElement.GetNode(i);
00376 
00377             if (NON_HEART)
00378             {
00379                 const c_vector<double, SPACE_DIM>& r_node_loc = p_node->rGetLocation();
00380                 // interpolate x
00381                 x.rGetLocation() += phi(i)*r_node_loc;
00382             }
00383 
00384             // interpolate u and grad u if a current solution or guess exists
00385             unsigned node_global_index = rElement.GetNodeGlobalIndex(i);
00386             if (mCurrentSolutionOrGuessReplicated.GetSize()>0)
00387             {
00388                 for (unsigned index_of_unknown=0; index_of_unknown<(NON_HEART ? PROBLEM_DIM : 1); index_of_unknown++)
00389                 {
00390                     // If we have a current solution (e.g. this is a dynamic problem)
00391                     // get the value in a usable form.rElement
00392 
00393                     // NOTE - currentSolutionOrGuess input is actually now redundant at this point -
00394 
00395                     // NOTE - following assumes that, if say there are two unknowns u and v, they
00396                     // are stored in the current solution vector as
00397                     // [U1 V1 U2 V2 ... U_n V_n]
00398                     double u_at_node=GetCurrentSolutionOrGuessValue(node_global_index, index_of_unknown);
00399                     u(index_of_unknown) += phi(i)*u_at_node;
00400 
00401                     if (this->ProblemIsNonlinear() ) // don't need to construct grad_phi or grad_u in that case
00402                     {
00403                         for (unsigned j=0; j<SPACE_DIM; j++)
00404                         {
00405                             grad_u(index_of_unknown,j) += grad_phi(j,i)*u_at_node;
00406                         }
00407                     }
00408                 }
00409             }
00410 
00411             // allow the concrete version of the assembler to interpolate any
00412             // desired quantities
00413             static_cast<typename AssemblerTraits<CONCRETE>::INTERPOLATE_CLASS *>(this)->IncrementInterpolatedQuantities(phi(i), p_node);
00414         }
00415 
00416         double wJ = jacobian_determinant * quad_rule.GetWeight(quad_index);
00417 
00419         // create rAElem and rBElem
00421         if (assembleMatrix)
00422         {
00423             noalias(rAElem) += static_cast<typename AssemblerTraits<CONCRETE>::CMT_CLASS *>(this)->ComputeMatrixTerm(phi, grad_phi, x, u, grad_u, &rElement) * wJ;
00424         }
00425 
00426         if (assembleVector)
00427         {
00428             noalias(rBElem) += static_cast<typename AssemblerTraits<CONCRETE>::CVT_CLASS *>(this)->ComputeVectorTerm(phi, grad_phi, x, u, grad_u, &rElement) * wJ;
00429         }
00430     }
00431 }
00432 
00433 
00434 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool NON_HEART, class CONCRETE>
00435 void AbstractStaticAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE>::AssembleOnSurfaceElement(const BoundaryElement<ELEMENT_DIM-1,SPACE_DIM>& rSurfaceElement,
00436                                       c_vector<double, PROBLEM_DIM*ELEMENT_DIM>& rBSurfElem)
00437 {
00438     GaussianQuadratureRule<ELEMENT_DIM-1>& quad_rule =
00439         *(AbstractStaticAssembler<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM, NON_HEART, CONCRETE>::mpSurfaceQuadRule);
00440 
00441     c_vector<double, SPACE_DIM> weighted_direction;
00442     double jacobian_determinant;
00443     mpMesh->GetWeightedDirectionForBoundaryElement(rSurfaceElement.GetIndex(), weighted_direction, jacobian_determinant);
00444 
00445     rBSurfElem.clear();
00446 
00447     // allocate memory for the basis function values
00448     c_vector<double, ELEMENT_DIM>  phi;
00449 
00450     // loop over Gauss points
00451     for (unsigned quad_index=0; quad_index<quad_rule.GetNumQuadPoints(); quad_index++)
00452     {
00453         const ChastePoint<ELEMENT_DIM-1>& quad_point = quad_rule.rGetQuadPoint(quad_index);
00454 
00455         SurfaceBasisFunction::ComputeBasisFunctions(quad_point, phi);
00456 
00458         // interpolation
00460 
00461         // Location of the gauss point in the original element will be
00462         // stored in x
00463         ChastePoint<SPACE_DIM> x(0,0,0);
00464 
00465         this->ResetInterpolatedQuantities();
00466         for (unsigned i=0; i<rSurfaceElement.GetNumNodes(); i++)
00467         {
00468             const c_vector<double, SPACE_DIM> node_loc = rSurfaceElement.GetNode(i)->rGetLocation();
00469             x.rGetLocation() += phi(i)*node_loc;
00470 
00471             // allow the concrete version of the assembler to interpolate any
00472             // desired quantities
00473             IncrementInterpolatedQuantities(phi(i), rSurfaceElement.GetNode(i));
00474 
00476         }
00477 
00478         double wJ = jacobian_determinant * quad_rule.GetWeight(quad_index);
00479 
00481         // create rAElem and rBElem
00484         noalias(rBSurfElem) += ComputeVectorSurfaceTerm(rSurfaceElement, phi, x) * wJ;
00485     }
00486 }
00487 
00488 
00489 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool NON_HEART, class CONCRETE>
00490 void AbstractStaticAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE>::AssembleSystem(bool assembleVector, bool assembleMatrix,
00491                             Vec currentSolutionOrGuess, double currentTime)
00492 {
00493     HeartEventHandler::EventType assemble_event;
00494     if (assembleMatrix)
00495     {
00496         assemble_event = HeartEventHandler::ASSEMBLE_SYSTEM;
00497     }
00498     else
00499     {
00500         assemble_event = HeartEventHandler::ASSEMBLE_RHS;
00501     }
00502 
00503     // Check we've actually been asked to do something!
00504     assert(assembleVector || assembleMatrix);
00505 
00506     // Check the linear system object has been set up correctly
00507     assert(mpLinearSystem != NULL);
00508     assert(mpLinearSystem->GetSize() == PROBLEM_DIM * this->mpMesh->GetNumNodes());
00509     assert(!assembleVector || mpLinearSystem->rGetRhsVector() != NULL);
00510     assert(!assembleMatrix || mpLinearSystem->rGetLhsMatrix() != NULL);
00511 
00512     // Replicate the current solution and store so can be used in
00513     // AssembleOnElement
00514     if (currentSolutionOrGuess != NULL)
00515     {
00516         HeartEventHandler::BeginEvent(HeartEventHandler::COMMUNICATION);
00517         this->mCurrentSolutionOrGuessReplicated.ReplicatePetscVector(currentSolutionOrGuess);
00518         HeartEventHandler::EndEvent(HeartEventHandler::COMMUNICATION);
00519     }
00520 
00521     // the AssembleOnElement type methods will determine if a current solution or
00522     // current guess exists by looking at the size of the replicated vector, so
00523     // check the size is zero if there isn't a current solution
00524     assert(    ( currentSolutionOrGuess && mCurrentSolutionOrGuessReplicated.GetSize()>0)
00525             || ( !currentSolutionOrGuess && mCurrentSolutionOrGuessReplicated.GetSize()==0));
00526 
00527     // the concrete class can override this following method if there is
00528     // work to be done before assembly
00529     this->PrepareForAssembleSystem(currentSolutionOrGuess, currentTime);
00530 
00531     // this has to be below PrepareForAssembleSystem as in that
00532     // method the odes are solved in cardiac problems
00533     HeartEventHandler::BeginEvent(assemble_event);
00534 
00535     // Zero the matrix/vector if it is to be assembled
00536     if (assembleVector)
00537     {
00538         mpLinearSystem->ZeroRhsVector();
00539     }
00540     if (assembleMatrix)
00541     {
00542         mpLinearSystem->ZeroLhsMatrix();
00543     }
00544 
00545     const size_t STENCIL_SIZE=PROBLEM_DIM*(ELEMENT_DIM+1);
00546     c_matrix<double, STENCIL_SIZE, STENCIL_SIZE> a_elem;
00547     c_vector<double, STENCIL_SIZE> b_elem;
00548 
00550     // loop over elements
00552     for (typename AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::ElementIterator iter = this->mpMesh->GetElementIteratorBegin();
00553          iter != this->mpMesh->GetElementIteratorEnd();
00554          ++iter)
00555     {
00556         Element<ELEMENT_DIM, SPACE_DIM>& element = *iter;
00557 
00558         if (element.GetOwnership() == true)
00559         {
00560             AssembleOnElement(element, a_elem, b_elem, assembleVector, assembleMatrix);
00561 
00562             unsigned p_indices[STENCIL_SIZE];
00563             element.GetStiffnessMatrixGlobalIndices(PROBLEM_DIM, p_indices);
00564 
00565             if (assembleMatrix)
00566             {
00567                 mpLinearSystem->AddLhsMultipleValues(p_indices, a_elem);
00568             }
00569 
00570             if (assembleVector)
00571             {
00572                 mpLinearSystem->AddRhsMultipleValues(p_indices, b_elem);
00573             }
00574         }
00575     }
00576 
00577     // add the integrals associated with Neumann boundary conditions to the linear system
00578     typename AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::BoundaryElementIterator
00579         surf_iter = this->mpMesh->GetBoundaryElementIteratorBegin();
00580 
00582     // Apply any Neumann boundary conditions
00583     //
00584     // NB. We assume that if an element has a boundary condition on any unknown there is a boundary condition
00585     // on unknown 0. This can be so for any problem by adding zero constant conditions where required
00586     // although this is a bit inefficient. Proper solution involves changing BCC to have a map of arrays
00587     // boundary conditions rather than an array of maps.
00589     if (assembleVector)
00590     {
00591         HeartEventHandler::EndEvent(assemble_event);
00592         this->ApplyNeummanBoundaryConditions();
00593         HeartEventHandler::BeginEvent(assemble_event);
00594     }
00595 
00596     if (assembleVector)
00597     {
00598         mpLinearSystem->AssembleRhsVector();
00599     }
00600 
00601     if (assembleMatrix)
00602     {
00603         mpLinearSystem->AssembleIntermediateLhsMatrix();
00604     }
00605 
00606     // Apply Dirichlet boundary conditions
00607     this->ApplyDirichletConditions(currentSolutionOrGuess, assembleMatrix);
00608 
00609     this->FinaliseLinearSystem(currentSolutionOrGuess, currentTime, assembleVector, assembleMatrix);
00610 
00611     if (assembleVector)
00612     {
00613         mpLinearSystem->AssembleRhsVector();
00614     }
00615     if (assembleMatrix)
00616     {
00617         mpLinearSystem->AssembleFinalLhsMatrix();
00618     }
00619 
00620     // overload this method if the assembler has to do anything else
00621     // required (like setting up a null basis (see BidomainDg0Assembler))
00622     this->FinaliseAssembleSystem(currentSolutionOrGuess, currentTime);
00623 
00624     HeartEventHandler::EndEvent(assemble_event);
00625 }
00626 
00627 
00628 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool NON_HEART, class CONCRETE>
00629 void AbstractStaticAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE>::PrepareForSolve()
00630 {
00631     assert(mpMesh != NULL);
00632     assert(this->mpBoundaryConditions != NULL);
00633 
00634     assert(mpMesh->GetNumNodes() == mpMesh->GetDistributedVectorFactory()->GetProblemSize());
00635 
00636     mpMesh->SetElementOwnerships(mpMesh->GetDistributedVectorFactory()->GetLow(),
00637                                  mpMesh->GetDistributedVectorFactory()->GetHigh());
00638 }
00639 
00640 
00641 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool NON_HEART, class CONCRETE>
00642 LinearSystem** AbstractStaticAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE>::GetLinearSystem()
00643 {
00644     return &mpLinearSystem;
00645 }
00646 
00647 
00648 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool NON_HEART, class CONCRETE>
00649 ReplicatableVector& AbstractStaticAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE>::rGetCurrentSolutionOrGuess()
00650 {
00651     return mCurrentSolutionOrGuessReplicated;
00652 }
00653 
00654 
00655 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool NON_HEART, class CONCRETE>
00656 double AbstractStaticAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE>::GetCurrentSolutionOrGuessValue(unsigned nodeIndex, unsigned indexOfUnknown)
00657 {
00658     return mCurrentSolutionOrGuessReplicated[ PROBLEM_DIM*nodeIndex + indexOfUnknown];
00659 }
00660 
00661 
00662 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool NON_HEART, class CONCRETE>
00663 AbstractStaticAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE>::AbstractStaticAssembler(unsigned numQuadPoints)
00664     : AbstractAssembler<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>()
00665 {
00666     // Initialise mesh and bcs to null, so we can check they
00667     // have been set before attempting to solve
00668     mpMesh = NULL;
00669 
00670     mpQuadRule = NULL;
00671     mpSurfaceQuadRule = NULL;
00672     SetNumberOfQuadraturePointsPerDimension(numQuadPoints);
00673 
00674     mpLinearSystem = NULL;
00675 }
00676 
00677 
00678 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool NON_HEART, class CONCRETE>
00679 void AbstractStaticAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE>::SetNumberOfQuadraturePointsPerDimension(unsigned numQuadPoints)
00680 {
00681     delete mpQuadRule;
00682     mpQuadRule = new GaussianQuadratureRule<ELEMENT_DIM>(numQuadPoints);
00683     delete mpSurfaceQuadRule;
00684     mpSurfaceQuadRule = new GaussianQuadratureRule<ELEMENT_DIM-1>(numQuadPoints);
00685 }
00686 
00687 
00688 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool NON_HEART, class CONCRETE>
00689 void AbstractStaticAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE>::SetMesh(AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh)
00690 {
00691     mpMesh = pMesh;
00692 }
00693 
00694 
00695 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool NON_HEART, class CONCRETE>
00696 AbstractStaticAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE>::~AbstractStaticAssembler()
00697 {
00698     delete mpQuadRule;
00699     delete mpSurfaceQuadRule;
00700     delete mpLinearSystem;
00701 }
00702 
00703 #endif //_ABSTRACTSTATICASSEMBLER_HPP_

Generated by  doxygen 1.6.2