AbstractFeObjectAssembler.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 
00029 #ifndef ABSTRACTFEOBJECTASSEMBLER_HPP_
00030 #define ABSTRACTFEOBJECTASSEMBLER_HPP_
00031 
00032 #include "LinearBasisFunction.hpp"
00033 #include "GaussianQuadratureRule.hpp"
00034 #include "BoundaryConditionsContainer.hpp"
00035 #include "ReplicatableVector.hpp"
00036 #include "DistributedVector.hpp"
00037 #include "HeartEventHandler.hpp"
00038 #include "AbstractTetrahedralMesh.hpp"
00039 #include "PetscTools.hpp"
00040 
00042 
00052 typedef enum InterpolationLevel_
00053 {
00054     CARDIAC = 0,
00055     NORMAL,
00056     NONLINEAR
00057 } InterpolationLevel;
00058 
00059 
00088 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
00089 class AbstractFeObjectAssembler : boost::noncopyable
00090 {
00091 protected:
00092 
00094     Vec mVectorToAssemble;
00095 
00097     Mat mMatrixToAssemble;
00098 
00103     bool mAssembleMatrix;
00104 
00106     bool mAssembleVector;
00107     
00109     bool mZeroVectorBeforeAssembly;
00110 
00112     bool mApplyNeummanBoundaryConditionsToVector;
00113     
00115     bool mOnlyAssembleOnSurfaceElements;
00116     
00118     PetscInt mOwnershipRangeLo;
00120     PetscInt mOwnershipRangeHi;
00121 
00123     GaussianQuadratureRule<ELEMENT_DIM>* mpQuadRule;
00124 
00126     GaussianQuadratureRule<ELEMENT_DIM-1>* mpSurfaceQuadRule;
00127 
00129     typedef LinearBasisFunction<ELEMENT_DIM> BasisFunction;
00131     typedef LinearBasisFunction<ELEMENT_DIM-1> SurfaceBasisFunction;
00132 
00133 
00138     ReplicatableVector mCurrentSolutionOrGuessReplicated;
00139 
00140 
00160     void ComputeTransformedBasisFunctionDerivatives(const ChastePoint<ELEMENT_DIM>& rPoint,
00161                                                     const c_matrix<double, ELEMENT_DIM, SPACE_DIM>& rInverseJacobian,
00162                                                     c_matrix<double, SPACE_DIM, ELEMENT_DIM+1>& rReturnValue);
00163 
00164 
00170     void DoAssemble();
00171 
00178     void ApplyNeummanBoundaryConditionsToVector();
00179     
00185     virtual double GetCurrentSolutionOrGuessValue(unsigned nodeIndex, unsigned indexOfUnknown)
00186     {
00187         return mCurrentSolutionOrGuessReplicated[ PROBLEM_DIM*nodeIndex + indexOfUnknown];
00188     }
00189 
00199     template<size_t SMALL_MATRIX_SIZE>
00200     void AddMultipleValuesToMatrix(unsigned* matrixRowAndColIndices, c_matrix<double, SMALL_MATRIX_SIZE, SMALL_MATRIX_SIZE>& rSmallMatrix);
00201 
00211     template<size_t SMALL_VECTOR_SIZE>
00212     void AddMultipleValuesToVector(unsigned* vectorIndices, c_vector<double, SMALL_VECTOR_SIZE>& smallVector);
00213     
00214 protected:
00216     AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>* mpMesh;
00217 
00222     BoundaryConditionsContainer<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>* mpBoundaryConditions;
00223 
00229     virtual void ResetInterpolatedQuantities()
00230     {}
00231 
00240     virtual void IncrementInterpolatedQuantities(double phiI, const Node<SPACE_DIM>* pNode)
00241     {}    
00242 
00243 
00265     virtual c_matrix<double,PROBLEM_DIM*(ELEMENT_DIM+1),PROBLEM_DIM*(ELEMENT_DIM+1)> ComputeMatrixTerm(
00266         c_vector<double, ELEMENT_DIM+1>& rPhi,
00267         c_matrix<double, SPACE_DIM, ELEMENT_DIM+1>& rGradPhi,
00268         ChastePoint<SPACE_DIM>& rX,
00269         c_vector<double,PROBLEM_DIM>& rU,
00270         c_matrix<double, PROBLEM_DIM, SPACE_DIM>& rGradU,
00271         Element<ELEMENT_DIM,SPACE_DIM>* pElement)
00272     {
00273         NEVER_REACHED;
00274         return zero_matrix<double>(PROBLEM_DIM*(ELEMENT_DIM+1),PROBLEM_DIM*(ELEMENT_DIM+1));
00275     }
00276 
00298     virtual c_vector<double,PROBLEM_DIM*(ELEMENT_DIM+1)> ComputeVectorTerm(
00299         c_vector<double, ELEMENT_DIM+1>& rPhi,
00300         c_matrix<double, SPACE_DIM, ELEMENT_DIM+1>& rGradPhi,
00301         ChastePoint<SPACE_DIM>& rX,
00302         c_vector<double,PROBLEM_DIM>& rU,
00303         c_matrix<double, PROBLEM_DIM, SPACE_DIM>& rGradU,
00304         Element<ELEMENT_DIM,SPACE_DIM>* pElement)
00305     {
00306         NEVER_REACHED;
00307         return zero_vector<double>(PROBLEM_DIM*(ELEMENT_DIM+1));
00308     }
00309 
00327     virtual c_vector<double, PROBLEM_DIM*ELEMENT_DIM> ComputeVectorSurfaceTerm(
00328         const BoundaryElement<ELEMENT_DIM-1,SPACE_DIM>& rSurfaceElement,
00329         c_vector<double, ELEMENT_DIM>& rPhi,
00330         ChastePoint<SPACE_DIM>& rX)
00331     {
00332         NEVER_REACHED;
00333         return zero_vector<double>(PROBLEM_DIM*ELEMENT_DIM);
00334     }
00335     
00342     virtual bool ElementAssemblyCriterion(Element<ELEMENT_DIM,SPACE_DIM>& rElement)
00343     {
00344         return true;
00345     }
00346 
00361     virtual void AssembleOnElement(Element<ELEMENT_DIM,SPACE_DIM>& rElement,
00362                                    c_matrix<double, PROBLEM_DIM*(ELEMENT_DIM+1), PROBLEM_DIM*(ELEMENT_DIM+1) >& rAElem,
00363                                    c_vector<double, PROBLEM_DIM*(ELEMENT_DIM+1)>& rBElem);
00364 
00374     virtual void AssembleOnSurfaceElement(const BoundaryElement<ELEMENT_DIM-1,SPACE_DIM>& rSurfaceElement,
00375                                           c_vector<double, PROBLEM_DIM*ELEMENT_DIM>& rBSurfElem);
00376 
00377 
00378 public:
00385     AbstractFeObjectAssembler(AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh, unsigned numQuadPoints=2);
00386 
00392     void SetApplyNeummanBoundaryConditionsToVector(BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>* pBcc);
00393 
00402     void OnlyAssembleOnSurfaceElements(bool onlyAssembleOnSurfaceElements = true)
00403     {
00404         assert(CAN_ASSEMBLE_VECTOR);
00405         mOnlyAssembleOnSurfaceElements = onlyAssembleOnSurfaceElements;
00406     }
00407     
00412     void SetMatrixToAssemble(Mat& rMatToAssemble);
00413 
00420     void SetVectorToAssemble(Vec& rVecToAssemble, bool zeroVectorBeforeAssembly);
00421 
00427     void SetCurrentSolution(Vec currentSolution);
00428 
00429 
00433     void Assemble()
00434     {
00435         mAssembleMatrix = CAN_ASSEMBLE_MATRIX;
00436         mAssembleVector = CAN_ASSEMBLE_VECTOR;
00437         DoAssemble();
00438     }
00439         
00443     void AssembleMatrix()
00444     {
00445         assert(CAN_ASSEMBLE_MATRIX);
00446         mAssembleMatrix = true;
00447         mAssembleVector = false;
00448         DoAssemble();
00449     }
00450 
00454     void AssembleVector()
00455     {
00456         assert(CAN_ASSEMBLE_VECTOR);
00457         mAssembleMatrix = false;
00458         mAssembleVector = true;
00459         DoAssemble();
00460     }
00461 
00465     virtual ~AbstractFeObjectAssembler()
00466     {
00467         delete mpQuadRule;
00468         if(mpSurfaceQuadRule)
00469         {
00470             delete mpSurfaceQuadRule;
00471         }
00472     }
00473 };
00474 
00475 
00476 
00477 
00478 
00479 
00480 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL> // class CONCRETE>
00481 AbstractFeObjectAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL>::AbstractFeObjectAssembler(AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh, unsigned numQuadPoints)
00482     : mVectorToAssemble(NULL),
00483       mMatrixToAssemble(NULL),
00484       mZeroVectorBeforeAssembly(true),
00485       mApplyNeummanBoundaryConditionsToVector(false),
00486       mOnlyAssembleOnSurfaceElements(false),
00487       mpMesh(pMesh),
00488       mpBoundaryConditions(NULL)
00489 {
00490     assert(pMesh);
00491     assert(numQuadPoints > 0);
00492     assert(CAN_ASSEMBLE_VECTOR || CAN_ASSEMBLE_MATRIX);
00493 
00494     mpQuadRule = new GaussianQuadratureRule<ELEMENT_DIM>(numQuadPoints);
00495 
00496     if(CAN_ASSEMBLE_VECTOR)
00497     {
00498         mpSurfaceQuadRule = new GaussianQuadratureRule<ELEMENT_DIM-1>(numQuadPoints);
00499     }
00500     else
00501     {
00502         mpSurfaceQuadRule = NULL;
00503     }
00504 }  
00505 
00506 
00507 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL> //, class CONCRETE>
00508 void AbstractFeObjectAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL>::SetMatrixToAssemble(Mat& rMatToAssemble)
00509 {
00510     assert(rMatToAssemble);
00511     MatGetOwnershipRange(rMatToAssemble, &mOwnershipRangeLo, &mOwnershipRangeHi);
00512 
00513     mMatrixToAssemble = rMatToAssemble;
00514 }
00515 
00516 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL> //, class CONCRETE>
00517 void AbstractFeObjectAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL>::SetVectorToAssemble(Vec& rVecToAssemble, bool zeroVectorBeforeAssembly)
00518 {
00519     assert(rVecToAssemble);
00520     VecGetOwnershipRange(rVecToAssemble, &mOwnershipRangeLo, &mOwnershipRangeHi);
00521 
00522     mVectorToAssemble = rVecToAssemble;
00523     mZeroVectorBeforeAssembly = zeroVectorBeforeAssembly;
00524 }
00525 
00526 
00527 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL> //, class CONCRETE>
00528 void AbstractFeObjectAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL>::SetCurrentSolution(Vec currentSolution)
00529 {
00530     assert(currentSolution != NULL);
00531     // Replicate the current solution and store so can be used in
00532     // AssembleOnElement
00533     
00534     HeartEventHandler::BeginEvent(HeartEventHandler::COMMUNICATION);
00535     mCurrentSolutionOrGuessReplicated.ReplicatePetscVector(currentSolution);
00536     HeartEventHandler::EndEvent(HeartEventHandler::COMMUNICATION);
00537     
00538     // the AssembleOnElement type methods will determine if a current solution or
00539     // current guess exists by looking at the size of the replicated vector, so
00540     // check the size is zero if there isn't a current solution
00541     assert( mCurrentSolutionOrGuessReplicated.GetSize()>0 );
00542 }
00543 
00544 
00545 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL> //, class CONCRETE>
00546 void AbstractFeObjectAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL>::DoAssemble()
00547 {
00548     HeartEventHandler::EventType assemble_event;
00549     if (mAssembleMatrix)
00550     {
00551         assemble_event = HeartEventHandler::ASSEMBLE_SYSTEM;
00552     }
00553     else
00554     {
00555         assemble_event = HeartEventHandler::ASSEMBLE_RHS;
00556     }
00557 
00558     if(mAssembleMatrix && mMatrixToAssemble==NULL)
00559     {
00560         EXCEPTION("Matrix to be assembled has not been set");
00561     }
00562     if(mAssembleVector && mVectorToAssemble==NULL)
00563     {
00564         EXCEPTION("Vector to be assembled has not been set");
00565     }
00566 
00567     // this has to be below PrepareForAssembleSystem as in that
00568     // method the odes are solved in cardiac problems
00569     HeartEventHandler::BeginEvent(assemble_event);
00570 
00571     // Zero the matrix/vector if it is to be assembled
00572     if (mAssembleVector && mZeroVectorBeforeAssembly)
00573     {
00574 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00575         PetscScalar zero = 0.0;
00576         VecSet(&zero, mVectorToAssemble);
00577 #else
00578         VecZeroEntries(mVectorToAssemble);
00579 #endif
00580     }
00581     if (mAssembleMatrix)
00582     {
00583         MatZeroEntries(mMatrixToAssemble);
00584     }
00585 
00586     const size_t STENCIL_SIZE=PROBLEM_DIM*(ELEMENT_DIM+1);
00587     c_matrix<double, STENCIL_SIZE, STENCIL_SIZE> a_elem;
00588     c_vector<double, STENCIL_SIZE> b_elem;
00589 
00591     // loop over elements
00593     if(mAssembleMatrix || (mAssembleVector && !mOnlyAssembleOnSurfaceElements))
00594     {
00595         for (typename AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::ElementIterator iter = mpMesh->GetElementIteratorBegin();
00596              iter != mpMesh->GetElementIteratorEnd();
00597              ++iter)
00598         {
00599             Element<ELEMENT_DIM, SPACE_DIM>& r_element = *iter;
00600     
00601             if (ElementAssemblyCriterion(r_element)==true && r_element.GetOwnership() == true)
00602             {
00603                 AssembleOnElement(r_element, a_elem, b_elem);
00604     
00605                 unsigned p_indices[STENCIL_SIZE];
00606                 r_element.GetStiffnessMatrixGlobalIndices(PROBLEM_DIM, p_indices);
00607     
00608                 if (mAssembleMatrix)
00609                 {
00610                     AddMultipleValuesToMatrix<STENCIL_SIZE>(p_indices, a_elem);
00611                 }
00612     
00613                 if (mAssembleVector)
00614                 {
00615                     AddMultipleValuesToVector<STENCIL_SIZE>(p_indices, b_elem);
00616                 }
00617             }
00618         }
00619     }
00620 
00622     // Apply any Neumann boundary conditions
00623     //
00624     // NB. We assume that if an element has a boundary condition on any unknown there is a boundary condition
00625     // on unknown 0. This can be so for any problem by adding zero constant conditions where required
00626     // although this is a bit inefficient. Proper solution involves changing BCC to have a map of arrays
00627     // boundary conditions rather than an array of maps.
00629 
00630 
00631     if (mApplyNeummanBoundaryConditionsToVector && mAssembleVector)
00632     {
00633         HeartEventHandler::EndEvent(assemble_event);
00634         ApplyNeummanBoundaryConditionsToVector();
00635         HeartEventHandler::BeginEvent(assemble_event);
00636     }
00637     
00638     HeartEventHandler::EndEvent(assemble_event);
00639 }
00640 
00641 
00642 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
00643 void AbstractFeObjectAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL>::SetApplyNeummanBoundaryConditionsToVector(BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>* pBcc)
00644 {
00645     assert(CAN_ASSEMBLE_VECTOR);
00646     mApplyNeummanBoundaryConditionsToVector = true;
00647     mpBoundaryConditions = pBcc;
00648 }
00649 
00650 
00651 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
00652 void AbstractFeObjectAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL>::ApplyNeummanBoundaryConditionsToVector()
00653 {
00654     assert(mAssembleVector);
00655     assert(mpBoundaryConditions != NULL);
00656     
00657     HeartEventHandler::BeginEvent(HeartEventHandler::NEUMANN_BCS);
00658     if (mpBoundaryConditions->AnyNonZeroNeumannConditions())
00659     {
00660         typename BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>::NeumannMapIterator
00661             neumann_iterator = mpBoundaryConditions->BeginNeumann();
00662         c_vector<double, PROBLEM_DIM*ELEMENT_DIM> b_surf_elem;
00663 
00664         // Iterate over defined conditions
00665         while (neumann_iterator != mpBoundaryConditions->EndNeumann())
00666         {
00667             const BoundaryElement<ELEMENT_DIM-1,SPACE_DIM>& r_surf_element = *(neumann_iterator->first);
00668             AssembleOnSurfaceElement(r_surf_element, b_surf_elem);
00669 
00670             const size_t STENCIL_SIZE=PROBLEM_DIM*ELEMENT_DIM; // problem_dim*num_nodes_on_surface_element
00671             unsigned p_indices[STENCIL_SIZE];
00672             r_surf_element.GetStiffnessMatrixGlobalIndices(PROBLEM_DIM, p_indices);
00673             AddMultipleValuesToVector<STENCIL_SIZE>(p_indices, b_surf_elem);
00674             ++neumann_iterator;
00675         }
00676     }
00677     HeartEventHandler::EndEvent(HeartEventHandler::NEUMANN_BCS);
00678 }
00679 
00680 
00681 
00682 
00683 
00685 // Implementation - AssembleOnElement and smaller
00687 
00688 
00689 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
00690 void AbstractFeObjectAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL>::ComputeTransformedBasisFunctionDerivatives(
00691         const ChastePoint<ELEMENT_DIM>& rPoint,
00692         const c_matrix<double, ELEMENT_DIM, SPACE_DIM>& rInverseJacobian,
00693         c_matrix<double, SPACE_DIM, ELEMENT_DIM+1>& rReturnValue)
00694 {
00695     assert(ELEMENT_DIM < 4 && ELEMENT_DIM > 0);
00696     static c_matrix<double, ELEMENT_DIM, ELEMENT_DIM+1> grad_phi;
00697 
00698     LinearBasisFunction<ELEMENT_DIM>::ComputeBasisFunctionDerivatives(rPoint, grad_phi);
00699     rReturnValue = prod(trans(rInverseJacobian), grad_phi);
00700 }
00701 
00702 
00703 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
00704 void AbstractFeObjectAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL>::AssembleOnElement(
00705     Element<ELEMENT_DIM,SPACE_DIM>& rElement,
00706     c_matrix<double, PROBLEM_DIM*(ELEMENT_DIM+1), PROBLEM_DIM*(ELEMENT_DIM+1) >& rAElem,
00707     c_vector<double, PROBLEM_DIM*(ELEMENT_DIM+1)>& rBElem)
00708 {
00714     c_matrix<double, SPACE_DIM, ELEMENT_DIM> jacobian;
00715     c_matrix<double, ELEMENT_DIM, SPACE_DIM> inverse_jacobian;
00716     double jacobian_determinant;
00717 
00718     mpMesh->GetInverseJacobianForElement(rElement.GetIndex(), jacobian, jacobian_determinant, inverse_jacobian);
00719 
00720     // With the new signature of GetInverseJacobianForElement, inverse and jacobian are returned at the same time
00721     //        // Initialise element contributions to zero
00722     //        if ( mAssembleMatrix || this->ProblemIsNonlinear() ) // don't need to construct grad_phi or grad_u in that case
00723     //        {
00724     //            this->mpMesh->GetInverseJacobianForElement(rElement.GetIndex(), inverse_jacobian);
00725     //        }
00726 
00727     if (mAssembleMatrix)
00728     {
00729         rAElem.clear();
00730     }
00731 
00732     if (mAssembleVector)
00733     {
00734         rBElem.clear();
00735     }
00736 
00737     const unsigned num_nodes = rElement.GetNumNodes();
00738 
00739     // allocate memory for the basis functions values and derivative values
00740     c_vector<double, ELEMENT_DIM+1> phi;
00741     c_matrix<double, SPACE_DIM, ELEMENT_DIM+1> grad_phi;
00742 
00743     // loop over Gauss points
00744     for (unsigned quad_index=0; quad_index < mpQuadRule->GetNumQuadPoints(); quad_index++)
00745     {
00746         const ChastePoint<ELEMENT_DIM>& quad_point = mpQuadRule->rGetQuadPoint(quad_index);
00747 
00748         BasisFunction::ComputeBasisFunctions(quad_point, phi);
00749 
00750         if ( mAssembleMatrix || INTERPOLATION_LEVEL==NONLINEAR )
00751         {
00752             ComputeTransformedBasisFunctionDerivatives(quad_point, inverse_jacobian, grad_phi);
00753         }
00754 
00755         // Location of the gauss point in the original element will be stored in x
00756         // Where applicable, u will be set to the value of the current solution at x
00757         ChastePoint<SPACE_DIM> x(0,0,0);
00758 
00759         c_vector<double,PROBLEM_DIM> u = zero_vector<double>(PROBLEM_DIM);
00760         c_matrix<double,PROBLEM_DIM,SPACE_DIM> grad_u = zero_matrix<double>(PROBLEM_DIM,SPACE_DIM);
00761 
00762         // allow the concrete version of the assembler to interpolate any
00763         // desired quantities
00764         ResetInterpolatedQuantities();
00765 
00767         // interpolation
00769         for (unsigned i=0; i<num_nodes; i++)
00770         {
00771             const Node<SPACE_DIM>* p_node = rElement.GetNode(i);
00772 
00773             if (INTERPOLATION_LEVEL != CARDIAC) // don't interpolate X if cardiac problem
00774             {
00775                 const c_vector<double, SPACE_DIM>& r_node_loc = p_node->rGetLocation();
00776                 // interpolate x
00777                 x.rGetLocation() += phi(i)*r_node_loc;
00778             }
00779 
00780             // interpolate u and grad u if a current solution or guess exists
00781             unsigned node_global_index = rElement.GetNodeGlobalIndex(i);
00782             if (mCurrentSolutionOrGuessReplicated.GetSize()>0)
00783             {
00784                 for (unsigned index_of_unknown=0; index_of_unknown<(INTERPOLATION_LEVEL!=CARDIAC ? PROBLEM_DIM : 1); index_of_unknown++)
00785                 {
00786                     // If we have a current solution (e.g. this is a dynamic problem)
00787                     // interpolate the value at this quad point
00788 
00789                     // NOTE - following assumes that, if say there are two unknowns u and v, they
00790                     // are stored in the current solution vector as
00791                     // [U1 V1 U2 V2 ... U_n V_n]
00792                     double u_at_node=GetCurrentSolutionOrGuessValue(node_global_index, index_of_unknown);
00793                     u(index_of_unknown) += phi(i)*u_at_node;
00794 
00795                     if (INTERPOLATION_LEVEL==NONLINEAR) // don't need to construct grad_phi or grad_u in other cases
00796                     {
00797                         for (unsigned j=0; j<SPACE_DIM; j++)
00798                         {
00799                             grad_u(index_of_unknown,j) += grad_phi(j,i)*u_at_node;
00800                         }
00801                     }
00802                 }
00803             }
00804 
00805             // allow the concrete version of the assembler to interpolate any
00806             // desired quantities
00807             IncrementInterpolatedQuantities(phi(i), p_node);
00808         }
00809 
00810         double wJ = jacobian_determinant * mpQuadRule->GetWeight(quad_index);
00811 
00813         // create rAElem and rBElem
00815         if (mAssembleMatrix)
00816         {
00817             noalias(rAElem) += ComputeMatrixTerm(phi, grad_phi, x, u, grad_u, &rElement) * wJ;
00818         }
00819 
00820         if (mAssembleVector)
00821         {
00822             noalias(rBElem) += ComputeVectorTerm(phi, grad_phi, x, u, grad_u, &rElement) * wJ;
00823         }
00824     }
00825 }
00826 
00827 
00828 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
00829 void AbstractFeObjectAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL>::AssembleOnSurfaceElement(
00830             const BoundaryElement<ELEMENT_DIM-1,SPACE_DIM>& rSurfaceElement,
00831             c_vector<double, PROBLEM_DIM*ELEMENT_DIM>& rBSurfElem)
00832 {
00833     c_vector<double, SPACE_DIM> weighted_direction;
00834     double jacobian_determinant;
00835     mpMesh->GetWeightedDirectionForBoundaryElement(rSurfaceElement.GetIndex(), weighted_direction, jacobian_determinant);
00836 
00837     rBSurfElem.clear();
00838 
00839     // allocate memory for the basis function values
00840     c_vector<double, ELEMENT_DIM>  phi;
00841 
00842     // loop over Gauss points
00843     for (unsigned quad_index=0; quad_index<mpSurfaceQuadRule->GetNumQuadPoints(); quad_index++)
00844     {
00845         const ChastePoint<ELEMENT_DIM-1>& quad_point = mpSurfaceQuadRule->rGetQuadPoint(quad_index);
00846 
00847         SurfaceBasisFunction::ComputeBasisFunctions(quad_point, phi);
00848 
00850         // interpolation
00852 
00853         // Location of the gauss point in the original element will be
00854         // stored in x
00855         ChastePoint<SPACE_DIM> x(0,0,0);
00856 
00857         this->ResetInterpolatedQuantities();
00858         for (unsigned i=0; i<rSurfaceElement.GetNumNodes(); i++)
00859         {
00860             const c_vector<double, SPACE_DIM> node_loc = rSurfaceElement.GetNode(i)->rGetLocation();
00861             x.rGetLocation() += phi(i)*node_loc;
00862 
00863             // allow the concrete version of the assembler to interpolate any
00864             // desired quantities
00865             IncrementInterpolatedQuantities(phi(i), rSurfaceElement.GetNode(i));
00866         }
00867 
00868         double wJ = jacobian_determinant * mpSurfaceQuadRule->GetWeight(quad_index);
00869 
00871         // create rAElem and rBElem
00874         noalias(rBSurfElem) += ComputeVectorSurfaceTerm(rSurfaceElement, phi, x) * wJ;
00875     }
00876 }
00877 
00879 //
00880 // see #1507
00881 //
00883 
00884 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
00885 template<size_t SMALL_MATRIX_SIZE>
00886 void AbstractFeObjectAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL>::AddMultipleValuesToMatrix(unsigned* matrixRowAndColIndices, c_matrix<double, SMALL_MATRIX_SIZE, SMALL_MATRIX_SIZE>& rSmallMatrix)
00887 {
00888     PetscInt matrix_row_indices[SMALL_MATRIX_SIZE];
00889     PetscInt num_rows_owned = 0;
00890     PetscInt global_row;
00891 
00892     for (unsigned row = 0; row<SMALL_MATRIX_SIZE; row++)
00893     {
00894         global_row = matrixRowAndColIndices[row];
00895 
00896         if (global_row >=mOwnershipRangeLo && global_row <mOwnershipRangeHi)
00897         {
00898             matrix_row_indices[num_rows_owned++] = global_row;
00899         }
00900     }
00901 
00902     if ( num_rows_owned == SMALL_MATRIX_SIZE)
00903     {
00904         MatSetValues(mMatrixToAssemble,
00905                      num_rows_owned,
00906                      matrix_row_indices,
00907                      SMALL_MATRIX_SIZE,
00908                      (PetscInt*) matrixRowAndColIndices,
00909                      rSmallMatrix.data(),
00910                      ADD_VALUES);
00911     }
00912     else
00913     {
00914         // We need continuous data, if some of the rows do not belong to the processor their values
00915         // are not passed to MatSetValues
00916         double values[SMALL_MATRIX_SIZE*SMALL_MATRIX_SIZE];
00917         unsigned num_values_owned = 0;
00918         for (unsigned row = 0; row<SMALL_MATRIX_SIZE; row++)
00919         {
00920             global_row = matrixRowAndColIndices[row];
00921 
00922             if (global_row >=mOwnershipRangeLo && global_row <mOwnershipRangeHi)
00923             {
00924                 for (unsigned col=0; col<SMALL_MATRIX_SIZE; col++)
00925                 {
00926                     values[num_values_owned++] = rSmallMatrix(row, col);
00927                 }
00928             }
00929         }
00930 
00931         MatSetValues(mMatrixToAssemble,
00932                      num_rows_owned,
00933                      matrix_row_indices,
00934                      SMALL_MATRIX_SIZE,
00935                      (PetscInt*) matrixRowAndColIndices,
00936                      values,
00937                      ADD_VALUES);
00938     }
00939 }
00940 
00941 
00942 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, bool CAN_ASSEMBLE_VECTOR, bool CAN_ASSEMBLE_MATRIX, InterpolationLevel INTERPOLATION_LEVEL>
00943 template<size_t SMALL_VECTOR_SIZE>
00944 void AbstractFeObjectAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL>::AddMultipleValuesToVector(unsigned* vectorIndices, c_vector<double, SMALL_VECTOR_SIZE>& smallVector)
00945 {
00946     PetscInt indices_owned[SMALL_VECTOR_SIZE];
00947     PetscInt num_indices_owned = 0;
00948     PetscInt global_row;
00949 
00950     for (unsigned row = 0; row<SMALL_VECTOR_SIZE; row++)
00951     {
00952         global_row = vectorIndices[row];
00953 
00954         if (global_row >=mOwnershipRangeLo && global_row <mOwnershipRangeHi)
00955         {
00956             indices_owned[num_indices_owned++] = global_row;
00957         }
00958     }
00959 
00960     if (num_indices_owned == SMALL_VECTOR_SIZE)
00961     {
00962         VecSetValues(mVectorToAssemble,
00963                      num_indices_owned,
00964                      indices_owned,
00965                      smallVector.data(),
00966                      ADD_VALUES);
00967     }
00968     else
00969     {
00970         // We need continuous data, if some of the rows do not belong to the processor their values
00971         // are not passed to MatSetValues
00972         double values[SMALL_VECTOR_SIZE];
00973         unsigned num_values_owned = 0;
00974 
00975         for (unsigned row = 0; row<SMALL_VECTOR_SIZE; row++)
00976         {
00977             global_row = vectorIndices[row];
00978 
00979             if (global_row >=mOwnershipRangeLo && global_row <mOwnershipRangeHi)
00980             {
00981                 values[num_values_owned++] = smallVector(row);
00982             }
00983         }
00984 
00985         VecSetValues(mVectorToAssemble,
00986                      num_indices_owned,
00987                      indices_owned,
00988                      values,
00989                      ADD_VALUES);
00990     }
00991 }
00992     
00993 
00994 #endif /*ABSTRACTFEOBJECTASSEMBLER_HPP_*/

Generated on Mon Nov 1 12:35:24 2010 for Chaste by  doxygen 1.5.5