Chaste Release::3.1
AbstractNonlinearElasticitySolver.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 ABSTRACTNONLINEARELASTICITYSOLVER_HPP_
00037 #define ABSTRACTNONLINEARELASTICITYSOLVER_HPP_
00038 
00039 #include <vector>
00040 #include <cmath>
00041 #include "AbstractContinuumMechanicsSolver.hpp"
00042 #include "LinearSystem.hpp"
00043 #include "LogFile.hpp"
00044 #include "MechanicsEventHandler.hpp"
00045 #include "ReplicatableVector.hpp"
00046 #include "FourthOrderTensor.hpp"
00047 #include "CmguiDeformedSolutionsWriter.hpp"
00048 #include "AbstractMaterialLaw.hpp"
00049 #include "QuadraticBasisFunction.hpp"
00050 #include "SolidMechanicsProblemDefinition.hpp"
00051 #include "Timer.hpp"
00052 #include "petscsnes.h"
00053 
00054 
00055 //#define MECH_USE_HYPRE    // uses HYPRE to solve linear systems, requires PETSc to be installed with HYPRE
00056 
00057 
00058 
00059 // Bizarrely PETSc 2.2 has this, but doesn't put it in the petscksp.h header...
00060 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00061 extern PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);
00062 #endif
00063 
00064 
00066 //  Globals functions used by the SNES solver
00068 
00078 template<unsigned DIM>
00079 PetscErrorCode AbstractNonlinearElasticitySolver_ComputeResidual(SNES snes,
00080                                                                  Vec currentGuess,
00081                                                                  Vec residualVector,
00082                                                                  void* pContext);
00083 
00094 template<unsigned DIM>
00095 PetscErrorCode AbstractNonlinearElasticitySolver_ComputeJacobian(SNES snes,
00096                                                                  Vec currentGuess,
00097                                                                  Mat* pGlobalJacobian,
00098                                                                  Mat* pPreconditioner,
00099                                                                  MatStructure* pMatStructure,
00100                                                                  void* pContext);
00101 
00102 
00114 template<unsigned DIM>
00115 class AbstractNonlinearElasticitySolver : public AbstractContinuumMechanicsSolver<DIM>
00116 {
00117 protected:
00118 
00120     static const size_t NUM_VERTICES_PER_ELEMENT = DIM+1;
00121 
00123     static const size_t NUM_NODES_PER_ELEMENT = (DIM+1)*(DIM+2)/2;      // assuming quadratic
00124 
00126     static const size_t NUM_NODES_PER_BOUNDARY_ELEMENT = DIM*(DIM+1)/2; // assuming quadratic
00127 
00132     static const size_t BOUNDARY_STENCIL_SIZE = DIM*NUM_NODES_PER_BOUNDARY_ELEMENT;
00133 
00141     static double MAX_NEWTON_ABS_TOL;
00142 
00144     static double MIN_NEWTON_ABS_TOL;
00145 
00147     static double NEWTON_REL_TOL;
00148 
00153     SolidMechanicsProblemDefinition<DIM>& mrProblemDefinition;
00154 
00159     Mat& mrJacobianMatrix;
00160 
00167     c_matrix<double,DIM,DIM> mChangeOfBasisMatrix;
00168 
00174     double mKspAbsoluteTol;
00175 
00181     bool mWriteOutputEachNewtonIteration;
00182 
00186     FourthOrderTensor<DIM,DIM,DIM,DIM> dTdE;
00187 
00189     unsigned mNumNewtonIterations;
00190 
00196     double mCurrentTime;
00197 
00198 
00202     bool mCheckedOutwardNormals;
00203 
00208     bool mUseSnesSolver;
00209 
00215     double mLastDampingValue;
00216 
00217 
00225     bool mIncludeActiveTension;
00226 
00231     bool mSetComputeAverageStressPerElement;
00232 
00240     std::vector<c_vector<double,DIM*(DIM+1)/2> > mAverageStressesPerElement;
00241 
00249     void AddStressToAverageStressPerElement(c_matrix<double,DIM,DIM>& rT, unsigned elementIndex);
00250 
00259     virtual void SetKspSolverAndPcType(KSP solver);
00260 
00261 
00274     virtual void AssembleSystem(bool assembleResidual, bool assembleLinearSystem)=0;
00275 
00276 
00277 
00285     virtual void FinishAssembleSystem(bool assembleResidual, bool assembleLinearSystem);
00286 
00287 
00289     //
00290     // Element level methods
00291     //
00293 
00300     void GetElementCentroidDeformationGradient(Element<DIM,DIM>& rElement,
00301                                                c_matrix<double,DIM,DIM>& rDeformationGradient);
00302 
00321     virtual void AddActiveStressAndStressDerivative(c_matrix<double,DIM,DIM>& rC,
00322                                                     unsigned elementIndex,
00323                                                     unsigned currentQuadPointGlobalIndex,
00324                                                     c_matrix<double,DIM,DIM>& rT,
00325                                                     FourthOrderTensor<DIM,DIM,DIM,DIM>& rDTdE,
00326                                                     bool addToDTdE)
00327     {
00328         // does nothing - subclass needs to overload this if there are active stresses
00329     }
00330 
00338     virtual void SetupChangeOfBasisMatrix(unsigned elementIndex, unsigned currentQuadPointGlobalIndex)
00339     {
00340     }
00341 
00342 
00343 
00362     void AssembleOnBoundaryElement(BoundaryElement<DIM-1, DIM>& rBoundaryElement,
00363                                    c_matrix<double, BOUNDARY_STENCIL_SIZE, BOUNDARY_STENCIL_SIZE>& rAelem,
00364                                    c_vector<double, BOUNDARY_STENCIL_SIZE>& rBelem,
00365                                    bool assembleResidual,
00366                                    bool assembleJacobian,
00367                                    unsigned boundaryConditionIndex);
00368 
00369 
00383     bool ShouldAssembleMatrixTermForPressureOnDeformedBc();
00384 
00404     void AssembleOnBoundaryElementForPressureOnDeformedBc(BoundaryElement<DIM-1,DIM>& rBoundaryElement,
00405                                                           c_matrix<double,BOUNDARY_STENCIL_SIZE,BOUNDARY_STENCIL_SIZE>& rAelem,
00406                                                           c_vector<double,BOUNDARY_STENCIL_SIZE>& rBelem,
00407                                                           bool assembleResidual,
00408                                                           bool assembleJacobian,
00409                                                           unsigned boundaryConditionIndex);
00410 
00412     //
00413     //    These methods form the non-SNES nonlinear solver
00414     //
00416 
00427     double ComputeResidualAndGetNorm(bool allowException);
00428 
00430     double CalculateResidualNorm();
00431 
00441     void VectorSum(std::vector<double>& rX, ReplicatableVector& rY, double a, std::vector<double>& rZ);
00442 
00450     void PrintLineSearchResult(double s, double residNorm);
00451 
00459     double TakeNewtonStep();
00460 
00469     double UpdateSolutionUsingLineSearch(Vec solution);
00470 
00478     virtual void PostNewtonStep(unsigned counter, double normResidual);
00479 
00485     void SolveNonSnes(double tol=-1.0);
00486 
00487 
00489     //
00490     //    These methods form the SNES nonlinear solver
00491     //
00493 public: // need to be public as are called by global functions
00500     void ComputeResidual(Vec currentGuess, Vec residualVector);
00501 
00509     void ComputeJacobian(Vec currentGuess, Mat* pJacobian, Mat* pPreconditioner);
00510 
00511 private:
00516     void SolveSnes();
00517 
00518 public:
00519 
00530     AbstractNonlinearElasticitySolver(QuadraticMesh<DIM>& rQuadMesh,
00531                                       SolidMechanicsProblemDefinition<DIM>& rProblemDefinition,
00532                                       std::string outputDirectory,
00533                                       CompressibilityType compressibilityType);
00534 
00538     virtual ~AbstractNonlinearElasticitySolver();
00539 
00545     void Solve(double tol=-1.0);
00546 
00556     void SetIncludeActiveTension(bool includeActiveTension = true)
00557     {
00558         mIncludeActiveTension = includeActiveTension;
00559     }
00560 
00564     unsigned GetNumNewtonIterations();
00565 
00566 
00574     void SetWriteOutputEachNewtonIteration(bool writeOutputEachNewtonIteration=true)
00575     {
00576         mWriteOutputEachNewtonIteration = writeOutputEachNewtonIteration;
00577     }
00578 
00585     void SetKspAbsoluteTolerance(double kspAbsoluteTolerance)
00586     {
00587         assert(kspAbsoluteTolerance > 0);
00588         mKspAbsoluteTol = kspAbsoluteTolerance;
00589     }
00590 
00597     void SetCurrentTime(double time)
00598     {
00599         mCurrentTime = time;
00600     }
00601 
00606     void CreateCmguiOutput();
00607 
00608 
00620     void WriteCurrentDeformationGradients(std::string fileName, int counterToAppend = -1);
00621 
00628     void SetComputeAverageStressPerElementDuringSolve(bool setComputeAverageStressPerElement = true);
00629 
00643     void WriteCurrentAverageElementStresses(std::string fileName, int counterToAppend = -1);
00644 
00649     std::vector<c_vector<double,DIM> >& rGetSpatialSolution();
00650 
00655     std::vector<c_vector<double,DIM> >& rGetDeformedPosition();
00656 
00664     c_matrix<double,DIM,DIM> GetAverageStressPerElement(unsigned elementIndex);
00665 };
00666 
00668 // Implementation: first, the non-nonlinear-solve methods
00670 
00671 template<unsigned DIM>
00672 AbstractNonlinearElasticitySolver<DIM>::AbstractNonlinearElasticitySolver(QuadraticMesh<DIM>& rQuadMesh,
00673                                                                           SolidMechanicsProblemDefinition<DIM>& rProblemDefinition,
00674                                                                           std::string outputDirectory,
00675                                                                           CompressibilityType compressibilityType)
00676     : AbstractContinuumMechanicsSolver<DIM>(rQuadMesh, rProblemDefinition, outputDirectory, compressibilityType),
00677       mrProblemDefinition(rProblemDefinition),
00678       mrJacobianMatrix(this->mSystemLhsMatrix),
00679       mKspAbsoluteTol(-1),
00680       mWriteOutputEachNewtonIteration(false),
00681       mNumNewtonIterations(0),
00682       mCurrentTime(0.0),
00683       mCheckedOutwardNormals(false),
00684       mLastDampingValue(0.0),
00685       mIncludeActiveTension(true),
00686       mSetComputeAverageStressPerElement(false)
00687 {
00688     mUseSnesSolver = (mrProblemDefinition.GetSolveUsingSnes() ||
00689                       CommandLineArguments::Instance()->OptionExists("-mech_use_snes") );
00690 
00691     mChangeOfBasisMatrix = identity_matrix<double>(DIM,DIM);
00692 }
00693 
00694 template<unsigned DIM>
00695 AbstractNonlinearElasticitySolver<DIM>::~AbstractNonlinearElasticitySolver()
00696 {
00697 }
00698 
00699 
00700 template<unsigned DIM>
00701 void AbstractNonlinearElasticitySolver<DIM>::FinishAssembleSystem(bool assembleResidual, bool assembleJacobian)
00702 {
00703     PetscVecTools::Finalise(this->mResidualVector);
00704 
00705     if (assembleJacobian)
00706     {
00707         PetscMatTools::SwitchWriteMode(mrJacobianMatrix);
00708         PetscMatTools::SwitchWriteMode(this->mPreconditionMatrix);
00709 
00710         VecCopy(this->mResidualVector, this->mLinearSystemRhsVector);
00711     }
00712 
00713     // Apply Dirichlet boundary conditions
00714     if(assembleJacobian)
00715     {
00716         this->ApplyDirichletBoundaryConditions(NONLINEAR_PROBLEM_APPLY_TO_EVERYTHING, this->mCompressibilityType==COMPRESSIBLE);
00717     }
00718     else if (assembleResidual)
00719     {
00720         this->ApplyDirichletBoundaryConditions(NONLINEAR_PROBLEM_APPLY_TO_RESIDUAL_ONLY, this->mCompressibilityType==COMPRESSIBLE);
00721     }
00722 
00723     if (assembleResidual)
00724     {
00725         PetscVecTools::Finalise(this->mResidualVector);
00726     }
00727     if (assembleJacobian)
00728     {
00729         PetscMatTools::Finalise(mrJacobianMatrix);
00730         PetscMatTools::Finalise(this->mPreconditionMatrix);
00731         PetscVecTools::Finalise(this->mLinearSystemRhsVector);
00732     }
00733 }
00734 
00735 
00736 template<unsigned DIM>
00737 std::vector<c_vector<double,DIM> >& AbstractNonlinearElasticitySolver<DIM>::rGetSpatialSolution()
00738 {
00739     this->mSpatialSolution.clear();
00740     this->mSpatialSolution.resize(this->mrQuadMesh.GetNumNodes(), zero_vector<double>(DIM));
00741     for (unsigned i=0; i<this->mrQuadMesh.GetNumNodes(); i++)
00742     {
00743         for (unsigned j=0; j<DIM; j++)
00744         {
00745             this->mSpatialSolution[i](j) = this->mrQuadMesh.GetNode(i)->rGetLocation()[j] + this->mCurrentSolution[this->mProblemDimension*i+j];
00746         }
00747     }
00748     return this->mSpatialSolution;
00749 }
00750 
00751 template<unsigned DIM>
00752 std::vector<c_vector<double,DIM> >& AbstractNonlinearElasticitySolver<DIM>::rGetDeformedPosition()
00753 {
00754     return rGetSpatialSolution();
00755 }
00756 
00757 
00758 template<unsigned DIM>
00759 void AbstractNonlinearElasticitySolver<DIM>::WriteCurrentDeformationGradients(std::string fileName, int counterToAppend)
00760 {
00761     if (!this->mWriteOutput)
00762     {
00763         return;
00764     }
00765 
00766     std::stringstream file_name;
00767     file_name << fileName;
00768     if (counterToAppend >= 0)
00769     {
00770         file_name << "_" << counterToAppend;
00771     }
00772     file_name << ".strain";
00773 
00774     out_stream p_file = this->mpOutputFileHandler->OpenOutputFile(file_name.str());
00775 
00776     c_matrix<double,DIM,DIM> deformation_gradient;
00777 
00778     for (typename AbstractTetrahedralMesh<DIM,DIM>::ElementIterator iter = this->mrQuadMesh.GetElementIteratorBegin();
00779          iter != this->mrQuadMesh.GetElementIteratorEnd();
00780          ++iter)
00781     {
00782         GetElementCentroidDeformationGradient(*iter, deformation_gradient);
00783         for(unsigned i=0; i<DIM; i++)
00784         {
00785             for(unsigned j=0; j<DIM; j++)
00786             {
00787                 *p_file << deformation_gradient(i,j) << " ";
00788             }
00789         }
00790         *p_file << "\n";
00791     }
00792     p_file->close();
00793 }
00794 
00795 
00796 template<unsigned DIM>
00797 void AbstractNonlinearElasticitySolver<DIM>::WriteCurrentAverageElementStresses(std::string fileName, int counterToAppend)
00798 {
00799     if (!this->mWriteOutput)
00800     {
00801         return;
00802     }
00803 
00804     if(!mSetComputeAverageStressPerElement)
00805     {
00806         EXCEPTION("Call SetComputeAverageStressPerElementDuringSolve() before solve if calling WriteCurrentAverageElementStresses()");
00807     }
00808 
00809     std::stringstream file_name;
00810     file_name << fileName;
00811     if (counterToAppend >= 0)
00812     {
00813         file_name << "_" << counterToAppend;
00814     }
00815     file_name << ".stress";
00816 
00817     out_stream p_file = this->mpOutputFileHandler->OpenOutputFile(file_name.str());
00818 
00819     assert(mAverageStressesPerElement.size()==this->mrQuadMesh.GetNumElements());
00820 
00821     for(unsigned i=0; i<mAverageStressesPerElement.size(); i++)
00822     {
00823         c_matrix<double,DIM,DIM> stress = GetAverageStressPerElement(i);
00824 //        c_vector<double,DIM> centroid = this->mrQuadMesh.GetElement(i)->CalculateCentroid();
00825 //
00826 //        for(unsigned j=0; j<DIM; j++)
00827 //        {
00828 //            *p_file << centroid(j) << " ";
00829 //        }
00830 
00831         for(unsigned j=0; j<DIM; j++)
00832         {
00833             for(unsigned k=0; k<DIM; k++)
00834             {
00835                 *p_file << stress(j,k) << " ";
00836             }
00837         }
00838         *p_file << "\n";
00839     }
00840     p_file->close();
00841 }
00842 
00843 
00844 template<unsigned DIM>
00845 void AbstractNonlinearElasticitySolver<DIM>::CreateCmguiOutput()
00846 {
00847     if (this->mOutputDirectory=="")
00848     {
00849         EXCEPTION("No output directory was given so no output was written, cannot convert to cmgui format");
00850     }
00851 
00852     CmguiDeformedSolutionsWriter<DIM> writer(this->mOutputDirectory + "/cmgui",
00853                                              "solution",
00854                                              this->mrQuadMesh,
00855                                              WRITE_QUADRATIC_MESH);
00856 
00857     std::vector<c_vector<double,DIM> >& r_deformed_positions = this->rGetDeformedPosition();
00858     writer.WriteInitialMesh(); // this writes solution_0.exnode and .exelem
00859     writer.WriteDeformationPositions(r_deformed_positions, 1); // this writes the final solution as solution_1.exnode
00860     writer.WriteCmguiScript(); // writes LoadSolutions.com
00861 }
00862 
00863 template<unsigned DIM>
00864 void AbstractNonlinearElasticitySolver<DIM>::SetComputeAverageStressPerElementDuringSolve(bool setComputeAverageStressPerElement)
00865 {
00866     if(setComputeAverageStressPerElement && PetscTools::IsParallel())
00867     {
00868         EXCEPTION("SetComputeAverageStressPerElementDuringSolve() is not yet implemented for parallel simulations");
00869     }
00870 
00871     mSetComputeAverageStressPerElement = setComputeAverageStressPerElement;
00872     if(setComputeAverageStressPerElement && mAverageStressesPerElement.size()==0)
00873     {
00874         mAverageStressesPerElement.resize(this->mrQuadMesh.GetNumElements(), zero_vector<double>(DIM*(DIM+1)/2));
00875     }
00876 }
00877 
00878 template<unsigned DIM>
00879 void AbstractNonlinearElasticitySolver<DIM>::AddStressToAverageStressPerElement(c_matrix<double,DIM,DIM>& rT, unsigned elemIndex)
00880 {
00881     assert(mSetComputeAverageStressPerElement);
00882     assert(elemIndex<this->mrQuadMesh.GetNumElements());
00883 
00884     // In 2d the matrix is
00885     // [T00 T01]
00886     // [T10 T11]
00887     // where T01 = T10. We store this as a vector
00888     // [T00 T01 T11]
00889     //
00890     // Similarly, for 3d we store
00891     // [T00 T01 T02 T11 T12 T22]
00892     for(unsigned i=0; i<DIM*(DIM+1)/2; i++)
00893     {
00894         unsigned row;
00895         unsigned col;
00896         if(DIM==2)
00897         {
00898             row = i<=1 ? 0 : 1;
00899             col = i==0 ? 0 : 1;
00900         }
00901         else // DIM==3
00902         {
00903             row = i<=2 ? 0 : (i<=4? 1 : 2);
00904             col = i==0 ? 0 : (i==1 || i==3? 1 : 2);
00905         }
00906 
00907         this->mAverageStressesPerElement[elemIndex](i) += rT(row,col);
00908     }
00909 }
00910 
00911 
00912 template<unsigned DIM>
00913 c_matrix<double,DIM,DIM> AbstractNonlinearElasticitySolver<DIM>::GetAverageStressPerElement(unsigned elementIndex)
00914 {
00915     if(!mSetComputeAverageStressPerElement)
00916     {
00917         EXCEPTION("Call SetComputeAverageStressPerElementDuringSolve() before solve if calling GetAverageStressesPerElement()");
00918     }
00919     assert(elementIndex<this->mrQuadMesh.GetNumElements());
00920 
00921     c_matrix<double,DIM,DIM> stress;
00922 
00923     // In 2d the matrix is
00924     // [T00 T01]
00925     // [T10 T11]
00926     // where T01 = T10, and was stored as
00927     // [T00 T01 T11]
00928     //
00929     // Similarly, for 3d the matrix was stored as
00930     // [T00 T01 T02 T11 T12 T22]
00931     if(DIM==2)
00932     {
00933         stress(0,0) = mAverageStressesPerElement[elementIndex](0);
00934         stress(1,0) = stress(0,1) = mAverageStressesPerElement[elementIndex](1);
00935         stress(1,1) = mAverageStressesPerElement[elementIndex](2);
00936     }
00937     else
00938     {
00939         stress(0,0) = mAverageStressesPerElement[elementIndex](0);
00940         stress(1,0) = stress(0,1) = mAverageStressesPerElement[elementIndex](1);
00941         stress(2,0) = stress(0,2) = mAverageStressesPerElement[elementIndex](2);
00942         stress(1,1) = mAverageStressesPerElement[elementIndex](3);
00943         stress(2,1) = stress(1,2) = mAverageStressesPerElement[elementIndex](4);
00944         stress(2,2) = mAverageStressesPerElement[elementIndex](5);
00945     }
00946 
00947     return stress;
00948 }
00949 
00950 
00952 // Methods at the 'element level'.
00954 
00955 template<unsigned DIM>
00956 void AbstractNonlinearElasticitySolver<DIM>::GetElementCentroidDeformationGradient(Element<DIM,DIM>& rElement,
00957                                                                                    c_matrix<double,DIM,DIM>& rDeformationGradient)
00958 {
00959     static c_matrix<double,DIM,DIM> jacobian;
00960     static c_matrix<double,DIM,DIM> inverse_jacobian;
00961     double jacobian_determinant;
00962 
00963     this->mrQuadMesh.GetInverseJacobianForElement(rElement.GetIndex(), jacobian, jacobian_determinant, inverse_jacobian);
00964 
00965     // Get the current displacement at the nodes
00966     static c_matrix<double,DIM,NUM_NODES_PER_ELEMENT> element_current_displacements;
00967     for (unsigned II=0; II<NUM_NODES_PER_ELEMENT; II++)
00968     {
00969         for (unsigned JJ=0; JJ<DIM; JJ++)
00970         {
00971             element_current_displacements(JJ,II) = this->mCurrentSolution[this->mProblemDimension*rElement.GetNodeGlobalIndex(II) + JJ];
00972         }
00973     }
00974 
00975     // Allocate memory for the basis functions values and derivative values
00976     static c_matrix<double, DIM, NUM_NODES_PER_ELEMENT> grad_quad_phi;
00977     static c_matrix<double,DIM,DIM> grad_u; // grad_u = (du_i/dX_M)
00978 
00979     // we need the point in the canonical element which corresponds to the centroid of the
00980     // version of the element in physical space. This point can be shown to be (1/3,1/3).
00981     ChastePoint<DIM> quadrature_point;
00982     if(DIM==2)
00983     {
00984         quadrature_point.rGetLocation()(0) = 1.0/3.0;
00985         quadrature_point.rGetLocation()(1) = 1.0/3.0;
00986     }
00987     else
00988     {
00989         assert(DIM==3);
00990         quadrature_point.rGetLocation()(0) = 1.0/4.0;
00991         quadrature_point.rGetLocation()(1) = 1.0/4.0;
00992         quadrature_point.rGetLocation()(2) = 1.0/4.0;
00993     }
00994 
00995     QuadraticBasisFunction<DIM>::ComputeTransformedBasisFunctionDerivatives(quadrature_point, inverse_jacobian, grad_quad_phi);
00996 
00997     // Interpolate grad_u
00998     grad_u = zero_matrix<double>(DIM,DIM);
00999     for (unsigned node_index=0; node_index<NUM_NODES_PER_ELEMENT; node_index++)
01000     {
01001         for (unsigned i=0; i<DIM; i++)
01002         {
01003             for (unsigned M=0; M<DIM; M++)
01004             {
01005                 grad_u(i,M) += grad_quad_phi(M,node_index)*element_current_displacements(i,node_index);
01006             }
01007         }
01008     }
01009 
01010     for (unsigned i=0; i<DIM; i++)
01011     {
01012         for (unsigned M=0; M<DIM; M++)
01013         {
01014             rDeformationGradient(i,M) = (i==M?1:0) + grad_u(i,M);
01015         }
01016     }
01017 }
01018 
01019 
01020 
01021 template<unsigned DIM>
01022 void AbstractNonlinearElasticitySolver<DIM>::AssembleOnBoundaryElement(
01023             BoundaryElement<DIM-1,DIM>& rBoundaryElement,
01024             c_matrix<double,BOUNDARY_STENCIL_SIZE,BOUNDARY_STENCIL_SIZE>& rAelem,
01025             c_vector<double,BOUNDARY_STENCIL_SIZE>& rBelem,
01026             bool assembleResidual,
01027             bool assembleJacobian,
01028             unsigned boundaryConditionIndex)
01029 {
01030     if(   this->mrProblemDefinition.GetTractionBoundaryConditionType() == PRESSURE_ON_DEFORMED
01031        || this->mrProblemDefinition.GetTractionBoundaryConditionType() == FUNCTIONAL_PRESSURE_ON_DEFORMED)
01032     {
01033         AssembleOnBoundaryElementForPressureOnDeformedBc(rBoundaryElement, rAelem, rBelem,
01034                                                          assembleResidual, assembleJacobian, boundaryConditionIndex);
01035         return;
01036     }
01037 
01038     rAelem.clear();
01039     rBelem.clear();
01040 
01041     if (assembleJacobian && !assembleResidual)
01042     {
01043         // Nothing to do
01044         return;
01045     }
01046 
01047     c_vector<double, DIM> weighted_direction;
01048     double jacobian_determinant;
01049     this->mrQuadMesh.GetWeightedDirectionForBoundaryElement(rBoundaryElement.GetIndex(), weighted_direction, jacobian_determinant);
01050 
01051     c_vector<double,NUM_NODES_PER_BOUNDARY_ELEMENT> phi;
01052 
01053     for (unsigned quad_index=0; quad_index<this->mpBoundaryQuadratureRule->GetNumQuadPoints(); quad_index++)
01054     {
01055         double wJ = jacobian_determinant * this->mpBoundaryQuadratureRule->GetWeight(quad_index);
01056 
01057         const ChastePoint<DIM-1>& quad_point = this->mpBoundaryQuadratureRule->rGetQuadPoint(quad_index);
01058 
01059         QuadraticBasisFunction<DIM-1>::ComputeBasisFunctions(quad_point, phi);
01060 
01061         // Get the required traction, interpolating X (slightly inefficiently,
01062         // as interpolating using quad bases) if necessary
01063         c_vector<double,DIM> traction = zero_vector<double>(DIM);
01064         switch (this->mrProblemDefinition.GetTractionBoundaryConditionType())
01065         {
01066             case FUNCTIONAL_TRACTION:
01067             {
01068                 c_vector<double,DIM> X = zero_vector<double>(DIM);
01069                 for (unsigned node_index=0; node_index<NUM_NODES_PER_BOUNDARY_ELEMENT; node_index++)
01070                 {
01071                     X += phi(node_index)*this->mrQuadMesh.GetNode( rBoundaryElement.GetNodeGlobalIndex(node_index) )->rGetLocation();
01072                 }
01073                 traction = this->mrProblemDefinition.EvaluateTractionFunction(X, this->mCurrentTime);
01074                 break;
01075             }
01076             case ELEMENTWISE_TRACTION:
01077             {
01078                 traction = this->mrProblemDefinition.rGetElementwiseTractions()[boundaryConditionIndex];
01079                 break;
01080             }
01081             default:
01082                 NEVER_REACHED;
01083         }
01084 
01085 
01086         for (unsigned index=0; index<NUM_NODES_PER_BOUNDARY_ELEMENT*DIM; index++)
01087         {
01088             unsigned spatial_dim = index%DIM;
01089             unsigned node_index = (index-spatial_dim)/DIM;
01090 
01091             assert(node_index < NUM_NODES_PER_BOUNDARY_ELEMENT);
01092 
01093             rBelem(index) -=   traction(spatial_dim)
01094                              * phi(node_index)
01095                              * wJ;
01096         }
01097     }
01098 }
01099 
01100 template<unsigned DIM>
01101 bool AbstractNonlinearElasticitySolver<DIM>::ShouldAssembleMatrixTermForPressureOnDeformedBc()
01102 {
01103     if(mUseSnesSolver)
01104     {
01105         // although not using this in the first few steps might be useful when the deformation
01106         // is large, the snes solver is more robust, so we have this on all the time. (Also because
01107         // for cardiac problems and in timesteps after the initial large deformation, we want this on
01108         // in the first step
01109         return true;
01110 
01111         // could do something like this, if make the snes a member variable
01112         //PetscInt iteration_number;
01113         //SNESGetIterationNumber(mSnes,&iteration_number);
01114         //return (iteration_number >= 3);
01115     }
01116     else
01117     {
01118         return (mLastDampingValue >= 0.5);
01119     }
01120 
01121 }
01122 
01123 template<unsigned DIM>
01124 void AbstractNonlinearElasticitySolver<DIM>::AssembleOnBoundaryElementForPressureOnDeformedBc(
01125             BoundaryElement<DIM-1,DIM>& rBoundaryElement,
01126             c_matrix<double,BOUNDARY_STENCIL_SIZE,BOUNDARY_STENCIL_SIZE>& rAelem,
01127             c_vector<double,BOUNDARY_STENCIL_SIZE>& rBelem,
01128             bool assembleResidual,
01129             bool assembleJacobian,
01130             unsigned boundaryConditionIndex)
01131 {
01132     assert(   this->mrProblemDefinition.GetTractionBoundaryConditionType()==PRESSURE_ON_DEFORMED
01133            || this->mrProblemDefinition.GetTractionBoundaryConditionType()==FUNCTIONAL_PRESSURE_ON_DEFORMED);
01134 
01135     rAelem.clear();
01136     rBelem.clear();
01137 
01138     c_vector<double, DIM> weighted_direction;
01139     double jacobian_determinant;
01140     // note: jacobian determinant may be over-written below
01141     this->mrQuadMesh.GetWeightedDirectionForBoundaryElement(rBoundaryElement.GetIndex(), weighted_direction, jacobian_determinant);
01142 
01143 
01145     // Find the volume element of the mesh which
01146     // contains this boundary element
01148 
01149     Element<DIM,DIM>* p_containing_vol_element = NULL;
01150 
01151     std::set<unsigned> potential_elements = rBoundaryElement.GetNode(0)->rGetContainingElementIndices();
01152     for(std::set<unsigned>::iterator iter = potential_elements.begin();
01153         iter != potential_elements.end();
01154         iter++)
01155     {
01156         p_containing_vol_element = this->mrQuadMesh.GetElement(*iter);
01157 
01158         bool this_vol_ele_contains_surf_ele = true;
01159         // loop over the nodes of boundary element and see if they are in the volume element
01160         for(unsigned i=1; i<NUM_NODES_PER_BOUNDARY_ELEMENT; i++) // don't need to start at 0, given looping over contain elems of node 0
01161         {
01162             unsigned surf_element_node_index = rBoundaryElement.GetNodeGlobalIndex(i);
01163             bool found_this_node = false;
01164             for(unsigned j=0; j<p_containing_vol_element->GetNumNodes(); j++)
01165             {
01166                 unsigned vol_element_node_index = p_containing_vol_element->GetNodeGlobalIndex(j);
01167                 if(surf_element_node_index == vol_element_node_index)
01168                 {
01169                     found_this_node = true;
01170                     break;
01171                 }
01172             }
01173             if(!found_this_node)
01174             {
01175                 this_vol_ele_contains_surf_ele = false;
01176                 break;
01177             }
01178         }
01179         if(this_vol_ele_contains_surf_ele)
01180         {
01181             break;
01182         }
01183     }
01184 
01185     // Find the local node index in the volume element for each node in the boundary element
01186     std::vector<unsigned> surf_to_vol_map(NUM_NODES_PER_BOUNDARY_ELEMENT);
01187     for(unsigned i=0; i<NUM_NODES_PER_BOUNDARY_ELEMENT; i++)
01188     {
01189         unsigned index = rBoundaryElement.GetNodeGlobalIndex(i);
01190         for(unsigned j=0; j<NUM_NODES_PER_ELEMENT; j++)
01191         {
01192             if(p_containing_vol_element->GetNodeGlobalIndex(j)==index)
01193             {
01194                 surf_to_vol_map[i] = j;
01195                 break;
01196             }
01197         }
01198     }
01199 
01200 
01201     // We require the volume element to compute F, which requires grad_phi on the volume element. For this we will
01202     // need the inverse jacobian for the volume element
01203     static c_matrix<double,DIM,DIM> jacobian_vol_element;
01204     static c_matrix<double,DIM,DIM> inverse_jacobian_vol_element;
01205     double jacobian_determinant_vol_element;
01206     this->mrQuadMesh.GetInverseJacobianForElement(p_containing_vol_element->GetIndex(), jacobian_vol_element, jacobian_determinant_vol_element, inverse_jacobian_vol_element);
01207 
01208     // Get the current displacements at each node of the volume element, to be used in computing F
01209     static c_matrix<double,DIM,NUM_NODES_PER_ELEMENT> element_current_displacements;
01210     for (unsigned II=0; II<NUM_NODES_PER_ELEMENT; II++)
01211     {
01212         for (unsigned JJ=0; JJ<DIM; JJ++)
01213         {
01214             element_current_displacements(JJ,II) = this->mCurrentSolution[this->mProblemDimension*p_containing_vol_element->GetNodeGlobalIndex(II) + JJ];
01215         }
01216     }
01217 
01218 
01219     // We will need both {grad phi_i} for the quadratic bases of the volume element, for computing F..
01220     static c_matrix<double, DIM, NUM_NODES_PER_ELEMENT> grad_quad_phi_vol_element;
01221     // ..the phi_i for each of the quadratic bases of the surface element, for the standard FE assembly part.
01222     c_vector<double,NUM_NODES_PER_BOUNDARY_ELEMENT> quad_phi_surf_element;
01223     // We need this too, which is obtained by taking a subset of grad_quad_phi_vol_element
01224     static c_matrix<double, DIM, NUM_NODES_PER_BOUNDARY_ELEMENT> grad_quad_phi_surf_element;
01225 
01226     c_matrix<double,DIM,DIM> F;
01227     c_matrix<double,DIM,DIM> invF;
01228 
01229     c_vector<double,DIM> normal = rBoundaryElement.CalculateNormal();
01230     c_matrix<double,1,DIM> normal_as_mat;
01231     for(unsigned i=0; i<DIM; i++)
01232     {
01233         normal_as_mat(0,i) = normal(i);
01234     }
01235 
01236     double normal_pressure;
01237     switch (this->mrProblemDefinition.GetTractionBoundaryConditionType())
01238     {
01239         case PRESSURE_ON_DEFORMED:
01240             normal_pressure = this->mrProblemDefinition.GetNormalPressure();
01241             break;
01242         case FUNCTIONAL_PRESSURE_ON_DEFORMED:
01243             normal_pressure = this->mrProblemDefinition.EvaluateNormalPressureFunction(this->mCurrentTime);
01244             break;
01245         default:
01246             NEVER_REACHED;
01247     }
01248 
01249 
01250 
01251     for (unsigned quad_index=0; quad_index<this->mpBoundaryQuadratureRule->GetNumQuadPoints(); quad_index++)
01252     {
01253         double wJ = jacobian_determinant * this->mpBoundaryQuadratureRule->GetWeight(quad_index);
01254 
01255         // Get the quadrature point on this surface element (in canonical space) - so eg, for a 2D problem,
01256         // the quad point is in 1D space
01257         const ChastePoint<DIM-1>& quadrature_point = this->mpBoundaryQuadratureRule->rGetQuadPoint(quad_index);
01258         QuadraticBasisFunction<DIM-1>::ComputeBasisFunctions(quadrature_point, quad_phi_surf_element);
01259 
01260         // We will need the xi coordinates of this quad point in the volume element. We could do this by figuring
01261         // out how the nodes of the surface element are ordered in the list of nodes in the volume element,
01262         // however it is less fiddly to compute directly. Firstly, compute the corresponding physical location
01263         // of the quad point, by interpolating
01264         c_vector<double,DIM> X = zero_vector<double>(DIM);
01265         for (unsigned node_index=0; node_index<NUM_NODES_PER_BOUNDARY_ELEMENT; node_index++)
01266         {
01267             X += quad_phi_surf_element(node_index)*rBoundaryElement.GetNode(node_index)->rGetLocation();
01268         }
01269 
01270 
01271         // Now compute the xi coordinates of the quad point in the volume element
01272         c_vector<double,DIM+1> weight = p_containing_vol_element->CalculateInterpolationWeights(X);
01273         c_vector<double,DIM> xi;
01274         for(unsigned i=0; i<DIM; i++)
01275         {
01276             xi(i) = weight(i+1); // Note, in 2d say, weights = [1-xi(0)-xi(1), xi(0), xi(1)]
01277         }
01278 
01279         // check one of the weights was zero, as the quad point is on the boundary of the volume element
01280         if(DIM==2)
01281         {
01282             assert( DIM!=2 || (fabs(weight(0))<1e-6) || (fabs(weight(1))<1e-6) || (fabs(weight(2))<1e-6) );
01283         }
01284         else
01285         {
01286             #define COVERAGE_IGNORE
01287             assert( DIM!=3 || (fabs(weight(0))<1e-6) || (fabs(weight(1))<1e-6) || (fabs(weight(2))<1e-6)  || (fabs(weight(3))<1e-6) );
01288             #undef COVERAGE_IGNORE
01289         }
01290 
01291         // Now we can compute the grad_phi and then interpolate F
01292         QuadraticBasisFunction<DIM>::ComputeTransformedBasisFunctionDerivatives(xi, inverse_jacobian_vol_element, grad_quad_phi_vol_element);
01293 
01294         F = identity_matrix<double>(DIM,DIM);
01295         for (unsigned node_index=0; node_index<NUM_NODES_PER_ELEMENT; node_index++)
01296         {
01297             for (unsigned i=0; i<DIM; i++)
01298             {
01299                 for (unsigned M=0; M<DIM; M++)
01300                 {
01301                     F(i,M) += grad_quad_phi_vol_element(M,node_index)*element_current_displacements(i,node_index);
01302                 }
01303             }
01304         }
01305 
01306         double detF = Determinant(F);
01307         invF = Inverse(F);
01308 
01309         if(assembleResidual)
01310         {
01311             c_vector<double,DIM> traction = detF*normal_pressure*prod(trans(invF),normal);
01312 
01313             // assemble
01314             for (unsigned index=0; index<NUM_NODES_PER_BOUNDARY_ELEMENT*DIM; index++)
01315             {
01316                 unsigned spatial_dim = index%DIM;
01317                 unsigned node_index = (index-spatial_dim)/DIM;
01318 
01319                 assert(node_index < NUM_NODES_PER_BOUNDARY_ELEMENT);
01320 
01321                 rBelem(index) -=   traction(spatial_dim)
01322                                  * quad_phi_surf_element(node_index)
01323                                  * wJ;
01324             }
01325         }
01326 
01327 
01328         // Sometimes we don't include the analytic jacobian for this integral
01329         // in the jacobian that is assembled - the ShouldAssembleMatrixTermForPressureOnDeformedBc()
01330         // bit below - see the documentation for this method to see why.
01331         if(assembleJacobian && ShouldAssembleMatrixTermForPressureOnDeformedBc())
01332         {
01333             for(unsigned II=0; II<NUM_NODES_PER_BOUNDARY_ELEMENT; II++)
01334             {
01335                 for(unsigned N=0; N<DIM; N++)
01336                 {
01337                     grad_quad_phi_surf_element(N,II) = grad_quad_phi_vol_element(N,surf_to_vol_map[II]);
01338                 }
01339             }
01340 
01341             static FourthOrderTensor<DIM,DIM,DIM,DIM> tensor1;
01342             for(unsigned N=0; N<DIM; N++)
01343             {
01344                 for(unsigned e=0; e<DIM; e++)
01345                 {
01346                     for(unsigned M=0; M<DIM; M++)
01347                     {
01348                         for(unsigned d=0; d<DIM; d++)
01349                         {
01350                             tensor1(N,e,M,d) = invF(N,e)*invF(M,d) - invF(M,e)*invF(N,d);
01351                         }
01352                     }
01353                 }
01354             }
01355 
01356             // tensor2(II,e,M,d) = tensor1(N,e,M,d)*grad_quad_phi_surf_element(N,II)
01357             static FourthOrderTensor<NUM_NODES_PER_BOUNDARY_ELEMENT,DIM,DIM,DIM> tensor2;
01358             tensor2.template SetAsContractionOnFirstDimension<DIM>( trans(grad_quad_phi_surf_element), tensor1);
01359 
01360             // tensor3 is really a third-order tensor
01361             // tensor3(II,e,0,d) = tensor2(II,e,M,d)*normal(M)
01362             static FourthOrderTensor<NUM_NODES_PER_BOUNDARY_ELEMENT,DIM,1,DIM> tensor3;
01363             tensor3.template SetAsContractionOnThirdDimension<DIM>( normal_as_mat, tensor2);
01364 
01365             for (unsigned index1=0; index1<NUM_NODES_PER_BOUNDARY_ELEMENT*DIM; index1++)
01366             {
01367                 unsigned spatial_dim1 = index1%DIM;
01368                 unsigned node_index1 = (index1-spatial_dim1)/DIM;
01369 
01370                 for (unsigned index2=0; index2<NUM_NODES_PER_BOUNDARY_ELEMENT*DIM; index2++)
01371                 {
01372                     unsigned spatial_dim2 = index2%DIM;
01373                     unsigned node_index2 = (index2-spatial_dim2)/DIM;
01374 
01375                     rAelem(index1,index2) -=    normal_pressure
01376                                               * detF
01377                                               * tensor3(node_index2,spatial_dim2,0,spatial_dim1)
01378                                               * quad_phi_surf_element(node_index1)
01379                                               * wJ;
01380                 }
01381             }
01382         }
01383     }
01384 }
01385 
01386 
01387 
01388 
01389 template<unsigned DIM>
01390 void AbstractNonlinearElasticitySolver<DIM>::Solve(double tol)
01391 {
01392 
01393     // check the problem definition is set up correctly (and fully).
01394     mrProblemDefinition.Validate();
01395 
01396     // If the problem includes specified pressures on deformed surfaces (as opposed
01397     // to specified tractions), the code needs to compute normals, and they need
01398     // to be consistently all facing outward (or all facing inward). Check the undeformed
01399     // mesh boundary elements has nodes that are ordered so that all normals are
01400     // outward-facing
01401     if(mrProblemDefinition.GetTractionBoundaryConditionType()==PRESSURE_ON_DEFORMED && mCheckedOutwardNormals==false)
01402     {
01403         this->mrQuadMesh.CheckOutwardNormals();
01404         mCheckedOutwardNormals = true;
01405     }
01406 
01407     // Write the initial solution
01408     this->WriteCurrentSpatialSolution("initial", "nodes");
01409 
01410     if(mUseSnesSolver)
01411     {
01412         SolveSnes();
01413     }
01414     else
01415     {
01416         SolveNonSnes(tol);
01417     }
01418 
01419     // Remove pressure dummy values (P=0 at internal nodes, which should have been
01420     // been the result of the solve above), by linear interpolating from vertices of
01421     // edges to the internal node of the edge
01422     if(this->mCompressibilityType==INCOMPRESSIBLE)
01423     {
01424         this->RemovePressureDummyValuesThroughLinearInterpolation();
01425     }
01426 
01427     // Write the final solution
01428     this->WriteCurrentSpatialSolution("solution", "nodes");
01429 
01430 }
01431 
01432 
01436 template<unsigned DIM>
01437 void AbstractNonlinearElasticitySolver<DIM>::SetKspSolverAndPcType(KSP solver)
01438 {
01439     // Three alternatives
01440     //   (a) Incompressible: GMRES with ILU preconditioner (or bjacobi=ILU on each process) [default]. Very poor on large problems.
01441     //   (b) Incompressible: GMRES with AMG preconditioner. Uncomment #define MECH_USE_HYPRE above. Requires Petsc3 with HYPRE installed.
01442     //   (c) Compressible: CG with ???
01443 
01444     PC pc;
01445     KSPGetPC(solver, &pc);
01446 
01447     if (this->mCompressibilityType==COMPRESSIBLE)
01448     {
01449         KSPSetType(solver,KSPCG);
01450         if(PetscTools::IsSequential())
01451         {
01452             PCSetType(pc, PCICC);
01453 //Note that PetscOptionsSetValue is dangerous because we can't easily do
01454 //regression testing.  If a name changes, then the behaviour of the code changes
01455 //because it won't recognise the old name.  However, it won't fail to compile/run.
01456 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 1) //PETSc 3.1 or later
01457             PetscOptionsSetValue("-pc_factor_shift_type", "positive_definite");
01458 #else
01459             PetscOptionsSetValue("-pc_factor_shift_positive_definite", "");
01460 #endif
01461         }
01462         else
01463         {
01464             PCSetType(pc, PCBJACOBI);
01465         }
01466     }
01467     else
01468     {
01469         unsigned num_restarts = 100;
01470         KSPSetType(solver,KSPGMRES);
01471         KSPGMRESSetRestart(solver,num_restarts);
01472 
01473         #ifndef MECH_USE_HYPRE
01474             PCSetType(pc, PCBJACOBI); // BJACOBI = ILU on each block (block = part of matrix on each process)
01475         #else
01476 
01477             // Speed up linear solve time massively for larger simulations (in fact GMRES may stagnate without
01478             // this for larger problems), by using a AMG preconditioner -- needs HYPRE installed
01480             PetscOptionsSetValue("-pc_hypre_type", "boomeramg");
01481             // PetscOptionsSetValue("-pc_hypre_boomeramg_max_iter", "1");
01482             // PetscOptionsSetValue("-pc_hypre_boomeramg_strong_threshold", "0.0");
01483 
01484             PCSetType(pc, PCHYPRE);
01485             KSPSetPreconditionerSide(solver, PC_RIGHT);
01486 
01487             // other possible preconditioners..
01488             //PCBlockDiagonalMechanics* p_custom_pc = new PCBlockDiagonalMechanics(solver, this->mPreconditionMatrix, mBlock1Size, mBlock2Size);
01489             //PCLDUFactorisationMechanics* p_custom_pc = new PCLDUFactorisationMechanics(solver, this->mPreconditionMatrix, mBlock1Size, mBlock2Size);
01490             //remember to delete memory..
01491             //KSPSetPreconditionerSide(solver, PC_RIGHT);
01492         #endif
01493     }
01494 }
01495 
01496 
01497 
01499 //  The code for the non-SNES solver - maybe remove all this
01500 //  as SNES solver appears better
01502 
01503 template<unsigned DIM>
01504 double AbstractNonlinearElasticitySolver<DIM>::ComputeResidualAndGetNorm(bool allowException)
01505 {
01506     if (!allowException)
01507     {
01508         // Assemble the residual
01509         AssembleSystem(true, false);
01510     }
01511     else
01512     {
01513         try
01514         {
01515             // Try to assemble the residual using this current solution
01516             AssembleSystem(true, false);
01517         }
01518         catch(Exception& e)
01519         {
01520             // If fail (because e.g. ODEs fail to solve, or strains are too large for material law), return infinity
01521             return DBL_MAX;
01522         }
01523     }
01524 
01525     // Return the scaled norm of the residual
01526     return CalculateResidualNorm();
01527 }
01528 
01529 template<unsigned DIM>
01530 double AbstractNonlinearElasticitySolver<DIM>::CalculateResidualNorm()
01531 {
01532     double norm;
01533     VecNorm(this->mResidualVector, NORM_2, &norm);
01534     return norm/this->mNumDofs;
01535 }
01536 
01537 template<unsigned DIM>
01538 void AbstractNonlinearElasticitySolver<DIM>::VectorSum(std::vector<double>& rX,
01539                                                        ReplicatableVector& rY,
01540                                                        double a,
01541                                                        std::vector<double>& rZ)
01542 {
01543     assert(rX.size()==rY.GetSize());
01544     assert(rY.GetSize()==rZ.size());
01545     for (unsigned i=0; i<rX.size(); i++)
01546     {
01547         rZ[i] = rX[i] + a*rY[i];
01548     }
01549 }
01550 
01551 
01552 template<unsigned DIM>
01553 double AbstractNonlinearElasticitySolver<DIM>::TakeNewtonStep()
01554 {
01555     if(this->mVerbose)
01556     {
01557         Timer::Reset();
01558     }
01559 
01561     // Assemble Jacobian (and preconditioner)
01563     MechanicsEventHandler::BeginEvent(MechanicsEventHandler::ASSEMBLE);
01564     AssembleSystem(true, true);
01565     MechanicsEventHandler::EndEvent(MechanicsEventHandler::ASSEMBLE);
01566     if(this->mVerbose)
01567     {
01568         Timer::PrintAndReset("AssembleSystem");
01569     }
01570 
01572     // Solve the linear system.
01574     MechanicsEventHandler::BeginEvent(MechanicsEventHandler::SOLVE);
01575 
01576     Vec solution;
01577     VecDuplicate(this->mResidualVector,&solution);
01578 
01579     KSP solver;
01580     KSPCreate(PETSC_COMM_WORLD,&solver);
01581 
01582     KSPSetOperators(solver, mrJacobianMatrix, this->mPreconditionMatrix, DIFFERENT_NONZERO_PATTERN /*in precond between successive solves*/);
01583 
01584     // Set the type of KSP solver (CG, GMRES etc) and preconditioner (ILU, HYPRE, etc)
01585     SetKspSolverAndPcType(solver);
01586 
01587     //PetscOptionsSetValue("-ksp_monitor","");
01588     //PetscOptionsSetValue("-ksp_norm_type","natural");
01589 
01590     KSPSetFromOptions(solver);
01591     KSPSetUp(solver);
01592 
01593 
01594     // Set the linear system absolute tolerance.
01595     // This is either the user provided value, or set to
01596     // max {1e-6 * initial_residual, 1e-12}
01597     if (mKspAbsoluteTol < 0)
01598     {
01599         Vec temp;
01600         VecDuplicate(this->mResidualVector, &temp);
01601         Vec temp2;
01602         VecDuplicate(this->mResidualVector, &temp2);
01603         Vec linsys_residual;
01604         VecDuplicate(this->mResidualVector, &linsys_residual);
01605 
01606         KSPInitialResidual(solver, solution, temp, temp2, linsys_residual, this->mLinearSystemRhsVector);
01607         double initial_resid_norm;
01608         VecNorm(linsys_residual, NORM_2, &initial_resid_norm);
01609 
01610         PetscTools::Destroy(temp);
01611         PetscTools::Destroy(temp2);
01612         PetscTools::Destroy(linsys_residual);
01613 
01614         double ksp_rel_tol = 1e-6;
01615         double absolute_tol = ksp_rel_tol * initial_resid_norm;
01616         if(absolute_tol < 1e-12)
01617         {
01618             absolute_tol = 1e-12;
01619         }
01620         KSPSetTolerances(solver, 1e-16, absolute_tol, PETSC_DEFAULT, 1000 /* max iters */); // Note: some machines - max iters seems to be 1000 whatever we give here
01621     }
01622     else
01623     {
01624         KSPSetTolerances(solver, 1e-16, mKspAbsoluteTol, PETSC_DEFAULT, 1000 /* max iters */); // Note: some machines - max iters seems to be 1000 whatever we give here
01625     }
01626 
01627     if(this->mVerbose)
01628     {
01629         Timer::PrintAndReset("KSP Setup");
01630     }
01631 
01632     KSPSolve(solver,this->mLinearSystemRhsVector,solution);
01633 
01634 //    ///// For printing matrix when debugging
01635 //    OutputFileHandler handler("TEMP",false);
01636 //    std::stringstream ss;
01637 //    static unsigned counter = 0;
01638 //    ss << "all_" << counter++ << ".txt";
01639 //    out_stream p_file = handler.OpenOutputFile(ss.str());
01640 //    *p_file << std::setprecision(10);
01641 //    for(unsigned i=0; i<this->mNumDofs; i++)
01642 //    {
01643 //        for(unsigned j=0; j<this->mNumDofs; j++)
01644 //        {
01645 //            *p_file << PetscMatTools::GetElement(mrJacobianMatrix, i, j) << " ";
01646 //        }
01647 //        *p_file << PetscVecTools::GetElement(this->mLinearSystemRhsVector, i) << " ";
01648 //        *p_file << PetscVecTools::GetElement(solution, i) << "\n";
01649 //    }
01650 //    p_file->close();
01651 
01652 
01654     // Error checking for linear solve
01656 
01657     // warn if ksp reports failure
01658     KSPConvergedReason reason;
01659     KSPGetConvergedReason(solver,&reason);
01660 
01661     if(reason != KSP_DIVERGED_ITS)
01662     {
01663         // Throw an exception if the solver failed for any reason other than DIVERGED_ITS.
01664         // This is not covered as would be difficult to cover - requires a bad matrix to
01665         // assembled, for example.
01666         #define COVERAGE_IGNORE
01667         KSPEXCEPT(reason);
01668         #undef COVERAGE_IGNORE
01669     }
01670     else
01671     {
01672         // DIVERGED_ITS just means it didn't converge in the given maximum number of iterations,
01673         // which is potentially not a problem, as the nonlinear solver may (and often will) still converge.
01674         // Just warn once.
01675         // (Very difficult to cover in normal tests, requires relative and absolute ksp tols to be very small, there
01676         // is no interface for setting both of these. Could be covered by setting up a problem the solver
01677         // finds difficult to solve, but this would be overkill.)
01678         #define COVERAGE_IGNORE
01679         WARN_ONCE_ONLY("Linear solve (within a mechanics solve) didn't converge, but this may not stop nonlinear solve converging")
01680         #undef COVERAGE_IGNORE
01681     }
01682 
01683     // quit if no ksp iterations were done
01684     int num_iters;
01685     KSPGetIterationNumber(solver, &num_iters);
01686     if (num_iters==0)
01687     {
01688         PetscTools::Destroy(solution);
01689         KSPDestroy(PETSC_DESTROY_PARAM(solver));
01690         EXCEPTION("KSP Absolute tolerance was too high, linear system wasn't solved - there will be no decrease in Newton residual. Decrease KspAbsoluteTolerance");
01691     }
01692 
01693 
01694     if(this->mVerbose)
01695     {
01696         Timer::PrintAndReset("KSP Solve");
01697         std::cout << "[" << PetscTools::GetMyRank() << "]: Num iterations = " << num_iters << "\n" << std::flush;
01698     }
01699 
01700     MechanicsEventHandler::EndEvent(MechanicsEventHandler::SOLVE);
01701 
01703     // Update the solution
01704     //  Newton method:       sol = sol - update, where update=Jac^{-1}*residual
01705     //  Newton with damping: sol = sol - s*update, where s is chosen
01706     //   such that |residual(sol)| is minimised. Damping is important to
01707     //   avoid initial divergence.
01708     //
01709     // Normally, finding the best s from say 0.05,0.1,0.2,..,1.0 is cheap,
01710     // but this is not the case in cardiac electromechanics calculations.
01711     // Therefore, we initially check s=1 (expected to be the best most of the
01712     // time, then s=0.9. If the norm of the residual increases, we assume
01713     // s=1 is the best. Otherwise, check s=0.8 to see if s=0.9 is a local min.
01715     MechanicsEventHandler::BeginEvent(MechanicsEventHandler::UPDATE);
01716     double new_norm_resid = UpdateSolutionUsingLineSearch(solution);
01717     MechanicsEventHandler::EndEvent(MechanicsEventHandler::UPDATE);
01718 
01719     PetscTools::Destroy(solution);
01720     KSPDestroy(PETSC_DESTROY_PARAM(solver));
01721 
01722     return new_norm_resid;
01723 }
01724 
01725 
01726 template<unsigned DIM>
01727 void AbstractNonlinearElasticitySolver<DIM>::PrintLineSearchResult(double s, double residNorm)
01728 {
01729     if(this->mVerbose)
01730     {
01731         std::cout << "\tTesting s = " << s << ", |f| = " << residNorm << "\n" << std::flush;
01732     }
01733 }
01734 
01735 template<unsigned DIM>
01736 double AbstractNonlinearElasticitySolver<DIM>::UpdateSolutionUsingLineSearch(Vec solution)
01737 {
01738     double initial_norm_resid = CalculateResidualNorm();
01739     if(this->mVerbose)
01740     {
01741         std::cout << "\tInitial |f| [corresponding to s=0] is " << initial_norm_resid << "\n"  << std::flush;
01742     }
01743 
01744     ReplicatableVector update(solution);
01745 
01746     std::vector<double> old_solution = this->mCurrentSolution;
01747 
01748     std::vector<double> damping_values; // = {1.0, 0.9, .., 0.2, 0.1, 0.05} ie size 11
01749     for (unsigned i=10; i>=1; i--)
01750     {
01751         damping_values.push_back((double)i/10.0);
01752     }
01753     damping_values.push_back(0.05);
01754     assert(damping_values.size()==11);
01755 
01757     // let mCurrentSolution = old_solution - damping_val[0]*update; and compute residual
01758     unsigned index = 0;
01759     VectorSum(old_solution, update, -damping_values[index], this->mCurrentSolution);
01760     double current_resid_norm = ComputeResidualAndGetNorm(true);
01761     PrintLineSearchResult(damping_values[index], current_resid_norm);
01762 
01764     // let mCurrentSolution = old_solution - damping_val[1]*update; and compute residual
01765     index = 1;
01766     VectorSum(old_solution, update, -damping_values[index], this->mCurrentSolution);
01767     double next_resid_norm = ComputeResidualAndGetNorm(true);
01768     PrintLineSearchResult(damping_values[index], next_resid_norm);
01769 
01770     index = 2;
01771     // While f(s_next) < f(s_current), [f = residnorm], keep trying new damping values,
01772     // ie exit thus loop when next norm of the residual first increases
01773     while (    (next_resid_norm==DBL_MAX) // the residual is returned as infinity if the deformation is so large to cause exceptions in the material law/EM contraction model
01774             || ( (next_resid_norm < current_resid_norm) && (index<damping_values.size()) ) )
01775     {
01776         current_resid_norm = next_resid_norm;
01777 
01778         // let mCurrentSolution = old_solution - damping_val*update; and compute residual
01779         VectorSum(old_solution, update, -damping_values[index], this->mCurrentSolution);
01780         next_resid_norm = ComputeResidualAndGetNorm(true);
01781         PrintLineSearchResult(damping_values[index], next_resid_norm);
01782 
01783         index++;
01784     }
01785 
01786     unsigned best_index;
01787 
01788     if (index==damping_values.size() && (next_resid_norm < current_resid_norm))
01789     {
01790         // Difficult to come up with large forces/tractions such that it had to
01791         // test right down to s=0.05, but overall doesn't fail.
01792         // The possible damping values have been manually temporarily altered to
01793         // get this code to be called, it appears to work correctly. Even if it
01794         // didn't tests wouldn't fail, they would just be v. slightly less efficient.
01795         #define COVERAGE_IGNORE
01796         // if we exited because we got to the end of the possible damping values, the
01797         // best one was the last one (excl the final index++ at the end)
01798         current_resid_norm = next_resid_norm;
01799         best_index = index-1;
01800         #undef COVERAGE_IGNORE
01801     }
01802     else
01803     {
01804         // else the best one must have been the second last one (excl the final index++ at the end)
01805         // (as we would have exited when the resid norm first increased)
01806         best_index = index-2;
01807     }
01808 
01809     // Check out best was better than the original residual-norm
01810     if (initial_norm_resid < current_resid_norm)
01811     {
01812         #define COVERAGE_IGNORE
01813         // Have to use an assert/exit here as the following exception causes a seg fault (in cardiac mech problems?)
01814         // Don't know why
01815         std::cout << "CHASTE ERROR: (AbstractNonlinearElasticitySolver.hpp): Residual does not appear to decrease in newton direction, quitting.\n" << std::flush;
01816         EXCEPTION("Residual does not appear to decrease in newton direction, quitting");
01817         #undef COVERAGE_IGNORE
01818     }
01819 
01820     if(this->mVerbose)
01821     {
01822         std::cout << "\tBest s = " << damping_values[best_index] << "\n"  << std::flush;
01823     }
01824 
01825     VectorSum(old_solution, update, -damping_values[best_index], this->mCurrentSolution);
01826 
01827     mLastDampingValue = damping_values[best_index];
01828     return current_resid_norm;
01829 }
01830 
01831 template<unsigned DIM>
01832 void AbstractNonlinearElasticitySolver<DIM>::PostNewtonStep(unsigned counter, double normResidual)
01833 {
01834 }
01835 
01836 
01837 template<unsigned DIM>
01838 void AbstractNonlinearElasticitySolver<DIM>::SolveNonSnes(double tol)
01839 {
01840     mLastDampingValue = 0;
01841 
01842     if (mWriteOutputEachNewtonIteration)
01843     {
01844         this->WriteCurrentSpatialSolution("newton_iteration", "nodes", 0);
01845     }
01846 
01847     // Compute residual
01848     double norm_resid = ComputeResidualAndGetNorm(false);
01849     if(this->mVerbose)
01850     {
01851         std::cout << "\nNorm of residual is " << norm_resid << "\n";
01852     }
01853 
01854     mNumNewtonIterations = 0;
01855     unsigned iteration_number = 1;
01856 
01857     if (tol < 0) // i.e. if wasn't passed in as a parameter
01858     {
01859         tol = NEWTON_REL_TOL*norm_resid;
01860 
01861         #define COVERAGE_IGNORE // not going to have tests in cts for everything
01862         if (tol > MAX_NEWTON_ABS_TOL)
01863         {
01864             tol = MAX_NEWTON_ABS_TOL;
01865         }
01866         if (tol < MIN_NEWTON_ABS_TOL)
01867         {
01868             tol = MIN_NEWTON_ABS_TOL;
01869         }
01870         #undef COVERAGE_IGNORE
01871     }
01872 
01873     if(this->mVerbose)
01874     {
01875         std::cout << "Solving with tolerance " << tol << "\n";
01876     }
01877 
01878     while (norm_resid > tol)
01879     {
01880         if(this->mVerbose)
01881         {
01882             std::cout <<  "\n-------------------\n"
01883                       <<   "Newton iteration " << iteration_number
01884                       << ":\n-------------------\n";
01885         }
01886 
01887         // take newton step (and get returned residual)
01888         norm_resid = TakeNewtonStep();
01889 
01890         if(this->mVerbose)
01891         {
01892             std::cout << "Norm of residual is " << norm_resid << "\n";
01893         }
01894 
01895         if (mWriteOutputEachNewtonIteration)
01896         {
01897             this->WriteCurrentSpatialSolution("newton_iteration", "nodes", iteration_number);
01898         }
01899 
01900         mNumNewtonIterations = iteration_number;
01901 
01902         PostNewtonStep(iteration_number,norm_resid);
01903 
01904         iteration_number++;
01905         if (iteration_number==20)
01906         {
01907             #define COVERAGE_IGNORE
01908             EXCEPTION("Not converged after 20 newton iterations, quitting");
01909             #undef COVERAGE_IGNORE
01910         }
01911     }
01912 
01913     if (norm_resid > tol)
01914     {
01915         #define COVERAGE_IGNORE
01916         EXCEPTION("Failed to converge");
01917         #undef COVERAGE_IGNORE
01918     }
01919 }
01920 
01921 
01922 
01923 template<unsigned DIM>
01924 unsigned AbstractNonlinearElasticitySolver<DIM>::GetNumNewtonIterations()
01925 {
01926     return mNumNewtonIterations;
01927 }
01928 
01929 
01930 
01932 //  SNES version of the nonlinear solver
01934 
01935 
01936 template<unsigned DIM>
01937 void AbstractNonlinearElasticitySolver<DIM>::SolveSnes()
01938 {
01939     // Set up solution guess for residuals
01940     Vec initial_guess;
01941     VecDuplicate(this->mResidualVector, &initial_guess);
01942     double* p_initial_guess;
01943     VecGetArray(initial_guess, &p_initial_guess);
01944     int lo, hi;
01945     VecGetOwnershipRange(initial_guess, &lo, &hi);
01946     for (int global_index=lo; global_index<hi; global_index++)
01947     {
01948         int local_index = global_index - lo;
01949         p_initial_guess[local_index] = this->mCurrentSolution[global_index];
01950     }
01951     VecRestoreArray(initial_guess, &p_initial_guess);
01952     PetscVecTools::Finalise(initial_guess);
01953 
01954     Vec snes_residual_vec;
01955     VecDuplicate(this->mResidualVector, &snes_residual_vec);
01956 
01957     SNES snes;
01958 
01959     SNESCreate(PETSC_COMM_WORLD, &snes);
01960     SNESSetFunction(snes, snes_residual_vec, &AbstractNonlinearElasticitySolver_ComputeResidual<DIM>, this);
01961     SNESSetJacobian(snes, mrJacobianMatrix, this->mPreconditionMatrix, &AbstractNonlinearElasticitySolver_ComputeJacobian<DIM>, this);
01962     SNESSetType(snes,SNESLS);
01963     SNESSetTolerances(snes,1e-5,1e-5,1e-5,PETSC_DEFAULT,PETSC_DEFAULT);
01964     SNESSetMaxLinearSolveFailures(snes,100);
01965 
01966     KSP ksp;
01967     SNESGetKSP(snes, &ksp);
01968 
01969     KSPSetTolerances(ksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 1000 /* max iters */); // Note: some machines - max iters seems to be 1000 whatever we give here
01970 
01971     // Set the type of KSP solver (CG, GMRES etc) and preconditioner (ILU, HYPRE, etc)
01972     SetKspSolverAndPcType(ksp);
01973 
01974     if(this->mVerbose)
01975     {
01976         PetscOptionsSetValue("-snes_monitor","");
01977     }
01978     SNESSetFromOptions(snes);
01979 
01980 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
01981     //int ierr =
01982     SNESSolve(snes, initial_guess);
01983 #else
01984     SNESSolve(snes, PETSC_NULL, initial_guess);
01985 #endif
01986 
01987     SNESConvergedReason reason;
01988     SNESGetConvergedReason(snes,&reason);
01989 #define COVERAGE_IGNORE
01990     if (reason < 0)
01991     {
01992         std::stringstream reason_stream;
01993         reason_stream << reason;
01994         PetscTools::Destroy(initial_guess);
01995         PetscTools::Destroy(snes_residual_vec);
01996         SNESDestroy(PETSC_DESTROY_PARAM(snes));
01997         EXCEPTION("Nonlinear Solver did not converge. PETSc reason code:"+reason_stream.str()+" .");
01998     }
01999 #undef COVERAGE_IGNORE
02000 
02001     PetscInt num_iters;
02002     SNESGetIterationNumber(snes,&num_iters);
02003     mNumNewtonIterations = num_iters;
02004 
02005     PetscTools::Destroy(initial_guess);
02006     PetscTools::Destroy(snes_residual_vec);
02007     SNESDestroy(PETSC_DESTROY_PARAM(snes));
02008 }
02009 
02010 
02011 template<unsigned DIM>
02012 void AbstractNonlinearElasticitySolver<DIM>::ComputeResidual(Vec currentGuess, Vec residualVector)
02013 {
02014     // Note: AssembleSystem() assumes the current solution is in this->mCurrentSolution and assembles
02015     // this->mResiduaVector and/or this->mrJacobianMatrix. Since PETSc wants us to use the input
02016     // currentGuess, and write the output to residualVector, we have to copy do some copies below.
02017 
02018     ReplicatableVector guess_repl(currentGuess);
02019     for(unsigned i=0; i<guess_repl.GetSize(); i++)
02020     {
02021         this->mCurrentSolution[i] = guess_repl[i];
02022     }
02023     AssembleSystem(true,false);
02024     VecCopy(this->mResidualVector, residualVector);
02025 }
02026 
02027 template<unsigned DIM>
02028 void AbstractNonlinearElasticitySolver<DIM>::ComputeJacobian(Vec currentGuess, Mat* pJacobian, Mat* pPreconditioner)
02029 {
02030     // Note: AssembleSystem() assumes the current solution is in this->mCurrentSolution and assembles
02031     // this->mResiduaVector and/or this->mrJacobianMatrix.
02032     // We need to copy the input currentGuess into the local mCurrentGuess.
02033     // We don't have to copy mrJacobianMatrix to pJacobian, which would be expensive, as they will
02034     // point to the same memory.
02035 
02036     // check Petsc data corresponds to internal Mats
02037     assert(mrJacobianMatrix==*pJacobian);
02038     assert(this->mPreconditionMatrix==*pPreconditioner);
02039 
02040     MechanicsEventHandler::BeginEvent(MechanicsEventHandler::ASSEMBLE);
02041     ReplicatableVector guess_repl(currentGuess);
02042     for(unsigned i=0; i<guess_repl.GetSize(); i++)
02043     {
02044         this->mCurrentSolution[i] = guess_repl[i];
02045     }
02046 
02047     AssembleSystem(false,true);
02048     MechanicsEventHandler::EndEvent(MechanicsEventHandler::ASSEMBLE);
02049 }
02050 
02051 
02052 
02053 template<unsigned DIM>
02054 PetscErrorCode AbstractNonlinearElasticitySolver_ComputeResidual(SNES snes,
02055                                                                  Vec currentGuess,
02056                                                                  Vec residualVector,
02057                                                                  void* pContext)
02058 {
02059     // Extract the solver from the void*
02060     AbstractNonlinearElasticitySolver<DIM>* p_solver = (AbstractNonlinearElasticitySolver<DIM>*)pContext;
02061     p_solver->ComputeResidual(currentGuess, residualVector);
02062     return 0;
02063 }
02064 
02065 template<unsigned DIM>
02066 PetscErrorCode AbstractNonlinearElasticitySolver_ComputeJacobian(SNES snes,
02067                                                                  Vec currentGuess,
02068                                                                  Mat* pGlobalJacobian,
02069                                                                  Mat* pPreconditioner,
02070                                                                  MatStructure* pMatStructure,
02071                                                                  void* pContext)
02072 {
02073     // Extract the solver from the void*
02074     AbstractNonlinearElasticitySolver<DIM>* p_solver = (AbstractNonlinearElasticitySolver<DIM>*) pContext;
02075     p_solver->ComputeJacobian(currentGuess, pGlobalJacobian, pPreconditioner);
02076     return 0;
02077 }
02078 
02079 
02080 // Constant setting definitions
02081 
02082 template<unsigned DIM>
02083 double AbstractNonlinearElasticitySolver<DIM>::MAX_NEWTON_ABS_TOL = 1e-7;
02084 
02085 template<unsigned DIM>
02086 double AbstractNonlinearElasticitySolver<DIM>::MIN_NEWTON_ABS_TOL = 1e-10;
02087 
02088 template<unsigned DIM>
02089 double AbstractNonlinearElasticitySolver<DIM>::NEWTON_REL_TOL = 1e-4;
02090 
02091 #endif /*ABSTRACTNONLINEARELASTICITYSOLVER_HPP_*/