AbstractNonlinearElasticitySolver.hpp

00001 /*
00002 
00003 Copyright (C) University of Oxford, 2005-2011
00004 
00005 University of Oxford means the Chancellor, Masters and Scholars of the
00006 University of Oxford, having an administrative office at Wellington
00007 Square, Oxford OX1 2JD, UK.
00008 
00009 This file is part of Chaste.
00010 
00011 Chaste is free software: you can redistribute it and/or modify it
00012 under the terms of the GNU Lesser General Public License as published
00013 by the Free Software Foundation, either version 2.1 of the License, or
00014 (at your option) any later version.
00015 
00016 Chaste is distributed in the hope that it will be useful, but WITHOUT
00017 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018 FITNESS FOR A PARTICULAR PURPOSE.  See the GNU Lesser General Public
00019 License for more details. The offer of Chaste under the terms of the
00020 License is subject to the License being interpreted in accordance with
00021 English Law and subject to any action against the University of Oxford
00022 being under the jurisdiction of the English Courts.
00023 
00024 You should have received a copy of the GNU Lesser General Public License
00025 along with Chaste. If not, see <http://www.gnu.org/licenses/>.
00026 
00027 */
00028 
00029 #ifndef ABSTRACTNONLINEARELASTICITYSOLVER_HPP_
00030 #define ABSTRACTNONLINEARELASTICITYSOLVER_HPP_
00031 
00032 #include <vector>
00033 #include <cmath>
00034 #include "LinearSystem.hpp"
00035 #include "OutputFileHandler.hpp"
00036 #include "LogFile.hpp"
00037 #include "PetscTools.hpp"
00038 #include "MechanicsEventHandler.hpp"
00039 #include "ReplicatableVector.hpp"
00040 #include "FourthOrderTensor.hpp"
00041 #include "QuadraticMesh.hpp"
00042 #include "GaussianQuadratureRule.hpp"
00043 #include "CmguiDeformedSolutionsWriter.hpp"
00044 
00045 //#include "PCBlockDiagonalMechanics.hpp"
00046 //#include "PCLDUFactorisationMechanics.hpp"
00047 
00048 
00049 //#define MECH_VERBOSE
00050 //#define MECH_VERY_VERBOSE
00051 
00053 //#define MECH_USE_MUMPS  // requires Petsc to be installed with MUMPS
00054 //#define MECH_USE_HYPRE  // requires Petsc to be installed with HYPRE
00055 
00056 
00057 
00058 #ifdef MECH_VERBOSE
00059 #include "Timer.hpp"
00060 #endif
00061 
00062 typedef enum CompressibilityType_
00063 {
00064     COMPRESSIBLE,
00065     INCOMPRESSIBLE
00066 } CompressibilityType;
00067 
00068 
00072 template<CompressibilityType COMPRESSIBILITY_TYPE, unsigned DIM>
00073 class AbstractNonlinearElasticitySolver
00074 {
00075 protected:
00076 
00082     static double MAX_NEWTON_ABS_TOL;
00083 
00085     static double MIN_NEWTON_ABS_TOL;
00086 
00088     static double NEWTON_REL_TOL;
00089 
00090 
00095     QuadraticMesh<DIM>* mpQuadMesh;
00096 
00098     std::vector<BoundaryElement<DIM-1,DIM>*> mBoundaryElements;
00099 
00101     GaussianQuadratureRule<DIM>* mpQuadratureRule;
00102 
00104     GaussianQuadratureRule<DIM-1>* mpBoundaryQuadratureRule;
00105     
00106 
00110     double mKspAbsoluteTol;
00111 
00117     unsigned mNumDofs;
00118 
00124     LinearSystem* mpLinearSystem;
00125 
00126 
00147     LinearSystem* mpPreconditionMatrixLinearSystem;
00148 
00150     c_vector<double,DIM> mBodyForce;
00151 
00153     double mDensity;
00154 
00156     std::vector<unsigned> mFixedNodes;
00157 
00159     std::vector<c_vector<double,DIM> > mFixedNodeDisplacements;
00160 
00162     bool mWriteOutput;
00163 
00165     std::string mOutputDirectory;
00166 
00168     OutputFileHandler* mpOutputFileHandler;
00169 
00173     bool mWriteOutputEachNewtonIteration;
00174 
00175 
00182     std::vector<double> mCurrentSolution;
00183 
00187     FourthOrderTensor<DIM,DIM,DIM,DIM> dTdE;
00188 
00190     unsigned mNumNewtonIterations;
00191 
00193     std::vector<c_vector<double,DIM> > mDeformedPosition;
00194 
00195 
00199     std::vector<c_vector<double,DIM> > mSurfaceTractions;
00200 
00202     c_vector<double,DIM> (*mpBodyForceFunction)(c_vector<double,DIM>& X, double t);
00203 
00208     c_vector<double,DIM> (*mpTractionBoundaryConditionFunction)(c_vector<double,DIM>& X, double t);
00209 
00211     bool mUsingBodyForceFunction;
00212 
00214     bool mUsingTractionBoundaryConditionFunction;
00215     
00220     double mCurrentTime;
00221 
00232     virtual void AssembleSystem(bool assembleResidual, bool assembleJacobian)=0;
00233 
00234 
00235 
00241     void Initialise(std::vector<c_vector<double,DIM> >* pFixedNodeLocations);
00242 
00243 
00247     void AllocateMatrixMemory();
00248 
00249 
00255     void ApplyBoundaryConditions(bool applyToMatrix);
00256 
00267     double ComputeResidualAndGetNorm(bool allowException);
00268 
00270     double CalculateResidualNorm();
00271 
00279     void VectorSum(std::vector<double>& rX, ReplicatableVector& rY, double a, std::vector<double>& rZ);
00280 
00287     void PrintLineSearchResult(double s, double residNorm);
00288 
00289 
00297     double TakeNewtonStep();
00298 
00306     double UpdateSolutionUsingLineSearch(Vec solution);
00307 
00308 
00316     virtual void PostNewtonStep(unsigned counter, double normResidual);
00317 
00318 public:
00319 
00329     AbstractNonlinearElasticitySolver(QuadraticMesh<DIM>* pQuadMesh,
00330                                       c_vector<double,DIM> bodyForce,
00331                                       double density,
00332                                       std::string outputDirectory,
00333                                       std::vector<unsigned>& fixedNodes);
00334 
00335 
00339     virtual ~AbstractNonlinearElasticitySolver();
00340 
00348     void Solve(double tol=-1.0,
00349                unsigned maxNumNewtonIterations=INT_MAX,
00350                bool quitIfNoConvergence=true);
00351 
00362     void WriteCurrentDeformation(std::string fileName, int counterToAppend=-1);
00363 
00367     unsigned GetNumNewtonIterations();
00368 
00376     void SetFunctionalBodyForce(c_vector<double,DIM> (*pFunction)(c_vector<double,DIM>& X, double t));
00377 
00383     void SetWriteOutput(bool writeOutput=true);
00384 
00391     void SetWriteOutputEachNewtonIteration(bool writeOutputEachNewtonIteration=true)
00392     {
00393         mWriteOutputEachNewtonIteration = writeOutputEachNewtonIteration;
00394     }
00395 
00401     void SetKspAbsoluteTolerance(double kspAbsoluteTolerance)
00402     {
00403         assert(kspAbsoluteTolerance>0);
00404         mKspAbsoluteTol = kspAbsoluteTolerance;
00405     }
00406 
00411     std::vector<double>& rGetCurrentSolution()
00412     {
00413         return mCurrentSolution;
00414     }
00415 
00424     void SetSurfaceTractionBoundaryConditions(std::vector<BoundaryElement<DIM-1,DIM>*>& rBoundaryElements,
00425                                               std::vector<c_vector<double,DIM> >& rSurfaceTractions);
00426 
00435     void SetFunctionalTractionBoundaryCondition(std::vector<BoundaryElement<DIM-1,DIM>*> rBoundaryElements,
00436                                                 c_vector<double,DIM> (*pFunction)(c_vector<double,DIM>& X, double t));
00437 
00441     std::vector<c_vector<double,DIM> >& rGetDeformedPosition();
00442 
00447     void SetCurrentTime(double time)
00448     {
00449         mCurrentTime = time;
00450     }
00451 
00455     void CreateCmguiOutput();
00456 };
00457 
00458 
00460 // Implementation
00462 
00463 
00464 template<CompressibilityType COMPRESSIBILITY_TYPE, unsigned DIM>
00465 void AbstractNonlinearElasticitySolver<COMPRESSIBILITY_TYPE,DIM>::Initialise(std::vector<c_vector<double,DIM> >* pFixedNodeLocations)
00466 {
00467     assert(mpQuadMesh);
00468 
00469     AllocateMatrixMemory();
00470 
00471     for (unsigned i=0; i<mFixedNodes.size(); i++)
00472     {
00473         assert(mFixedNodes[i] < mpQuadMesh->GetNumNodes());
00474     }
00475 
00476     mpQuadratureRule = new GaussianQuadratureRule<DIM>(3);
00477     mpBoundaryQuadratureRule = new GaussianQuadratureRule<DIM-1>(3);
00478 
00479     mCurrentSolution.resize(mNumDofs, 0.0);
00480 
00481     // compute the displacements at each of the fixed nodes, given the
00482     // fixed nodes locations.
00483     if (pFixedNodeLocations == NULL)
00484     {
00485         mFixedNodeDisplacements.clear();
00486         for (unsigned i=0; i<mFixedNodes.size(); i++)
00487         {
00488             mFixedNodeDisplacements.push_back(zero_vector<double>(DIM));
00489         }
00490     }
00491     else
00492     {
00493         assert(pFixedNodeLocations->size()==mFixedNodes.size());
00494         for (unsigned i=0; i<mFixedNodes.size(); i++)
00495         {
00496             unsigned index = mFixedNodes[i];
00497             c_vector<double,DIM> displacement = (*pFixedNodeLocations)[i] - mpQuadMesh->GetNode(index)->rGetLocation();
00498             mFixedNodeDisplacements.push_back(displacement);
00499         }
00500     }
00501     assert(mFixedNodeDisplacements.size()==mFixedNodes.size());
00502 }
00503 
00504 
00505 
00506 template<CompressibilityType COMPRESSIBILITY_TYPE, unsigned DIM>
00507 void AbstractNonlinearElasticitySolver<COMPRESSIBILITY_TYPE,DIM>::AllocateMatrixMemory()
00508 {
00509     if(DIM==2)
00510     {
00511         // 2D: N elements around a point => 7N+3 non-zeros in that row? Assume N<=10 (structured mesh would have N_max=6) => 73.
00512         unsigned num_non_zeros = std::min(75u, mNumDofs);
00513 
00514         mpLinearSystem = new LinearSystem(mNumDofs, num_non_zeros);
00515         mpPreconditionMatrixLinearSystem = new LinearSystem(mNumDofs, num_non_zeros);
00516     }
00517     else
00518     {
00519         assert(DIM==3);
00520 
00521         // pass in 0 as the preallocation number
00522 
00523         // Note: linear system always tries to use MAT_IGNORE_OFF_PROC_ENTRIES,
00524         // which PetscTools::SetupMat() can do with num_dofs = 0, hence the warnings
00525         // that get output. To fix this: replace linear systems with matrices, calling
00526         // PetscTools::SetupMat(), and in the future, when parallelising, remember to
00527         // think about MAT_IGNORE_OFF_PROC_ENTRIES
00528 
00529         mpLinearSystem = new LinearSystem(mNumDofs, 0);
00530         mpPreconditionMatrixLinearSystem = new LinearSystem(mNumDofs, 0);
00531 
00532         // 3D: N elements around a point. nz < (3*10+6)N (lazy estimate). Better estimate is 23N+4?. Assume N<20 => 500ish
00533 
00534         // in 3d we get the number of containing elements for each node and use that to obtain an upper bound
00535         // for the number of non-zeros for each DOF associated with that node.
00536 
00537         int* num_non_zeros_each_row = new int[mNumDofs];
00538         for(unsigned i=0; i<mNumDofs; i++)
00539         {
00540             num_non_zeros_each_row[i] = 0;
00541         }
00542 
00543         for(unsigned i=0; i<mpQuadMesh->GetNumNodes(); i++)
00544         {
00545             // this upper bound neglects the fact that two containing elements will share the same nodes..
00546             // 4 = max num dofs associated with this node
00547             // 30 = 3*9+3 = 3 dimensions x 9 other nodes on this element   +  3 vertices with a pressure unknown
00548             unsigned num_non_zeros_upper_bound = 4 + 30*mpQuadMesh->GetNode(i)->GetNumContainingElements();
00549 
00550             num_non_zeros_upper_bound = std::min(num_non_zeros_upper_bound, mNumDofs);
00551 
00552             num_non_zeros_each_row[DIM*i + 0] = num_non_zeros_upper_bound;
00553             num_non_zeros_each_row[DIM*i + 1] = num_non_zeros_upper_bound;
00554             num_non_zeros_each_row[DIM*i + 2] = num_non_zeros_upper_bound;
00555 
00556             if(COMPRESSIBILITY_TYPE==INCOMPRESSIBLE)
00557             {
00558                 //Could do !mpQuadMesh->GetNode(i)->IsInternal()
00559                 if(i<mpQuadMesh->GetNumVertices()) // then this is a vertex
00560                 {
00561                     num_non_zeros_each_row[DIM*mpQuadMesh->GetNumNodes() + i] = num_non_zeros_upper_bound;
00562                 }
00563             }
00564         }
00565 
00566         // NOTE: PetscTools::SetupMat() creates a MATAIJ matrix, which means the matrix will
00567         // be of type MATSEQAIJ if num_procs=1 and MATMPIAIJ otherwise. In the former case
00568         // MatSeqAIJSetPreallocation MUST be called [MatMPIAIJSetPreallocation will have
00569         // no effect (silently)], and vice versa in the latter case
00570 
00571         if(PetscTools::IsSequential())
00572         {
00573             MatSeqAIJSetPreallocation(mpLinearSystem->rGetLhsMatrix(),                   PETSC_NULL, num_non_zeros_each_row);
00574             MatSeqAIJSetPreallocation(mpPreconditionMatrixLinearSystem->rGetLhsMatrix(), PETSC_NULL, num_non_zeros_each_row);
00575         }
00576         else
00577         {
00578             PetscInt lo, hi;
00579             mpLinearSystem->GetOwnershipRange(lo, hi);
00580             int* num_non_zeros_each_row_this_proc = new int[hi-lo];
00581             int* zero = new int[hi-lo];
00582             for(unsigned i=0; i<unsigned(hi-lo); i++)
00583             {
00584                 num_non_zeros_each_row_this_proc[i] = num_non_zeros_each_row[lo+i];
00585                 zero[i] = 0;
00586             }
00587 
00588             MatMPIAIJSetPreallocation(mpLinearSystem->rGetLhsMatrix(), PETSC_NULL, num_non_zeros_each_row_this_proc, PETSC_NULL, num_non_zeros_each_row_this_proc);
00589             MatMPIAIJSetPreallocation(mpPreconditionMatrixLinearSystem->rGetLhsMatrix(), PETSC_NULL, num_non_zeros_each_row_this_proc, PETSC_NULL, num_non_zeros_each_row_this_proc);
00590         }
00591 
00592         //unsigned total_non_zeros = 0;
00593         //for(unsigned i=0; i<mNumDofs; i++)
00594         //{
00595         //   total_non_zeros += num_non_zeros_each_row[i];
00596         //}
00597         //std::cout << total_non_zeros << " versus " << 500*mNumDofs << "\n" << std::flush;
00598 
00599         delete [] num_non_zeros_each_row;
00600     }
00601 }
00602 
00603 
00604 
00605 template<CompressibilityType COMPRESSIBILITY_TYPE, unsigned DIM>
00606 void AbstractNonlinearElasticitySolver<COMPRESSIBILITY_TYPE,DIM>::ApplyBoundaryConditions(bool applyToMatrix)
00607 {
00608     assert(mFixedNodeDisplacements.size()==mFixedNodes.size());
00609 
00610     // The boundary conditions on the NONLINEAR SYSTEM are x=boundary_values
00611     // on the boundary nodes. However:
00612     // The boundary conditions on the LINEAR SYSTEM  Ju=f, where J is the
00613     // u the negative update vector and f is the residual is
00614     // u=current_soln-boundary_values on the boundary nodes
00615 
00616     std::vector<unsigned> rows;
00617     if(applyToMatrix)
00618     {
00619         rows.resize(DIM*mFixedNodes.size());
00620     }
00621 
00622     for (unsigned i=0; i<mFixedNodes.size(); i++)
00623     {
00624         unsigned node_index = mFixedNodes[i];
00625         for (unsigned j=0; j<DIM; j++)
00626         {
00627             unsigned dof_index = DIM*node_index+j;
00628 
00629             if(applyToMatrix)
00630             {
00631                 rows[DIM*i+j] = dof_index;
00632             }
00633 
00634             double value = mCurrentSolution[dof_index] - mFixedNodeDisplacements[i](j);
00635             mpLinearSystem->SetRhsVectorElement(dof_index, value);
00636         }
00637     }
00638 
00639     if(applyToMatrix)
00640     {
00641         mpLinearSystem->ZeroMatrixRowsWithValueOnDiagonal(rows, 1.0);
00642         mpPreconditionMatrixLinearSystem->ZeroMatrixRowsWithValueOnDiagonal(rows, 1.0);
00643     }
00644 }
00645 
00646 template<CompressibilityType COMPRESSIBILITY_TYPE, unsigned DIM>
00647 double AbstractNonlinearElasticitySolver<COMPRESSIBILITY_TYPE,DIM>::ComputeResidualAndGetNorm(bool allowException)
00648 {
00649     if(!allowException)
00650     {
00651         // assemble the residual
00652         AssembleSystem(true, false);
00653     }
00654     else
00655     {
00656         try
00657         {
00658             // try to assemble the residual using this current solution
00659             AssembleSystem(true, false);
00660         }
00661         catch(Exception& e)
00662         {
00663             // if fail (because eg ODEs fail to solve, or strains are too large for material law), 
00664             // return infinity
00665             return DBL_MAX;
00666         }
00667     }
00668 
00669     // return the scaled norm of the residual
00670     return CalculateResidualNorm();
00671 }
00672 
00673 template<CompressibilityType COMPRESSIBILITY_TYPE, unsigned DIM>
00674 double AbstractNonlinearElasticitySolver<COMPRESSIBILITY_TYPE,DIM>::CalculateResidualNorm()
00675 {
00676     double norm;
00677     VecNorm(mpLinearSystem->rGetRhsVector(), NORM_2, &norm);
00678     return norm/mNumDofs;
00679 }
00680 
00681 template<CompressibilityType COMPRESSIBILITY_TYPE, unsigned DIM>
00682 void AbstractNonlinearElasticitySolver<COMPRESSIBILITY_TYPE,DIM>::VectorSum(std::vector<double>& rX,
00683                                                        ReplicatableVector& rY,
00684                                                        double a,
00685                                                        std::vector<double>& rZ)
00686 {
00687     assert(rX.size()==rY.GetSize());
00688     assert(rY.GetSize()==rZ.size());
00689     for(unsigned i=0; i<rX.size(); i++)
00690     {
00691         rZ[i] = rX[i] + a*rY[i];
00692     }
00693 }
00694 
00695 
00696 template<CompressibilityType COMPRESSIBILITY_TYPE, unsigned DIM>
00697 double AbstractNonlinearElasticitySolver<COMPRESSIBILITY_TYPE,DIM>::TakeNewtonStep()
00698 {
00699     #ifdef MECH_VERBOSE
00700     Timer::Reset();
00701     #endif
00702 
00704     // Assemble Jacobian (and preconditioner)
00706     MechanicsEventHandler::BeginEvent(MechanicsEventHandler::ASSEMBLE);
00707     AssembleSystem(true, true);
00708     MechanicsEventHandler::EndEvent(MechanicsEventHandler::ASSEMBLE);
00709     #ifdef MECH_VERBOSE
00710     Timer::PrintAndReset("AssembleSystem");
00711     #endif
00712 
00713 
00715     //
00716     // Solve the linear system.
00717     // Three alternatives
00718     //   (a) GMRES with ILU preconditioner (or bjacobi=ILU on each process) [default]. Very poor on large problems.
00719     //   (b) GMRES with AMG preconditioner. Uncomment #define MECH_USE_HYPRE above. Requires Petsc3 with HYPRE installed.
00720     //   (c) Solve using MUMPS direct solver Uncomment #define MECH_USE_MUMPS above. Requires Petsc 3 with MUMPS installed.
00721     //
00723     MechanicsEventHandler::BeginEvent(MechanicsEventHandler::SOLVE);
00724 
00725     Vec solution;
00726     VecDuplicate(mpLinearSystem->rGetRhsVector(),&solution);
00727 
00728     Mat& r_jac = mpLinearSystem->rGetLhsMatrix();
00729 
00730     KSP solver;
00731     KSPCreate(PETSC_COMM_WORLD,&solver);
00732     PC pc;
00733     KSPGetPC(solver, &pc);
00734     
00735     #ifdef MECH_USE_MUMPS
00736         // no need for the precondition matrix
00737         KSPSetOperators(solver, r_jac, r_jac, DIFFERENT_NONZERO_PATTERN /*in precond between successive solves*/);
00738 
00739         KSPSetType(solver, KSPPREONLY);
00740         PCSetType(pc, PCLU);
00741         PCFactorSetMatSolverPackage(pc, MAT_SOLVER_MUMPS);
00742 
00743         // Allow matrices with zero diagonals to be solved
00744         PCFactorSetShiftNonzero(pc, PETSC_DECIDE);
00745         // "Do LU factorization in-place (saves memory)"
00746         PCASMSetUseInPlace(pc);
00747     #else
00748         Mat& r_precond_jac = mpPreconditionMatrixLinearSystem->rGetLhsMatrix();
00749         KSPSetOperators(solver, r_jac, r_precond_jac, DIFFERENT_NONZERO_PATTERN /*in precond between successive solves*/);
00750 
00751         unsigned num_restarts = 100;
00752         KSPSetType(solver,KSPGMRES);
00753         KSPGMRESSetRestart(solver,num_restarts); // gmres num restarts
00754     
00755         if(mKspAbsoluteTol < 0)
00756         {
00757             double ksp_rel_tol = 1e-6;
00758             KSPSetTolerances(solver, ksp_rel_tol, PETSC_DEFAULT, PETSC_DEFAULT, 10000 /*max iter*/); //hopefully with the preconditioner this max is way too high
00759         }
00760         else
00761         {
00762             KSPSetTolerances(solver, 1e-16, mKspAbsoluteTol, PETSC_DEFAULT, 10000 /*max iter*/); //hopefully with the preconditioner this max is way too high
00763         }
00764 
00765         #ifndef MECH_USE_HYPRE 
00766             PCSetType(pc, PCBJACOBI); // BJACOBI = ILU on each block (block = part of matrix on each process)
00767         #else
00769             // Speed up linear solve time massively for larger simulations (in fact GMRES may stagnate without
00770             // this for larger problems), by using a AMG preconditioner -- needs HYPRE installed 
00772             PetscOptionsSetValue("-pc_hypre_type", "boomeramg");
00773             // PetscOptionsSetValue("-pc_hypre_boomeramg_max_iter", "1");
00774             // PetscOptionsSetValue("-pc_hypre_boomeramg_strong_threshold", "0.0");
00775         
00776             PCSetType(pc, PCHYPRE);
00777             KSPSetPreconditionerSide(solver, PC_RIGHT);
00778             
00779             // other possible preconditioners..
00780             //PCBlockDiagonalMechanics* p_custom_pc = new PCBlockDiagonalMechanics(solver, r_precond_jac, mBlock1Size, mBlock2Size);
00781             //PCLDUFactorisationMechanics* p_custom_pc = new PCLDUFactorisationMechanics(solver, r_precond_jac, mBlock1Size, mBlock2Size);
00782             //remember to delete memory..
00783             //KSPSetPreconditionerSide(solver, PC_RIGHT);
00784         #endif
00785 
00786         // uncomment to get convergence information
00787         //#ifdef MECH_VERBOSE
00788         //PetscOptionsSetValue("-ksp_monitor","");
00789         //#endif
00790 
00791         KSPSetFromOptions(solver);
00792         KSPSetUp(solver);
00793     #endif     
00794 
00795     #ifdef MECH_VERBOSE
00796     Timer::PrintAndReset("KSP Setup");
00797     #endif
00798     
00799     KSPSolve(solver,mpLinearSystem->rGetRhsVector(),solution);
00800     
00801     int num_iters;
00802     KSPGetIterationNumber(solver, &num_iters);
00803     if(num_iters==0)
00804     {
00805         VecDestroy(solution);
00806         KSPDestroy(solver);
00807         EXCEPTION("KSP Absolute tolerance was too high, linear system wasn't solved - there will be no decrease in Newton residual. Decrease KspAbsoluteTolerance");
00808     }    
00809     
00810     #ifdef MECH_VERBOSE
00811     Timer::PrintAndReset("KSP Solve");
00812 
00813     std::cout << "[" << PetscTools::GetMyRank() << "]: Num iterations = " << num_iters << "\n" << std::flush;
00814     #endif
00815 
00816     MechanicsEventHandler::EndEvent(MechanicsEventHandler::SOLVE);
00818     // Update the solution
00819     //  Newton method:       sol = sol - update, where update=Jac^{-1}*residual
00820     //  Newton with damping: sol = sol - s*update, where s is chosen
00821     //   such that |residual(sol)| is minimised. Damping is important to
00822     //   avoid initial divergence.
00823     //
00824     // Normally, finding the best s from say 0.05,0.1,0.2,..,1.0 is cheap,
00825     // but this is not the case in cardiac electromechanics calculations.
00826     // Therefore, we initially check s=1 (expected to be the best most of the
00827     // time, then s=0.9. If the norm of the residual increases, we assume
00828     // s=1 is the best. Otherwise, check s=0.8 to see if s=0.9 is a local min.
00830     MechanicsEventHandler::BeginEvent(MechanicsEventHandler::UPDATE);
00831     double new_norm_resid = UpdateSolutionUsingLineSearch(solution);
00832     MechanicsEventHandler::EndEvent(MechanicsEventHandler::UPDATE);
00833 
00834     VecDestroy(solution);
00835     KSPDestroy(solver);
00836 
00837     return new_norm_resid;
00838 }
00839 
00840 template<CompressibilityType COMPRESSIBILITY_TYPE, unsigned DIM>
00841 void AbstractNonlinearElasticitySolver<COMPRESSIBILITY_TYPE,DIM>::PrintLineSearchResult(double s, double residNorm)
00842 {
00843     #ifdef MECH_VERBOSE
00844     std::cout << "\tTesting s = " << s << ", |f| = " << residNorm << "\n" << std::flush;
00845     #endif
00846 }
00847 
00848 template<CompressibilityType COMPRESSIBILITY_TYPE, unsigned DIM>
00849 double AbstractNonlinearElasticitySolver<COMPRESSIBILITY_TYPE,DIM>::UpdateSolutionUsingLineSearch(Vec solution)
00850 {
00851     double initial_norm_resid = CalculateResidualNorm();
00852     #ifdef MECH_VERBOSE
00853     std::cout << "\tInitial |f| [corresponding to s=0] is " << initial_norm_resid << "\n"  << std::flush;
00854     #endif
00855 
00856     ReplicatableVector update(solution);
00857 
00858     std::vector<double> old_solution = mCurrentSolution;
00859 
00860     std::vector<double> damping_values; // = {1.0, 0.9, .., 0.2, 0.1, 0.05} ie size 11
00861     for (unsigned i=10; i>=1; i--)
00862     {
00863         damping_values.push_back((double)i/10.0);
00864     }
00865     damping_values.push_back(0.05);
00866     assert(damping_values.size()==11);
00867 
00868 
00870     // let mCurrentSolution = old_solution - damping_val[0]*update; and compute residual
00871     unsigned index = 0;
00872     VectorSum(old_solution, update, -damping_values[index], mCurrentSolution);
00873     double current_resid_norm = ComputeResidualAndGetNorm(true);
00874     PrintLineSearchResult(damping_values[index], current_resid_norm);
00875 
00877     // let mCurrentSolution = old_solution - damping_val[1]*update; and compute residual
00878     index = 1;
00879     VectorSum(old_solution, update, -damping_values[index], mCurrentSolution);
00880     double next_resid_norm = ComputeResidualAndGetNorm(true);
00881     PrintLineSearchResult(damping_values[index], next_resid_norm);
00882 
00883     index = 2;
00884     // While f(s_next) < f(s_current), [f = residnorm], keep trying new damping values,
00885     // ie exit thus loop when next norm of the residual first increases
00886     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
00887             || ( (next_resid_norm < current_resid_norm) && (index<damping_values.size()) ) )
00888     {
00889         current_resid_norm = next_resid_norm;
00890 
00891         // let mCurrentSolution = old_solution - damping_val*update; and compute residual
00892         VectorSum(old_solution, update, -damping_values[index], mCurrentSolution);
00893         next_resid_norm = ComputeResidualAndGetNorm(true);
00894         PrintLineSearchResult(damping_values[index], next_resid_norm);
00895 
00896         index++;
00897     }
00898 
00899     unsigned best_index;
00900 
00901     if(index==damping_values.size() && (next_resid_norm < current_resid_norm))
00902     {
00903         // Difficult to come up with large forces/tractions such that it had to
00904         // test right down to s=0.05, but overall doesn't fail.
00905         // The possible damping values have been manually temporarily altered to
00906         // get this code to be called, it appears to work correctly. Even if it
00907         // didn't tests wouldn't fail, they would just be v. slightly less efficient.
00908         #define COVERAGE_IGNORE
00909         // if we exited because we got to the end of the possible damping values, the
00910         // best one was the last one (excl the final index++ at the end)
00911         current_resid_norm = next_resid_norm;
00912         best_index = index-1;
00913         #undef COVERAGE_IGNORE
00914     }
00915     else
00916     {
00917         // else the best one must have been the second last one (excl the final index++ at the end)
00918         // (as we would have exited when the resid norm first increased)
00919         best_index = index-2;
00920     }
00921 
00922     // check out best was better than the original residual-norm
00923     if (initial_norm_resid < current_resid_norm)
00924     {
00925         #define COVERAGE_IGNORE
00926         // Have to use an assert/exit here as the following exception causes a seg fault (in cardiac mech problems?)
00927         // Don't know why
00928         std::cout << "CHASTE ERROR: (AbstractNonlinearElasticitySolver.hpp): Residual does not appear to decrease in newton direction, quitting.\n" << std::flush;
00929         exit(0);
00930         //EXCEPTION("Residual does not appear to decrease in newton direction, quitting");
00931         #undef COVERAGE_IGNORE
00932     }
00933 
00934     #ifdef MECH_VERBOSE
00935     std::cout << "\tBest s = " << damping_values[best_index] << "\n"  << std::flush;
00936     #endif
00937     VectorSum(old_solution, update, -damping_values[best_index], mCurrentSolution);
00938 
00939     return current_resid_norm;
00940 
00942 //        double best_norm_resid = DBL_MAX;
00943 //        double best_damping_value = 0.0;
00944 //
00945 //        std::vector<double> damping_values;
00946 //        damping_values.reserve(12);
00947 //        damping_values.push_back(0.0);
00948 //        damping_values.push_back(0.05);
00949 //        for (unsigned i=1; i<=10; i++)
00950 //        {
00951 //            damping_values.push_back((double)i/10.0);
00952 //        }
00953 //
00954 //        for (unsigned i=0; i<damping_values.size(); i++)
00955 //        {
00956 //            for (unsigned j=0; j<mNumDofs; j++)
00957 //            {
00958 //                mCurrentSolution[j] = old_solution[j] - damping_values[i]*update[j];
00959 //            }
00960 //
00961 //            // compute residual
00962 //            double norm_resid = ComputeResidualAndGetNorm();
00963 //
00964 //            std::cout << "\tTesting s = " << damping_values[i] << ", |f| = " << norm_resid << "\n" << std::flush;
00965 //            if (norm_resid < best_norm_resid)
00966 //            {
00967 //                best_norm_resid = norm_resid;
00968 //                best_damping_value = damping_values[i];
00969 //            }
00970 //        }
00971 //
00972 //        if (best_damping_value == 0.0)
00973 //        {
00974 //            #define COVERAGE_IGNORE
00975 //            assert(0);
00976 //            EXCEPTION("Residual does not decrease in newton direction, quitting");
00977 //            #undef COVERAGE_IGNORE
00978 //        }
00979 //        else
00980 //        {
00981 //            std::cout << "\tBest s = " << best_damping_value << "\n"  << std::flush;
00982 //        }
00983 //        //Timer::PrintAndReset("Find best damping");
00984 //
00985 //        // implement best update and recalculate residual
00986 //        for (unsigned j=0; j<mNumDofs; j++)
00987 //        {
00988 //            mCurrentSolution[j] = old_solution[j] - best_damping_value*update[j];
00989 //        }
00990 }
00991 
00992 
00993 
00994 template<CompressibilityType COMPRESSIBILITY_TYPE, unsigned DIM>
00995 void AbstractNonlinearElasticitySolver<COMPRESSIBILITY_TYPE,DIM>::PostNewtonStep(unsigned counter, double normResidual)
00996 {
00997 }
00998 
00999 
01000 template<CompressibilityType COMPRESSIBILITY_TYPE, unsigned DIM>
01001 AbstractNonlinearElasticitySolver<COMPRESSIBILITY_TYPE,DIM>::AbstractNonlinearElasticitySolver(QuadraticMesh<DIM>* pQuadMesh,
01002                                                                                                c_vector<double,DIM> bodyForce,
01003                                                                                                double density,
01004                                                                                                std::string outputDirectory,
01005                                                                                                std::vector<unsigned>& fixedNodes)
01006     : mpQuadMesh(pQuadMesh),
01007       mpQuadratureRule(NULL),
01008       mpBoundaryQuadratureRule(NULL),
01009       mKspAbsoluteTol(-1),
01010       mBodyForce(bodyForce),
01011       mDensity(density),
01012       mFixedNodes(fixedNodes),
01013       mOutputDirectory(outputDirectory),
01014       mpOutputFileHandler(NULL),
01015       mWriteOutputEachNewtonIteration(false),
01016       mNumNewtonIterations(0),
01017       mUsingBodyForceFunction(false),
01018       mUsingTractionBoundaryConditionFunction(false),
01019       mCurrentTime(0.0)
01020 {
01021     assert(DIM==2 || DIM==3);
01022     assert(density > 0);
01023     assert(fixedNodes.size() > 0);
01024     assert(pQuadMesh != NULL);
01025 
01026     if(COMPRESSIBILITY_TYPE==COMPRESSIBLE)
01027     {
01028         mNumDofs = DIM*mpQuadMesh->GetNumNodes();
01029     }
01030     else
01031     {
01032         mNumDofs = DIM*mpQuadMesh->GetNumNodes() + mpQuadMesh->GetNumVertices();
01033     }
01034 
01035     mWriteOutput = (mOutputDirectory != "");
01036     if(mWriteOutput)
01037     {
01038         mpOutputFileHandler = new OutputFileHandler(mOutputDirectory);
01039     }
01040 }
01041 
01042 
01043 
01044 template<CompressibilityType COMPRESSIBILITY_TYPE, unsigned DIM>
01045 AbstractNonlinearElasticitySolver<COMPRESSIBILITY_TYPE,DIM>::~AbstractNonlinearElasticitySolver()
01046 {
01047     delete mpLinearSystem;
01048     delete mpPreconditionMatrixLinearSystem;
01049     delete mpQuadratureRule;
01050     delete mpBoundaryQuadratureRule;
01051     if(mpOutputFileHandler)
01052     {
01053         delete mpOutputFileHandler;
01054     }
01055 }
01056 
01057 
01058 template<CompressibilityType COMPRESSIBILITY_TYPE, unsigned DIM>
01059 void AbstractNonlinearElasticitySolver<COMPRESSIBILITY_TYPE,DIM>::Solve(double tol,
01060                                                                         unsigned maxNumNewtonIterations,
01061                                                                         bool quitIfNoConvergence)
01062 {
01063     WriteCurrentDeformation("initial");
01064 
01065 
01066     if(mWriteOutputEachNewtonIteration)
01067     {
01068         WriteCurrentDeformation("newton_iteration", 0);
01069     }
01070 
01071     // compute residual
01072     double norm_resid = ComputeResidualAndGetNorm(false);
01073     #ifdef MECH_VERBOSE
01074     std::cout << "\nNorm of residual is " << norm_resid << "\n";
01075     #endif
01076 
01077     mNumNewtonIterations = 0;
01078     unsigned iteration_number = 1;
01079 
01080     if (tol < 0) // ie if wasn't passed in as a parameter
01081     {
01082         tol = NEWTON_REL_TOL*norm_resid;
01083 
01084         #define COVERAGE_IGNORE // not going to have tests in cts for everything
01085         if (tol > MAX_NEWTON_ABS_TOL)
01086         {
01087             tol = MAX_NEWTON_ABS_TOL;
01088         }
01089         if (tol < MIN_NEWTON_ABS_TOL)
01090         {
01091             tol = MIN_NEWTON_ABS_TOL;
01092         }
01093         #undef COVERAGE_IGNORE
01094     }
01095 
01096     #ifdef MECH_VERBOSE
01097     std::cout << "Solving with tolerance " << tol << "\n";
01098     #endif
01099 
01100     while (norm_resid > tol && iteration_number<=maxNumNewtonIterations)
01101     {
01102         #ifdef MECH_VERBOSE
01103         std::cout <<  "\n-------------------\n"
01104                   <<   "Newton iteration " << counter
01105                   << ":\n-------------------\n";
01106         #endif
01107 
01108         // take newton step (and get returned residual)
01109         norm_resid = TakeNewtonStep();
01110 
01111         #ifdef MECH_VERBOSE
01112         std::cout << "Norm of residual is " << norm_resid << "\n";
01113         #endif
01114         if (mWriteOutputEachNewtonIteration)
01115         {
01116             WriteCurrentDeformation("newton_iteration", iteration_number);
01117         }
01118 
01119         mNumNewtonIterations = iteration_number;
01120 
01121         PostNewtonStep(iteration_number,norm_resid);
01122 
01123         iteration_number++;
01124         if (iteration_number==20)
01125         {
01126             #define COVERAGE_IGNORE
01127             EXCEPTION("Not converged after 20 newton iterations, quitting");
01128             #undef COVERAGE_IGNORE
01129         }
01130     }
01131 
01132     if ((norm_resid > tol) && quitIfNoConvergence)
01133     {
01134         #define COVERAGE_IGNORE
01135         EXCEPTION("Failed to converge");
01136         #undef COVERAGE_IGNORE
01137     }
01138 
01139     // write the final solution
01140     WriteCurrentDeformation("solution");
01141 }
01142 
01143 
01144 template<CompressibilityType COMPRESSIBILITY_TYPE, unsigned DIM>
01145 void AbstractNonlinearElasticitySolver<COMPRESSIBILITY_TYPE,DIM>::WriteCurrentDeformation(std::string fileName, int counterToAppend)
01146 {
01147     // only write output if the flag mWriteOutput has been set
01148     if (!mWriteOutput)
01149     {
01150         return;
01151     }
01152 
01153     std::stringstream file_name;
01154     file_name << fileName;
01155     if(counterToAppend >= 0)
01156     {
01157         file_name << "_" << counterToAppend;
01158     }
01159     file_name << ".nodes";
01160 
01161     out_stream p_file = mpOutputFileHandler->OpenOutputFile(file_name.str());
01162 
01163     std::vector<c_vector<double,DIM> >& r_deformed_position = rGetDeformedPosition();
01164     for (unsigned i=0; i<r_deformed_position.size(); i++)
01165     {
01166         for (unsigned j=0; j<DIM; j++)
01167         {
01168            * p_file << r_deformed_position[i](j) << " ";
01169         }
01170        * p_file << "\n";
01171     }
01172     p_file->close();
01173 }
01174 
01175 
01176 template<CompressibilityType COMPRESSIBILITY_TYPE, unsigned DIM>
01177 unsigned AbstractNonlinearElasticitySolver<COMPRESSIBILITY_TYPE,DIM>::GetNumNewtonIterations()
01178 {
01179     return mNumNewtonIterations;
01180 }
01181 
01182 
01183 template<CompressibilityType COMPRESSIBILITY_TYPE, unsigned DIM>
01184 void AbstractNonlinearElasticitySolver<COMPRESSIBILITY_TYPE,DIM>::SetFunctionalBodyForce(c_vector<double,DIM> (*pFunction)(c_vector<double,DIM>& X, double t))
01185 {
01186     mUsingBodyForceFunction = true;
01187     mpBodyForceFunction = pFunction;
01188 }
01189 
01190 
01191 template<CompressibilityType COMPRESSIBILITY_TYPE, unsigned DIM>
01192 void AbstractNonlinearElasticitySolver<COMPRESSIBILITY_TYPE,DIM>::SetWriteOutput(bool writeOutput)
01193 {
01194     if (writeOutput && (mOutputDirectory==""))
01195     {
01196         EXCEPTION("Can't write output if no output directory was given in constructor");
01197     }
01198     mWriteOutput = writeOutput;
01199 }
01200 
01201 template<CompressibilityType COMPRESSIBILITY_TYPE, unsigned DIM>
01202 void AbstractNonlinearElasticitySolver<COMPRESSIBILITY_TYPE,DIM>::SetSurfaceTractionBoundaryConditions(
01203             std::vector<BoundaryElement<DIM-1,DIM>*>& rBoundaryElements,
01204             std::vector<c_vector<double,DIM> >& rSurfaceTractions)
01205 {
01206     assert(rBoundaryElements.size()==rSurfaceTractions.size());
01207     mBoundaryElements = rBoundaryElements;
01208     mSurfaceTractions = rSurfaceTractions;
01209 }
01210 
01211 template<CompressibilityType COMPRESSIBILITY_TYPE, unsigned DIM>
01212 void AbstractNonlinearElasticitySolver<COMPRESSIBILITY_TYPE,DIM>::SetFunctionalTractionBoundaryCondition(
01213             std::vector<BoundaryElement<DIM-1,DIM>*> rBoundaryElements,
01214             c_vector<double,DIM> (*pFunction)(c_vector<double,DIM>& X, double t))
01215 {
01216     mBoundaryElements = rBoundaryElements;
01217     mUsingTractionBoundaryConditionFunction = true;
01218     mpTractionBoundaryConditionFunction = pFunction;
01219 }
01220 
01221 
01222 template<CompressibilityType COMPRESSIBILITY_TYPE, unsigned DIM>
01223 std::vector<c_vector<double,DIM> >& AbstractNonlinearElasticitySolver<COMPRESSIBILITY_TYPE,DIM>::rGetDeformedPosition()
01224 {
01225     mDeformedPosition.resize(mpQuadMesh->GetNumNodes(), zero_vector<double>(DIM));
01226     for (unsigned i=0; i<mpQuadMesh->GetNumNodes(); i++)
01227     {
01228         for (unsigned j=0; j<DIM; j++)
01229         {
01230             mDeformedPosition[i](j) = mpQuadMesh->GetNode(i)->rGetLocation()[j] + mCurrentSolution[DIM*i+j];
01231         }
01232     }
01233     return mDeformedPosition;
01234 }
01235 
01236 
01237 template<CompressibilityType COMPRESSIBILITY_TYPE, unsigned DIM>
01238 void AbstractNonlinearElasticitySolver<COMPRESSIBILITY_TYPE,DIM>::CreateCmguiOutput()
01239 {
01240     if(mOutputDirectory=="")
01241     {
01242         EXCEPTION("No output directory was given so no output was written, cannot convert to cmgui format");
01243     }
01244 
01245     CmguiDeformedSolutionsWriter<DIM> writer(mOutputDirectory + "/cmgui",
01246                                              "solution",
01247                                              *mpQuadMesh,
01248                                              WRITE_QUADRATIC_MESH);
01249 
01250     std::vector<c_vector<double,DIM> >& r_deformed_positions = rGetDeformedPosition();
01251     writer.WriteInitialMesh(); // this writes solution_0.exnode and .exelem
01252     writer.WriteDeformationPositions(r_deformed_positions, 1); // this writes the final solution as solution_1.exnode
01253     writer.WriteCmguiScript(); // writes LoadSolutions.com
01254 }
01255 
01256 //
01257 // Constant setting definitions
01258 //
01259 template<CompressibilityType COMPRESSIBILITY_TYPE, unsigned DIM>
01260 double AbstractNonlinearElasticitySolver<COMPRESSIBILITY_TYPE,DIM>::MAX_NEWTON_ABS_TOL = 1e-7;
01261 
01262 template<CompressibilityType COMPRESSIBILITY_TYPE, unsigned DIM>
01263 double AbstractNonlinearElasticitySolver<COMPRESSIBILITY_TYPE,DIM>::MIN_NEWTON_ABS_TOL = 1e-10;
01264 
01265 template<CompressibilityType COMPRESSIBILITY_TYPE, unsigned DIM>
01266 double AbstractNonlinearElasticitySolver<COMPRESSIBILITY_TYPE,DIM>::NEWTON_REL_TOL = 1e-4;
01267 
01268 
01269 #endif /*ABSTRACTNONLINEARELASTICITYSOLVER_HPP_*/

Generated on Mon Apr 18 11:35:36 2011 for Chaste by  doxygen 1.5.5