00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
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 "AbstractPerElementWriter.hpp"
00053 #include "petscsnes.h"
00054
00055
00056
00057
00058
00059
00064 typedef enum StrainType_
00065 {
00066 DEFORMATION_GRADIENT_F = 0,
00067 DEFORMATION_TENSOR_C,
00068 LAGRANGE_STRAIN_E
00069 } StrainType;
00070
00071
00072
00073
00074 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00075 extern PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);
00076 #endif
00077
00078
00080
00082
00092 template<unsigned DIM>
00093 PetscErrorCode AbstractNonlinearElasticitySolver_ComputeResidual(SNES snes,
00094 Vec currentGuess,
00095 Vec residualVector,
00096 void* pContext);
00097
00098 #if ( PETSC_VERSION_MAJOR==3 && PETSC_VERSION_MINOR>=5 )
00099
00111 template<unsigned DIM>
00112 PetscErrorCode AbstractNonlinearElasticitySolver_ComputeJacobian(SNES snes,
00113 Vec currentGuess,
00114 Mat globalJacobian,
00115 Mat preconditioner,
00116 void* pContext);
00117 #else
00118
00129 template<unsigned DIM>
00130 PetscErrorCode AbstractNonlinearElasticitySolver_ComputeJacobian(SNES snes,
00131 Vec currentGuess,
00132 Mat* pGlobalJacobian,
00133 Mat* pPreconditioner,
00134 MatStructure* pMatStructure,
00135 void* pContext);
00136 #endif
00137
00138 template <unsigned DIM>
00139 class AbstractNonlinearElasticitySolver;
00140
00147 template<unsigned DIM>
00148 class StressPerElementWriter : public AbstractPerElementWriter<DIM, DIM, DIM*DIM>
00149 {
00150 private:
00151 AbstractNonlinearElasticitySolver<DIM>* mpSolver;
00152 public:
00153
00158 StressPerElementWriter(AbstractTetrahedralMesh<DIM, DIM>* pMesh,
00159 AbstractNonlinearElasticitySolver<DIM>* pSolver)
00160 : AbstractPerElementWriter<DIM, DIM, DIM*DIM>(pMesh),
00161 mpSolver(pSolver)
00162 {
00163 }
00164
00172 void Visit(Element<DIM, DIM>* pElement, unsigned localIndex, c_vector<double, DIM*DIM>& rData)
00173 {
00177 assert(localIndex == pElement->GetIndex());
00178
00179 c_matrix<double, DIM, DIM> data = mpSolver->GetAverageStressPerElement(localIndex);
00180 for (unsigned i=0; i<DIM; i++)
00181 {
00182 for (unsigned j=0; j<DIM; j++)
00183 {
00184 rData[i*DIM+j] = data(i,j);
00185 }
00186 }
00187 }
00188 };
00189
00201 template<unsigned DIM>
00202 class AbstractNonlinearElasticitySolver : public AbstractContinuumMechanicsSolver<DIM>
00203 {
00204 friend class StressRecoveror<DIM>;
00205 friend class VtkNonlinearElasticitySolutionWriter<DIM>;
00206
00207
00208 protected:
00209
00211 static const size_t NUM_VERTICES_PER_ELEMENT = DIM+1;
00212
00214 static const size_t NUM_NODES_PER_ELEMENT = (DIM+1)*(DIM+2)/2;
00215
00217 static const size_t NUM_NODES_PER_BOUNDARY_ELEMENT = DIM*(DIM+1)/2;
00218
00223 static const size_t BOUNDARY_STENCIL_SIZE = DIM*NUM_NODES_PER_BOUNDARY_ELEMENT;
00224
00232 static double MAX_NEWTON_ABS_TOL;
00233
00235 static double MIN_NEWTON_ABS_TOL;
00236
00238 static double NEWTON_REL_TOL;
00239
00244 SolidMechanicsProblemDefinition<DIM>& mrProblemDefinition;
00245
00250 Mat& mrJacobianMatrix;
00251
00258 c_matrix<double,DIM,DIM> mChangeOfBasisMatrix;
00259
00265 double mKspAbsoluteTol;
00266
00272 bool mWriteOutputEachNewtonIteration;
00273
00277 FourthOrderTensor<DIM,DIM,DIM,DIM> dTdE;
00278
00280 unsigned mNumNewtonIterations;
00281
00287 double mCurrentTime;
00288
00289
00293 bool mCheckedOutwardNormals;
00294
00299 bool mUseSnesSolver;
00300
00306 double mLastDampingValue;
00307
00311 bool mFirstStep;
00312
00317 bool mTakeFullFirstNewtonStep;
00318
00322 bool mPetscDirectSolve;
00323
00331 bool mIncludeActiveTension;
00332
00337 bool mSetComputeAverageStressPerElement;
00338
00346 std::vector<c_vector<double,DIM*(DIM+1)/2> > mAverageStressesPerElement;
00347
00355 void AddStressToAverageStressPerElement(c_matrix<double,DIM,DIM>& rT, unsigned elementIndex);
00356
00365 virtual void SetKspSolverAndPcType(KSP solver);
00366
00367
00380 virtual void AssembleSystem(bool assembleResidual, bool assembleLinearSystem)=0;
00381
00382
00383
00391 virtual void FinishAssembleSystem(bool assembleResidual, bool assembleLinearSystem);
00392
00393
00395
00396
00397
00399
00407 void GetElementCentroidStrain(StrainType strainType,
00408 Element<DIM,DIM>& rElement,
00409 c_matrix<double,DIM,DIM>& rDeformationGradient);
00410
00429 virtual void AddActiveStressAndStressDerivative(c_matrix<double,DIM,DIM>& rC,
00430 unsigned elementIndex,
00431 unsigned currentQuadPointGlobalIndex,
00432 c_matrix<double,DIM,DIM>& rT,
00433 FourthOrderTensor<DIM,DIM,DIM,DIM>& rDTdE,
00434 bool addToDTdE)
00435 {
00436
00437 }
00438
00446 virtual void SetupChangeOfBasisMatrix(unsigned elementIndex, unsigned currentQuadPointGlobalIndex)
00447 {
00448 }
00449
00450
00451
00470 void AssembleOnBoundaryElement(BoundaryElement<DIM-1, DIM>& rBoundaryElement,
00471 c_matrix<double, BOUNDARY_STENCIL_SIZE, BOUNDARY_STENCIL_SIZE>& rAelem,
00472 c_vector<double, BOUNDARY_STENCIL_SIZE>& rBelem,
00473 bool assembleResidual,
00474 bool assembleJacobian,
00475 unsigned boundaryConditionIndex);
00476
00477
00492 bool ShouldAssembleMatrixTermForPressureOnDeformedBc();
00493
00513 void AssembleOnBoundaryElementForPressureOnDeformedBc(BoundaryElement<DIM-1,DIM>& rBoundaryElement,
00514 c_matrix<double,BOUNDARY_STENCIL_SIZE,BOUNDARY_STENCIL_SIZE>& rAelem,
00515 c_vector<double,BOUNDARY_STENCIL_SIZE>& rBelem,
00516 bool assembleResidual,
00517 bool assembleJacobian,
00518 unsigned boundaryConditionIndex);
00519
00521
00522
00523
00525
00537 double ComputeResidualAndGetNorm(bool allowException);
00538
00540 double CalculateResidualNorm();
00541
00551 void VectorSum(std::vector<double>& rX, ReplicatableVector& rY, double a, std::vector<double>& rZ);
00552
00560 void PrintLineSearchResult(double s, double residNorm);
00561
00569 double TakeNewtonStep();
00570
00580 double UpdateSolutionUsingLineSearch(Vec solution);
00581
00589 virtual void PostNewtonStep(unsigned counter, double normResidual);
00590
00596 void SolveNonSnes(double tol=-1.0);
00597
00598
00600
00601
00602
00604 public:
00611 void ComputeResidual(Vec currentGuess, Vec residualVector);
00612
00620 void ComputeJacobian(Vec currentGuess, Mat* pJacobian, Mat* pPreconditioner);
00621
00622 private:
00627 void SolveSnes();
00628
00629 public:
00630
00641 AbstractNonlinearElasticitySolver(AbstractTetrahedralMesh<DIM,DIM>& rQuadMesh,
00642 SolidMechanicsProblemDefinition<DIM>& rProblemDefinition,
00643 std::string outputDirectory,
00644 CompressibilityType compressibilityType);
00645
00649 virtual ~AbstractNonlinearElasticitySolver();
00650
00656 void Solve(double tol=-1.0);
00657
00667 void SetIncludeActiveTension(bool includeActiveTension = true)
00668 {
00669 mIncludeActiveTension = includeActiveTension;
00670 }
00671
00675 unsigned GetNumNewtonIterations();
00676
00677
00685 void SetWriteOutputEachNewtonIteration(bool writeOutputEachNewtonIteration=true)
00686 {
00687 mWriteOutputEachNewtonIteration = writeOutputEachNewtonIteration;
00688 }
00689
00696 void SetKspAbsoluteTolerance(double kspAbsoluteTolerance)
00697 {
00698 assert(kspAbsoluteTolerance > 0);
00699 mKspAbsoluteTol = kspAbsoluteTolerance;
00700 }
00701
00715 void SetTakeFullFirstNewtonStep(bool takeFullFirstStep = true)
00716 {
00717 mTakeFullFirstNewtonStep = takeFullFirstStep;
00718 }
00719
00733 void SetUsePetscDirectSolve(bool usePetscDirectSolve = true)
00734 {
00735 mPetscDirectSolve = usePetscDirectSolve;
00736 }
00737
00738
00745 void SetCurrentTime(double time)
00746 {
00747 mCurrentTime = time;
00748 }
00749
00754 void CreateCmguiOutput();
00755
00756
00770 void WriteCurrentStrains(StrainType strainType, std::string fileName, int counterToAppend = -1);
00771
00778 void SetComputeAverageStressPerElementDuringSolve(bool setComputeAverageStressPerElement = true);
00779
00793 void WriteCurrentAverageElementStresses(std::string fileName, int counterToAppend = -1);
00794
00799 std::vector<c_vector<double,DIM> >& rGetSpatialSolution();
00800
00805 std::vector<c_vector<double,DIM> >& rGetDeformedPosition();
00806
00815 c_matrix<double,DIM,DIM> GetAverageStressPerElement(unsigned elementIndex);
00816 };
00817
00819
00821
00822 template<unsigned DIM>
00823 AbstractNonlinearElasticitySolver<DIM>::AbstractNonlinearElasticitySolver(AbstractTetrahedralMesh<DIM,DIM>& rQuadMesh,
00824 SolidMechanicsProblemDefinition<DIM>& rProblemDefinition,
00825 std::string outputDirectory,
00826 CompressibilityType compressibilityType)
00827 : AbstractContinuumMechanicsSolver<DIM>(rQuadMesh, rProblemDefinition, outputDirectory, compressibilityType),
00828 mrProblemDefinition(rProblemDefinition),
00829 mrJacobianMatrix(this->mSystemLhsMatrix),
00830 mKspAbsoluteTol(-1),
00831 mWriteOutputEachNewtonIteration(false),
00832 mNumNewtonIterations(0),
00833 mCurrentTime(0.0),
00834 mCheckedOutwardNormals(false),
00835 mLastDampingValue(0.0),
00836 mIncludeActiveTension(true),
00837 mSetComputeAverageStressPerElement(false)
00838 {
00839 mUseSnesSolver = (mrProblemDefinition.GetSolveUsingSnes() ||
00840 CommandLineArguments::Instance()->OptionExists("-mech_use_snes") );
00841
00842 mChangeOfBasisMatrix = identity_matrix<double>(DIM,DIM);
00843
00844 mTakeFullFirstNewtonStep = CommandLineArguments::Instance()->OptionExists("-mech_full_first_newton_step");
00845 mPetscDirectSolve = CommandLineArguments::Instance()->OptionExists("-mech_petsc_direct_solve");
00846 }
00847
00848 template<unsigned DIM>
00849 AbstractNonlinearElasticitySolver<DIM>::~AbstractNonlinearElasticitySolver()
00850 {
00851 }
00852
00853
00854 template<unsigned DIM>
00855 void AbstractNonlinearElasticitySolver<DIM>::FinishAssembleSystem(bool assembleResidual, bool assembleJacobian)
00856 {
00857 PetscVecTools::Finalise(this->mResidualVector);
00858
00859 if (assembleJacobian)
00860 {
00861 PetscMatTools::SwitchWriteMode(mrJacobianMatrix);
00862 PetscMatTools::SwitchWriteMode(this->mPreconditionMatrix);
00863
00864 VecCopy(this->mResidualVector, this->mLinearSystemRhsVector);
00865 }
00866
00867
00868 if(assembleJacobian)
00869 {
00870 this->ApplyDirichletBoundaryConditions(NONLINEAR_PROBLEM_APPLY_TO_EVERYTHING, this->mCompressibilityType==COMPRESSIBLE);
00871 }
00872 else if (assembleResidual)
00873 {
00874 this->ApplyDirichletBoundaryConditions(NONLINEAR_PROBLEM_APPLY_TO_RESIDUAL_ONLY, this->mCompressibilityType==COMPRESSIBLE);
00875 }
00876
00877 if (assembleResidual)
00878 {
00879 PetscVecTools::Finalise(this->mResidualVector);
00880 }
00881 if (assembleJacobian)
00882 {
00883 PetscMatTools::Finalise(mrJacobianMatrix);
00884 PetscMatTools::Finalise(this->mPreconditionMatrix);
00885 PetscVecTools::Finalise(this->mLinearSystemRhsVector);
00886 }
00887 }
00888
00889
00890 template<unsigned DIM>
00891 std::vector<c_vector<double,DIM> >& AbstractNonlinearElasticitySolver<DIM>::rGetSpatialSolution()
00892 {
00893 this->mSpatialSolution.clear();
00894 this->mSpatialSolution.resize(this->mrQuadMesh.GetNumNodes(), zero_vector<double>(DIM));
00895 for (unsigned i=0; i<this->mrQuadMesh.GetNumNodes(); i++)
00896 {
00897 for (unsigned j=0; j<DIM; j++)
00898 {
00899 this->mSpatialSolution[i](j) = this->mrQuadMesh.GetNode(i)->rGetLocation()[j] + this->mCurrentSolution[this->mProblemDimension*i+j];
00900 }
00901 }
00902 return this->mSpatialSolution;
00903 }
00904
00905 template<unsigned DIM>
00906 std::vector<c_vector<double,DIM> >& AbstractNonlinearElasticitySolver<DIM>::rGetDeformedPosition()
00907 {
00908 return rGetSpatialSolution();
00909 }
00910
00911
00912 template<unsigned DIM>
00913 void AbstractNonlinearElasticitySolver<DIM>::WriteCurrentStrains(StrainType strainType, std::string fileName, int counterToAppend)
00914 {
00915 if (!this->mWriteOutput)
00916 {
00917 return;
00918 }
00919
00920 std::stringstream file_name;
00921 file_name << fileName;
00922 if (counterToAppend >= 0)
00923 {
00924 file_name << "_" << counterToAppend;
00925 }
00926 file_name << ".strain";
00927
00928 out_stream p_file = this->mpOutputFileHandler->OpenOutputFile(file_name.str());
00929
00930 c_matrix<double,DIM,DIM> strain;
00931
00932 for (typename AbstractTetrahedralMesh<DIM,DIM>::ElementIterator iter = this->mrQuadMesh.GetElementIteratorBegin();
00933 iter != this->mrQuadMesh.GetElementIteratorEnd();
00934 ++iter)
00935 {
00936 GetElementCentroidStrain(strainType, *iter, strain);
00937 for(unsigned i=0; i<DIM; i++)
00938 {
00939 for(unsigned j=0; j<DIM; j++)
00940 {
00941 *p_file << strain(i,j) << " ";
00942 }
00943 }
00944 *p_file << "\n";
00945 }
00946 p_file->close();
00947 }
00948
00949
00950
00951 template<unsigned DIM>
00952 void AbstractNonlinearElasticitySolver<DIM>::WriteCurrentAverageElementStresses(std::string fileName, int counterToAppend)
00953 {
00954 if (!this->mWriteOutput)
00955 {
00956 return;
00957 }
00958
00959 if(!mSetComputeAverageStressPerElement)
00960 {
00961 EXCEPTION("Call SetComputeAverageStressPerElementDuringSolve() before solve if calling WriteCurrentAverageElementStresses()");
00962 }
00963
00964 std::stringstream file_name;
00965 file_name << fileName;
00966 if (counterToAppend >= 0)
00967 {
00968 file_name << "_" << counterToAppend;
00969 }
00970 file_name << ".stress";
00971 assert(mAverageStressesPerElement.size()==this->mrQuadMesh.GetNumElements());
00972
00973 StressPerElementWriter<DIM> stress_writer(&(this->mrQuadMesh),this);
00974 stress_writer.WriteData(*(this->mpOutputFileHandler), file_name.str());
00975 }
00976
00977
00978 template<unsigned DIM>
00979 void AbstractNonlinearElasticitySolver<DIM>::CreateCmguiOutput()
00980 {
00981 if (this->mOutputDirectory=="")
00982 {
00983 EXCEPTION("No output directory was given so no output was written, cannot convert to cmgui format");
00984 }
00985
00986 CmguiDeformedSolutionsWriter<DIM> writer(this->mOutputDirectory + "/cmgui",
00987 "solution",
00988 this->mrQuadMesh,
00989 WRITE_QUADRATIC_MESH);
00990
00991 std::vector<c_vector<double,DIM> >& r_deformed_positions = this->rGetDeformedPosition();
00992 writer.WriteInitialMesh();
00993 writer.WriteDeformationPositions(r_deformed_positions, 1);
00994 writer.WriteCmguiScript();
00995 }
00996
00997 template<unsigned DIM>
00998 void AbstractNonlinearElasticitySolver<DIM>::SetComputeAverageStressPerElementDuringSolve(bool setComputeAverageStressPerElement)
00999 {
01000 mSetComputeAverageStressPerElement = setComputeAverageStressPerElement;
01001 if(setComputeAverageStressPerElement && mAverageStressesPerElement.size()==0)
01002 {
01003 mAverageStressesPerElement.resize(this->mrQuadMesh.GetNumElements(), zero_vector<double>(DIM*(DIM+1)/2));
01004 }
01005 }
01006
01007 template<unsigned DIM>
01008 void AbstractNonlinearElasticitySolver<DIM>::AddStressToAverageStressPerElement(c_matrix<double,DIM,DIM>& rT, unsigned elemIndex)
01009 {
01010 assert(mSetComputeAverageStressPerElement);
01011 assert(elemIndex<this->mrQuadMesh.GetNumElements());
01012
01013
01014
01015
01016
01017
01018
01019
01020
01021 for(unsigned i=0; i<DIM*(DIM+1)/2; i++)
01022 {
01023 unsigned row;
01024 unsigned col;
01025 if(DIM==2)
01026 {
01027 row = i<=1 ? 0 : 1;
01028 col = i==0 ? 0 : 1;
01029 }
01030 else
01031 {
01032 row = i<=2 ? 0 : (i<=4? 1 : 2);
01033 col = i==0 ? 0 : (i==1 || i==3? 1 : 2);
01034 }
01035
01036 this->mAverageStressesPerElement[elemIndex](i) += rT(row,col);
01037 }
01038 }
01039
01040
01041 template<unsigned DIM>
01042 c_matrix<double,DIM,DIM> AbstractNonlinearElasticitySolver<DIM>::GetAverageStressPerElement(unsigned elementIndex)
01043 {
01044 if(!mSetComputeAverageStressPerElement)
01045 {
01046 EXCEPTION("Call SetComputeAverageStressPerElementDuringSolve() before solve if calling GetAverageStressesPerElement()");
01047 }
01048 assert(elementIndex<this->mrQuadMesh.GetNumElements());
01049
01050 c_matrix<double,DIM,DIM> stress;
01051
01052
01053
01054
01055
01056
01057
01058
01059
01060 if(DIM==2)
01061 {
01062 stress(0,0) = mAverageStressesPerElement[elementIndex](0);
01063 stress(1,0) = stress(0,1) = mAverageStressesPerElement[elementIndex](1);
01064 stress(1,1) = mAverageStressesPerElement[elementIndex](2);
01065 }
01066 else
01067 {
01068 stress(0,0) = mAverageStressesPerElement[elementIndex](0);
01069 stress(1,0) = stress(0,1) = mAverageStressesPerElement[elementIndex](1);
01070 stress(2,0) = stress(0,2) = mAverageStressesPerElement[elementIndex](2);
01071 stress(1,1) = mAverageStressesPerElement[elementIndex](3);
01072 stress(2,1) = stress(1,2) = mAverageStressesPerElement[elementIndex](4);
01073 stress(2,2) = mAverageStressesPerElement[elementIndex](5);
01074 }
01075
01076 return stress;
01077 }
01078
01079
01081
01083
01084 template<unsigned DIM>
01085 void AbstractNonlinearElasticitySolver<DIM>::GetElementCentroidStrain(StrainType strainType,
01086 Element<DIM,DIM>& rElement,
01087 c_matrix<double,DIM,DIM>& rStrain)
01088 {
01089 static c_matrix<double,DIM,DIM> jacobian;
01090 static c_matrix<double,DIM,DIM> inverse_jacobian;
01091 double jacobian_determinant;
01092
01093 this->mrQuadMesh.GetInverseJacobianForElement(rElement.GetIndex(), jacobian, jacobian_determinant, inverse_jacobian);
01094
01095
01096 static c_matrix<double,DIM,NUM_NODES_PER_ELEMENT> element_current_displacements;
01097 for (unsigned II=0; II<NUM_NODES_PER_ELEMENT; II++)
01098 {
01099 for (unsigned JJ=0; JJ<DIM; JJ++)
01100 {
01101 element_current_displacements(JJ,II) = this->mCurrentSolution[this->mProblemDimension*rElement.GetNodeGlobalIndex(II) + JJ];
01102 }
01103 }
01104
01105
01106 static c_matrix<double, DIM, NUM_NODES_PER_ELEMENT> grad_quad_phi;
01107 static c_matrix<double,DIM,DIM> grad_u;
01108
01109
01110
01111 ChastePoint<DIM> quadrature_point;
01112 if(DIM==2)
01113 {
01114 quadrature_point.rGetLocation()(0) = 1.0/3.0;
01115 quadrature_point.rGetLocation()(1) = 1.0/3.0;
01116 }
01117 else
01118 {
01119 assert(DIM==3);
01120 quadrature_point.rGetLocation()(0) = 1.0/4.0;
01121 quadrature_point.rGetLocation()(1) = 1.0/4.0;
01122 quadrature_point.rGetLocation()(2) = 1.0/4.0;
01123 }
01124
01125 QuadraticBasisFunction<DIM>::ComputeTransformedBasisFunctionDerivatives(quadrature_point, inverse_jacobian, grad_quad_phi);
01126
01127
01128 grad_u = zero_matrix<double>(DIM,DIM);
01129 for (unsigned node_index=0; node_index<NUM_NODES_PER_ELEMENT; node_index++)
01130 {
01131 for (unsigned i=0; i<DIM; i++)
01132 {
01133 for (unsigned M=0; M<DIM; M++)
01134 {
01135 grad_u(i,M) += grad_quad_phi(M,node_index)*element_current_displacements(i,node_index);
01136 }
01137 }
01138 }
01139
01140 c_matrix<double,DIM,DIM> deformation_gradient;
01141
01142 for (unsigned i=0; i<DIM; i++)
01143 {
01144 for (unsigned M=0; M<DIM; M++)
01145 {
01146 deformation_gradient(i,M) = (i==M?1:0) + grad_u(i,M);
01147 }
01148 }
01149
01150 switch(strainType)
01151 {
01152 case DEFORMATION_GRADIENT_F:
01153 {
01154 rStrain = deformation_gradient;
01155 break;
01156 }
01157 case DEFORMATION_TENSOR_C:
01158 {
01159 rStrain = prod(trans(deformation_gradient),deformation_gradient);
01160 break;
01161 }
01162 case LAGRANGE_STRAIN_E:
01163 {
01164 c_matrix<double,DIM,DIM> C = prod(trans(deformation_gradient),deformation_gradient);
01165 for (unsigned M=0; M<DIM; M++)
01166 {
01167 for (unsigned N=0; N<DIM; N++)
01168 {
01169 rStrain(M,N) = 0.5* ( C(M,N)-(M==N?1:0) );
01170 }
01171 }
01172 break;
01173 }
01174 default:
01175 {
01176 NEVER_REACHED;
01177 break;
01178 }
01179 }
01180 }
01181
01182
01183
01184 template<unsigned DIM>
01185 void AbstractNonlinearElasticitySolver<DIM>::AssembleOnBoundaryElement(
01186 BoundaryElement<DIM-1,DIM>& rBoundaryElement,
01187 c_matrix<double,BOUNDARY_STENCIL_SIZE,BOUNDARY_STENCIL_SIZE>& rAelem,
01188 c_vector<double,BOUNDARY_STENCIL_SIZE>& rBelem,
01189 bool assembleResidual,
01190 bool assembleJacobian,
01191 unsigned boundaryConditionIndex)
01192 {
01193 if( this->mrProblemDefinition.GetTractionBoundaryConditionType() == PRESSURE_ON_DEFORMED
01194 || this->mrProblemDefinition.GetTractionBoundaryConditionType() == FUNCTIONAL_PRESSURE_ON_DEFORMED)
01195 {
01196 AssembleOnBoundaryElementForPressureOnDeformedBc(rBoundaryElement, rAelem, rBelem,
01197 assembleResidual, assembleJacobian, boundaryConditionIndex);
01198 return;
01199 }
01200
01201 rAelem.clear();
01202 rBelem.clear();
01203
01204 if (assembleJacobian && !assembleResidual)
01205 {
01206
01207 return;
01208 }
01209
01210 c_vector<double, DIM> weighted_direction;
01211 double jacobian_determinant;
01212 this->mrQuadMesh.GetWeightedDirectionForBoundaryElement(rBoundaryElement.GetIndex(), weighted_direction, jacobian_determinant);
01213
01214 c_vector<double,NUM_NODES_PER_BOUNDARY_ELEMENT> phi;
01215
01216 for (unsigned quad_index=0; quad_index<this->mpBoundaryQuadratureRule->GetNumQuadPoints(); quad_index++)
01217 {
01218 double wJ = jacobian_determinant * this->mpBoundaryQuadratureRule->GetWeight(quad_index);
01219
01220 const ChastePoint<DIM-1>& quad_point = this->mpBoundaryQuadratureRule->rGetQuadPoint(quad_index);
01221
01222 QuadraticBasisFunction<DIM-1>::ComputeBasisFunctions(quad_point, phi);
01223
01224
01225
01226 c_vector<double,DIM> traction = zero_vector<double>(DIM);
01227 switch (this->mrProblemDefinition.GetTractionBoundaryConditionType())
01228 {
01229 case FUNCTIONAL_TRACTION:
01230 {
01231 c_vector<double,DIM> X = zero_vector<double>(DIM);
01232 for (unsigned node_index=0; node_index<NUM_NODES_PER_BOUNDARY_ELEMENT; node_index++)
01233 {
01234 X += phi(node_index)*this->mrQuadMesh.GetNode( rBoundaryElement.GetNodeGlobalIndex(node_index) )->rGetLocation();
01235 }
01236 traction = this->mrProblemDefinition.EvaluateTractionFunction(X, this->mCurrentTime);
01237 break;
01238 }
01239 case ELEMENTWISE_TRACTION:
01240 {
01241 traction = this->mrProblemDefinition.rGetElementwiseTractions()[boundaryConditionIndex];
01242 break;
01243 }
01244 default:
01245 NEVER_REACHED;
01246 }
01247
01248
01249 for (unsigned index=0; index<NUM_NODES_PER_BOUNDARY_ELEMENT*DIM; index++)
01250 {
01251 unsigned spatial_dim = index%DIM;
01252 unsigned node_index = (index-spatial_dim)/DIM;
01253
01254 assert(node_index < NUM_NODES_PER_BOUNDARY_ELEMENT);
01255
01256 rBelem(index) -= traction(spatial_dim)
01257 * phi(node_index)
01258 * wJ;
01259 }
01260 }
01261 }
01262
01263 template<unsigned DIM>
01264 bool AbstractNonlinearElasticitySolver<DIM>::ShouldAssembleMatrixTermForPressureOnDeformedBc()
01265 {
01266 if(mUseSnesSolver)
01267 {
01268
01269
01270
01271
01272 return true;
01273
01274
01275
01276
01277
01278 }
01279 else
01280 {
01281 return (mLastDampingValue >= 0.5);
01282 }
01283
01284 }
01285
01286 template<unsigned DIM>
01287 void AbstractNonlinearElasticitySolver<DIM>::AssembleOnBoundaryElementForPressureOnDeformedBc(
01288 BoundaryElement<DIM-1,DIM>& rBoundaryElement,
01289 c_matrix<double,BOUNDARY_STENCIL_SIZE,BOUNDARY_STENCIL_SIZE>& rAelem,
01290 c_vector<double,BOUNDARY_STENCIL_SIZE>& rBelem,
01291 bool assembleResidual,
01292 bool assembleJacobian,
01293 unsigned boundaryConditionIndex)
01294 {
01295 assert( this->mrProblemDefinition.GetTractionBoundaryConditionType()==PRESSURE_ON_DEFORMED
01296 || this->mrProblemDefinition.GetTractionBoundaryConditionType()==FUNCTIONAL_PRESSURE_ON_DEFORMED);
01297
01298 rAelem.clear();
01299 rBelem.clear();
01300
01301 c_vector<double, DIM> weighted_direction;
01302 double jacobian_determinant;
01303
01304 this->mrQuadMesh.GetWeightedDirectionForBoundaryElement(rBoundaryElement.GetIndex(), weighted_direction, jacobian_determinant);
01305
01306
01308
01309
01311
01312 Element<DIM,DIM>* p_containing_vol_element = NULL;
01313
01314 std::set<unsigned> potential_elements = rBoundaryElement.GetNode(0)->rGetContainingElementIndices();
01315 for(std::set<unsigned>::iterator iter = potential_elements.begin();
01316 iter != potential_elements.end();
01317 iter++)
01318 {
01319 p_containing_vol_element = this->mrQuadMesh.GetElement(*iter);
01320
01321 bool this_vol_ele_contains_surf_ele = true;
01322
01323 for(unsigned i=1; i<NUM_NODES_PER_BOUNDARY_ELEMENT; i++)
01324 {
01325 unsigned surf_element_node_index = rBoundaryElement.GetNodeGlobalIndex(i);
01326 bool found_this_node = false;
01327 for(unsigned j=0; j<p_containing_vol_element->GetNumNodes(); j++)
01328 {
01329 unsigned vol_element_node_index = p_containing_vol_element->GetNodeGlobalIndex(j);
01330 if(surf_element_node_index == vol_element_node_index)
01331 {
01332 found_this_node = true;
01333 break;
01334 }
01335 }
01336 if(!found_this_node)
01337 {
01338 this_vol_ele_contains_surf_ele = false;
01339 break;
01340 }
01341 }
01342 if(this_vol_ele_contains_surf_ele)
01343 {
01344 break;
01345 }
01346 }
01347
01348
01349 std::vector<unsigned> surf_to_vol_map(NUM_NODES_PER_BOUNDARY_ELEMENT);
01350 for(unsigned i=0; i<NUM_NODES_PER_BOUNDARY_ELEMENT; i++)
01351 {
01352 unsigned index = rBoundaryElement.GetNodeGlobalIndex(i);
01353 for(unsigned j=0; j<NUM_NODES_PER_ELEMENT; j++)
01354 {
01355 if(p_containing_vol_element->GetNodeGlobalIndex(j)==index)
01356 {
01357 surf_to_vol_map[i] = j;
01358 break;
01359 }
01360 }
01361 }
01362
01363
01364
01365
01366 static c_matrix<double,DIM,DIM> jacobian_vol_element;
01367 static c_matrix<double,DIM,DIM> inverse_jacobian_vol_element;
01368 double jacobian_determinant_vol_element;
01369 this->mrQuadMesh.GetInverseJacobianForElement(p_containing_vol_element->GetIndex(), jacobian_vol_element, jacobian_determinant_vol_element, inverse_jacobian_vol_element);
01370
01371
01372 static c_matrix<double,DIM,NUM_NODES_PER_ELEMENT> element_current_displacements;
01373 for (unsigned II=0; II<NUM_NODES_PER_ELEMENT; II++)
01374 {
01375 for (unsigned JJ=0; JJ<DIM; JJ++)
01376 {
01377 element_current_displacements(JJ,II) = this->mCurrentSolution[this->mProblemDimension*p_containing_vol_element->GetNodeGlobalIndex(II) + JJ];
01378 }
01379 }
01380
01381
01382
01383 static c_matrix<double, DIM, NUM_NODES_PER_ELEMENT> grad_quad_phi_vol_element;
01384
01385 c_vector<double,NUM_NODES_PER_BOUNDARY_ELEMENT> quad_phi_surf_element;
01386
01387 static c_matrix<double, DIM, NUM_NODES_PER_BOUNDARY_ELEMENT> grad_quad_phi_surf_element;
01388
01389 c_matrix<double,DIM,DIM> F;
01390 c_matrix<double,DIM,DIM> invF;
01391
01392 c_vector<double,DIM> normal = rBoundaryElement.CalculateNormal();
01393 c_matrix<double,1,DIM> normal_as_mat;
01394 for(unsigned i=0; i<DIM; i++)
01395 {
01396 normal_as_mat(0,i) = normal(i);
01397 }
01398
01399 double normal_pressure;
01400 switch (this->mrProblemDefinition.GetTractionBoundaryConditionType())
01401 {
01402 case PRESSURE_ON_DEFORMED:
01403 normal_pressure = this->mrProblemDefinition.GetNormalPressure();
01404 break;
01405 case FUNCTIONAL_PRESSURE_ON_DEFORMED:
01406 normal_pressure = this->mrProblemDefinition.EvaluateNormalPressureFunction(this->mCurrentTime);
01407 break;
01408 default:
01409 NEVER_REACHED;
01410 }
01411
01412
01413
01414 for (unsigned quad_index=0; quad_index<this->mpBoundaryQuadratureRule->GetNumQuadPoints(); quad_index++)
01415 {
01416 double wJ = jacobian_determinant * this->mpBoundaryQuadratureRule->GetWeight(quad_index);
01417
01418
01419
01420 const ChastePoint<DIM-1>& quadrature_point = this->mpBoundaryQuadratureRule->rGetQuadPoint(quad_index);
01421 QuadraticBasisFunction<DIM-1>::ComputeBasisFunctions(quadrature_point, quad_phi_surf_element);
01422
01423
01424
01425
01426
01427 c_vector<double,DIM> X = zero_vector<double>(DIM);
01428 for (unsigned node_index=0; node_index<NUM_NODES_PER_BOUNDARY_ELEMENT; node_index++)
01429 {
01430 X += quad_phi_surf_element(node_index)*rBoundaryElement.GetNode(node_index)->rGetLocation();
01431 }
01432
01433
01434
01435 c_vector<double,DIM+1> weight = p_containing_vol_element->CalculateInterpolationWeights(X);
01436 c_vector<double,DIM> xi;
01437 for(unsigned i=0; i<DIM; i++)
01438 {
01439 xi(i) = weight(i+1);
01440 }
01441
01442
01443 if(DIM==2)
01444 {
01445 assert( DIM!=2 || (fabs(weight(0))<1e-6) || (fabs(weight(1))<1e-6) || (fabs(weight(2))<1e-6) );
01446 }
01447 else
01448 {
01449 #define COVERAGE_IGNORE
01450 assert( DIM!=3 || (fabs(weight(0))<1e-6) || (fabs(weight(1))<1e-6) || (fabs(weight(2))<1e-6) || (fabs(weight(3))<1e-6) );
01451 #undef COVERAGE_IGNORE
01452 }
01453
01454
01455 QuadraticBasisFunction<DIM>::ComputeTransformedBasisFunctionDerivatives(xi, inverse_jacobian_vol_element, grad_quad_phi_vol_element);
01456
01457 F = identity_matrix<double>(DIM,DIM);
01458 for (unsigned node_index=0; node_index<NUM_NODES_PER_ELEMENT; node_index++)
01459 {
01460 for (unsigned i=0; i<DIM; i++)
01461 {
01462 for (unsigned M=0; M<DIM; M++)
01463 {
01464 F(i,M) += grad_quad_phi_vol_element(M,node_index)*element_current_displacements(i,node_index);
01465 }
01466 }
01467 }
01468
01469 double detF = Determinant(F);
01470 invF = Inverse(F);
01471
01472 if(assembleResidual)
01473 {
01474 c_vector<double,DIM> traction = detF*normal_pressure*prod(trans(invF),normal);
01475
01476
01477 for (unsigned index=0; index<NUM_NODES_PER_BOUNDARY_ELEMENT*DIM; index++)
01478 {
01479 unsigned spatial_dim = index%DIM;
01480 unsigned node_index = (index-spatial_dim)/DIM;
01481
01482 assert(node_index < NUM_NODES_PER_BOUNDARY_ELEMENT);
01483
01484 rBelem(index) -= traction(spatial_dim)
01485 * quad_phi_surf_element(node_index)
01486 * wJ;
01487 }
01488 }
01489
01490
01491
01492
01493
01494 if(assembleJacobian && ShouldAssembleMatrixTermForPressureOnDeformedBc())
01495 {
01496 for(unsigned II=0; II<NUM_NODES_PER_BOUNDARY_ELEMENT; II++)
01497 {
01498 for(unsigned N=0; N<DIM; N++)
01499 {
01500 grad_quad_phi_surf_element(N,II) = grad_quad_phi_vol_element(N,surf_to_vol_map[II]);
01501 }
01502 }
01503
01504 static FourthOrderTensor<DIM,DIM,DIM,DIM> tensor1;
01505 for(unsigned N=0; N<DIM; N++)
01506 {
01507 for(unsigned e=0; e<DIM; e++)
01508 {
01509 for(unsigned M=0; M<DIM; M++)
01510 {
01511 for(unsigned d=0; d<DIM; d++)
01512 {
01513 tensor1(N,e,M,d) = invF(N,e)*invF(M,d) - invF(M,e)*invF(N,d);
01514 }
01515 }
01516 }
01517 }
01518
01519
01520 static FourthOrderTensor<NUM_NODES_PER_BOUNDARY_ELEMENT,DIM,DIM,DIM> tensor2;
01521 tensor2.template SetAsContractionOnFirstDimension<DIM>( trans(grad_quad_phi_surf_element), tensor1);
01522
01523
01524
01525 static FourthOrderTensor<NUM_NODES_PER_BOUNDARY_ELEMENT,DIM,1,DIM> tensor3;
01526 tensor3.template SetAsContractionOnThirdDimension<DIM>( normal_as_mat, tensor2);
01527
01528 for (unsigned index1=0; index1<NUM_NODES_PER_BOUNDARY_ELEMENT*DIM; index1++)
01529 {
01530 unsigned spatial_dim1 = index1%DIM;
01531 unsigned node_index1 = (index1-spatial_dim1)/DIM;
01532
01533 for (unsigned index2=0; index2<NUM_NODES_PER_BOUNDARY_ELEMENT*DIM; index2++)
01534 {
01535 unsigned spatial_dim2 = index2%DIM;
01536 unsigned node_index2 = (index2-spatial_dim2)/DIM;
01537
01538 rAelem(index1,index2) -= normal_pressure
01539 * detF
01540 * tensor3(node_index2,spatial_dim2,0,spatial_dim1)
01541 * quad_phi_surf_element(node_index1)
01542 * wJ;
01543 }
01544 }
01545 }
01546 }
01547 }
01548
01549
01550
01551
01552 template<unsigned DIM>
01553 void AbstractNonlinearElasticitySolver<DIM>::Solve(double tol)
01554 {
01555
01556
01557 mrProblemDefinition.Validate();
01558
01559
01560
01561
01562
01563
01564 if(mrProblemDefinition.GetTractionBoundaryConditionType()==PRESSURE_ON_DEFORMED && mCheckedOutwardNormals==false)
01565 {
01566 this->mrQuadMesh.CheckOutwardNormals();
01567 mCheckedOutwardNormals = true;
01568 }
01569
01570
01571 this->WriteCurrentSpatialSolution("initial", "nodes");
01572
01573 if(mUseSnesSolver)
01574 {
01575 SolveSnes();
01576 }
01577 else
01578 {
01579 SolveNonSnes(tol);
01580 }
01581
01582
01583
01584
01585 if(this->mCompressibilityType==INCOMPRESSIBLE)
01586 {
01587 this->RemovePressureDummyValuesThroughLinearInterpolation();
01588 }
01589
01590
01591 this->WriteCurrentSpatialSolution("solution", "nodes");
01592
01593 }
01594
01595
01599 template<unsigned DIM>
01600 void AbstractNonlinearElasticitySolver<DIM>::SetKspSolverAndPcType(KSP solver)
01601 {
01602
01603
01604
01605
01606
01607
01608
01609 PC pc;
01610 KSPGetPC(solver, &pc);
01611
01612
01613 if (mPetscDirectSolve)
01614 {
01615 if(this->mCompressibilityType==COMPRESSIBLE)
01616 {
01617 KSPSetType(solver,KSPPREONLY);
01618
01619 }
01620 PCSetType(pc, PCLU);
01621
01622
01623
01624 }
01625 else
01626 {
01627 if (this->mCompressibilityType==COMPRESSIBLE)
01628 {
01629 KSPSetType(solver,KSPCG);
01630 if(PetscTools::IsSequential())
01631 {
01632 PCSetType(pc, PCICC);
01633
01634
01635
01636 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 1) //PETSc 3.1 or later
01637 PetscOptionsSetValue("-pc_factor_shift_type", "positive_definite");
01638 #else
01639 PetscOptionsSetValue("-pc_factor_shift_positive_definite", "");
01640 #endif
01641 }
01642 else
01643 {
01644 PCSetType(pc, PCBJACOBI);
01645 }
01646 }
01647 else
01648 {
01649 unsigned num_restarts = 100;
01650 KSPSetType(solver,KSPGMRES);
01651 KSPGMRESSetRestart(solver,num_restarts);
01652
01653 #ifndef MECH_USE_HYPRE
01654 PCSetType(pc, PCBJACOBI);
01655 #else
01657 // Speed up linear solve time massively for larger simulations (in fact GMRES may stagnate without
01658
01660 PetscOptionsSetValue("-pc_hypre_type", "boomeramg");
01661
01662
01663
01664 PCSetType(pc, PCHYPRE);
01665 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >=2) //PETSc 3.2 or later
01666 KSPSetPCSide(solver, PC_RIGHT);
01667 #else
01668 KSPSetPreconditionerSide(solver, PC_RIGHT);
01669 #endif
01670
01671
01672
01673
01674
01675
01676 #endif
01677 }
01678 }
01679 }
01680
01681
01682
01684
01685
01687
01688 template<unsigned DIM>
01689 double AbstractNonlinearElasticitySolver<DIM>::ComputeResidualAndGetNorm(bool allowException)
01690 {
01691 if (!allowException)
01692 {
01693
01694 AssembleSystem(true, false);
01695 }
01696 else
01697 {
01698 try
01699 {
01700
01701 AssembleSystem(true, false);
01702 }
01703 catch(Exception&)
01704 {
01705
01706 return DBL_MAX;
01707 }
01708 }
01709
01710
01711 return CalculateResidualNorm();
01712 }
01713
01714 template<unsigned DIM>
01715 double AbstractNonlinearElasticitySolver<DIM>::CalculateResidualNorm()
01716 {
01717 double norm;
01718
01719
01720 VecNorm(this->mResidualVector, NORM_2, &norm);
01721 return norm/this->mNumDofs;
01722 }
01723
01724 template<unsigned DIM>
01725 void AbstractNonlinearElasticitySolver<DIM>::VectorSum(std::vector<double>& rX,
01726 ReplicatableVector& rY,
01727 double a,
01728 std::vector<double>& rZ)
01729 {
01730 assert(rX.size()==rY.GetSize());
01731 assert(rY.GetSize()==rZ.size());
01732 for (unsigned i=0; i<rX.size(); i++)
01733 {
01734 rZ[i] = rX[i] + a*rY[i];
01735 }
01736 }
01737
01738
01739 template<unsigned DIM>
01740 double AbstractNonlinearElasticitySolver<DIM>::TakeNewtonStep()
01741 {
01742 if(this->mVerbose)
01743 {
01744 Timer::Reset();
01745 }
01746
01748
01750 MechanicsEventHandler::BeginEvent(MechanicsEventHandler::ASSEMBLE);
01751 AssembleSystem(true, true);
01752 MechanicsEventHandler::EndEvent(MechanicsEventHandler::ASSEMBLE);
01753 if(this->mVerbose)
01754 {
01755 Timer::PrintAndReset("AssembleSystem");
01756 }
01757
01759
01761 MechanicsEventHandler::BeginEvent(MechanicsEventHandler::SOLVE);
01762
01763 Vec solution;
01764 VecDuplicate(this->mResidualVector,&solution);
01765
01766 KSP solver;
01767 KSPCreate(PETSC_COMM_WORLD,&solver);
01768
01769 #if ( PETSC_VERSION_MAJOR==3 && PETSC_VERSION_MINOR>=5 )
01770 KSPSetOperators(solver, mrJacobianMatrix, this->mPreconditionMatrix);
01771 #else
01772 KSPSetOperators(solver, mrJacobianMatrix, this->mPreconditionMatrix, DIFFERENT_NONZERO_PATTERN );
01773 #endif
01774
01775
01776 SetKspSolverAndPcType(solver);
01777
01778
01779
01780
01781 KSPSetFromOptions(solver);
01782 KSPSetUp(solver);
01783
01784
01785
01786
01787
01788 if (mKspAbsoluteTol < 0)
01789 {
01790 Vec temp;
01791 VecDuplicate(this->mResidualVector, &temp);
01792 Vec temp2;
01793 VecDuplicate(this->mResidualVector, &temp2);
01794 Vec linsys_residual;
01795 VecDuplicate(this->mResidualVector, &linsys_residual);
01796
01797 KSPInitialResidual(solver, solution, temp, temp2, linsys_residual, this->mLinearSystemRhsVector);
01798 double initial_resid_norm;
01799 VecNorm(linsys_residual, NORM_2, &initial_resid_norm);
01800
01801 PetscTools::Destroy(temp);
01802 PetscTools::Destroy(temp2);
01803 PetscTools::Destroy(linsys_residual);
01804
01805 double ksp_rel_tol = 1e-6;
01806 double absolute_tol = ksp_rel_tol * initial_resid_norm;
01807 if(absolute_tol < 1e-12)
01808 {
01809 absolute_tol = 1e-12;
01810 }
01811 KSPSetTolerances(solver, 1e-16, absolute_tol, PETSC_DEFAULT, 1000 );
01812 }
01813 else
01814 {
01815 KSPSetTolerances(solver, 1e-16, mKspAbsoluteTol, PETSC_DEFAULT, 1000 );
01816 }
01817
01818 if(this->mVerbose)
01819 {
01820 Timer::PrintAndReset("KSP Setup");
01821 }
01822
01823 KSPSolve(solver,this->mLinearSystemRhsVector,solution);
01824
01825
01826
01827
01828
01829
01830
01831
01832
01833
01834
01835
01836
01837
01838
01839
01840
01841
01842
01843
01845
01847
01848
01849 KSPConvergedReason reason;
01850 KSPGetConvergedReason(solver,&reason);
01851
01852 if(reason != KSP_DIVERGED_ITS)
01853 {
01854
01855
01856
01857 #define COVERAGE_IGNORE
01858 KSPEXCEPT(reason);
01859 #undef COVERAGE_IGNORE
01860 }
01861 else
01862 {
01863
01864
01865
01866
01867
01868
01869 #define COVERAGE_IGNORE
01870 WARN_ONCE_ONLY("Linear solve (within a mechanics solve) didn't converge, but this may not stop nonlinear solve converging")
01871 #undef COVERAGE_IGNORE
01872 }
01873
01874
01875 int num_iters;
01876 KSPGetIterationNumber(solver, &num_iters);
01877 if (num_iters==0)
01878 {
01879 PetscTools::Destroy(solution);
01880 KSPDestroy(PETSC_DESTROY_PARAM(solver));
01881 EXCEPTION("KSP Absolute tolerance was too high, linear system wasn't solved - there will be no decrease in Newton residual. Decrease KspAbsoluteTolerance");
01882 }
01883
01884
01885 if(this->mVerbose)
01886 {
01887 Timer::PrintAndReset("KSP Solve");
01888 std::cout << "[" << PetscTools::GetMyRank() << "]: Num iterations = " << num_iters << "\n" << std::flush;
01889 }
01890
01891 MechanicsEventHandler::EndEvent(MechanicsEventHandler::SOLVE);
01892
01894
01895
01896
01897
01898
01899
01900
01901
01902
01903
01904
01906 MechanicsEventHandler::BeginEvent(MechanicsEventHandler::UPDATE);
01907 double new_norm_resid = UpdateSolutionUsingLineSearch(solution);
01908 MechanicsEventHandler::EndEvent(MechanicsEventHandler::UPDATE);
01909
01910 PetscTools::Destroy(solution);
01911 KSPDestroy(PETSC_DESTROY_PARAM(solver));
01912
01913 return new_norm_resid;
01914 }
01915
01916
01917 template<unsigned DIM>
01918 void AbstractNonlinearElasticitySolver<DIM>::PrintLineSearchResult(double s, double residNorm)
01919 {
01920 if(this->mVerbose)
01921 {
01922 std::cout << "\tTesting s = " << s << ", |f| = " << residNorm << "\n" << std::flush;
01923 }
01924 }
01925
01926 template<unsigned DIM>
01927 double AbstractNonlinearElasticitySolver<DIM>::UpdateSolutionUsingLineSearch(Vec solution)
01928 {
01929 double initial_norm_resid = CalculateResidualNorm();
01930 if(this->mVerbose)
01931 {
01932 std::cout << "\tInitial |f| [corresponding to s=0] is " << initial_norm_resid << "\n" << std::flush;
01933 }
01934
01935 ReplicatableVector update(solution);
01936
01937 std::vector<double> old_solution = this->mCurrentSolution;
01938
01939 std::vector<double> damping_values;
01940 for (unsigned i=10; i>=1; i--)
01941 {
01942 damping_values.push_back((double)i/10.0);
01943 }
01944 damping_values.push_back(0.05);
01945 assert(damping_values.size()==11);
01946
01948
01949 unsigned index = 0;
01950 VectorSum(old_solution, update, -damping_values[index], this->mCurrentSolution);
01951 double current_resid_norm = ComputeResidualAndGetNorm(true);
01952 PrintLineSearchResult(damping_values[index], current_resid_norm);
01953
01955
01956 index = 1;
01957 VectorSum(old_solution, update, -damping_values[index], this->mCurrentSolution);
01958 double next_resid_norm = ComputeResidualAndGetNorm(true);
01959 PrintLineSearchResult(damping_values[index], next_resid_norm);
01960
01961 index = 2;
01962
01963
01964 while ( (next_resid_norm==DBL_MAX)
01965 || ( (next_resid_norm < current_resid_norm) && (index<damping_values.size()) ) )
01966 {
01967 current_resid_norm = next_resid_norm;
01968
01969
01970 VectorSum(old_solution, update, -damping_values[index], this->mCurrentSolution);
01971 next_resid_norm = ComputeResidualAndGetNorm(true);
01972 PrintLineSearchResult(damping_values[index], next_resid_norm);
01973
01974 index++;
01975 }
01976
01977 unsigned best_index;
01978
01979 if (index==damping_values.size() && (next_resid_norm < current_resid_norm))
01980 {
01981
01982
01983
01984
01985
01986 #define COVERAGE_IGNORE
01987
01988
01989 current_resid_norm = next_resid_norm;
01990 best_index = index-1;
01991 #undef COVERAGE_IGNORE
01992 }
01993 else
01994 {
01995
01996
01997 best_index = index-2;
01998 }
01999
02000
02001 bool full_first_step = mTakeFullFirstNewtonStep && mFirstStep;
02002
02003
02004
02005 if (initial_norm_resid < current_resid_norm && !full_first_step)
02006 {
02007 #define COVERAGE_IGNORE
02008 EXCEPTION("Residual does not appear to decrease in newton direction, quitting");
02009 #undef COVERAGE_IGNORE
02010 }
02011
02012
02013 if(full_first_step)
02014 {
02015 if(this->mVerbose)
02016 {
02017 std::cout << "\tTaking full first Newton step...\n";
02018 }
02019 best_index = 0;
02020 }
02021
02022 if(this->mVerbose)
02023 {
02024 std::cout << "\tChoosing s = " << damping_values[best_index] << "\n" << std::flush;
02025 }
02026
02027
02029
02030
02031
02032
02033
02034
02035
02036
02037
02038
02039
02040
02041
02042
02043
02044
02045
02046
02047
02048
02049
02050
02051
02052
02053
02054
02055
02056
02057
02058
02059
02060
02061 VectorSum(old_solution, update, -damping_values[best_index], this->mCurrentSolution);
02062
02063 mLastDampingValue = damping_values[best_index];
02064 return current_resid_norm;
02065 }
02066
02067 template<unsigned DIM>
02068 void AbstractNonlinearElasticitySolver<DIM>::PostNewtonStep(unsigned counter, double normResidual)
02069 {
02070 }
02071
02072
02073 template<unsigned DIM>
02074 void AbstractNonlinearElasticitySolver<DIM>::SolveNonSnes(double tol)
02075 {
02076 mLastDampingValue = 0;
02077
02078 if (mWriteOutputEachNewtonIteration)
02079 {
02080 this->WriteCurrentSpatialSolution("newton_iteration", "nodes", 0);
02081 }
02082
02083
02084 double norm_resid = ComputeResidualAndGetNorm(false);
02085 if(this->mVerbose)
02086 {
02087 std::cout << "\nNorm of residual is " << norm_resid << "\n";
02088 }
02089
02090 mNumNewtonIterations = 0;
02091 unsigned iteration_number = 1;
02092
02093 if (tol < 0)
02094 {
02095 tol = NEWTON_REL_TOL*norm_resid;
02096
02097 #define COVERAGE_IGNORE // not going to have tests in cts for everything
02098 if (tol > MAX_NEWTON_ABS_TOL)
02099 {
02100 tol = MAX_NEWTON_ABS_TOL;
02101 }
02102 if (tol < MIN_NEWTON_ABS_TOL)
02103 {
02104 tol = MIN_NEWTON_ABS_TOL;
02105 }
02106 #undef COVERAGE_IGNORE
02107 }
02108
02109 if(this->mVerbose)
02110 {
02111 std::cout << "Solving with tolerance " << tol << "\n";
02112 }
02113
02114 while (norm_resid > tol)
02115 {
02116 if(this->mVerbose)
02117 {
02118 std::cout << "\n-------------------\n"
02119 << "Newton iteration " << iteration_number
02120 << ":\n-------------------\n";
02121 }
02122
02123
02124 mFirstStep = (iteration_number==1);
02125 norm_resid = TakeNewtonStep();
02126
02127 if(this->mVerbose)
02128 {
02129 std::cout << "Norm of residual is " << norm_resid << "\n";
02130 }
02131
02132 if (mWriteOutputEachNewtonIteration)
02133 {
02134 this->WriteCurrentSpatialSolution("newton_iteration", "nodes", iteration_number);
02135 }
02136
02137 mNumNewtonIterations = iteration_number;
02138
02139 PostNewtonStep(iteration_number,norm_resid);
02140
02141 iteration_number++;
02142 if (iteration_number==20)
02143 {
02144 #define COVERAGE_IGNORE
02145 EXCEPTION("Not converged after 20 newton iterations, quitting");
02146 #undef COVERAGE_IGNORE
02147 }
02148 }
02149
02150 if (norm_resid > tol)
02151 {
02152 #define COVERAGE_IGNORE
02153 EXCEPTION("Failed to converge");
02154 #undef COVERAGE_IGNORE
02155 }
02156 }
02157
02158
02159
02160 template<unsigned DIM>
02161 unsigned AbstractNonlinearElasticitySolver<DIM>::GetNumNewtonIterations()
02162 {
02163 return mNumNewtonIterations;
02164 }
02165
02166
02167
02169
02171
02172
02173 template<unsigned DIM>
02174 void AbstractNonlinearElasticitySolver<DIM>::SolveSnes()
02175 {
02176
02177 Vec initial_guess;
02178 VecDuplicate(this->mResidualVector, &initial_guess);
02179 double* p_initial_guess;
02180 VecGetArray(initial_guess, &p_initial_guess);
02181 int lo, hi;
02182 VecGetOwnershipRange(initial_guess, &lo, &hi);
02183 for (int global_index=lo; global_index<hi; global_index++)
02184 {
02185 int local_index = global_index - lo;
02186 p_initial_guess[local_index] = this->mCurrentSolution[global_index];
02187 }
02188 VecRestoreArray(initial_guess, &p_initial_guess);
02189 PetscVecTools::Finalise(initial_guess);
02190
02191 Vec snes_residual_vec;
02192 VecDuplicate(this->mResidualVector, &snes_residual_vec);
02193
02194 SNES snes;
02195
02196 SNESCreate(PETSC_COMM_WORLD, &snes);
02197 SNESSetFunction(snes, snes_residual_vec, &AbstractNonlinearElasticitySolver_ComputeResidual<DIM>, this);
02198 SNESSetJacobian(snes, mrJacobianMatrix, this->mPreconditionMatrix, &AbstractNonlinearElasticitySolver_ComputeJacobian<DIM>, this);
02199 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 4) //PETSc 3.4 or later
02200 SNESSetType(snes, SNESNEWTONLS);
02201 #else
02202 SNESSetType(snes, SNESLS);
02203 #endif
02204 SNESSetTolerances(snes,1e-5,1e-5,1e-5,PETSC_DEFAULT,PETSC_DEFAULT);
02205 SNESSetMaxLinearSolveFailures(snes,100);
02206
02207 KSP ksp;
02208 SNESGetKSP(snes, &ksp);
02209
02210 KSPSetTolerances(ksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 1000 );
02211
02212
02213 SetKspSolverAndPcType(ksp);
02214
02215 if(this->mVerbose)
02216 {
02217 PetscOptionsSetValue("-snes_monitor","");
02218 }
02219 SNESSetFromOptions(snes);
02220
02221 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
02222 PetscErrorCode err = SNESSolve(snes, initial_guess);
02223 #else
02224 PetscErrorCode err = SNESSolve(snes, PETSC_NULL, initial_guess);
02225 #endif
02226
02227 SNESConvergedReason reason;
02228 SNESGetConvergedReason(snes,&reason);
02229
02230 #define COVERAGE_IGNORE
02231 if (err != 0)
02232 {
02233 std::stringstream err_stream;
02234 err_stream << err;
02235 PetscTools::Destroy(initial_guess);
02236 PetscTools::Destroy(snes_residual_vec);
02237 SNESDestroy(PETSC_DESTROY_PARAM(snes));
02238 EXCEPTION("Nonlinear Solver failed. PETSc error code: "+err_stream.str()+" .");
02239 }
02240
02241 if (reason < 0)
02242 {
02243 std::stringstream reason_stream;
02244 reason_stream << reason;
02245 PetscTools::Destroy(initial_guess);
02246 PetscTools::Destroy(snes_residual_vec);
02247 SNESDestroy(PETSC_DESTROY_PARAM(snes));
02248 EXCEPTION("Nonlinear Solver did not converge. PETSc reason code: "+reason_stream.str()+" .");
02249 }
02250 #undef COVERAGE_IGNORE
02251
02252 PetscInt num_iters;
02253 SNESGetIterationNumber(snes,&num_iters);
02254 mNumNewtonIterations = num_iters;
02255
02256 PetscTools::Destroy(initial_guess);
02257 PetscTools::Destroy(snes_residual_vec);
02258 SNESDestroy(PETSC_DESTROY_PARAM(snes));
02259 }
02260
02261
02262 template<unsigned DIM>
02263 void AbstractNonlinearElasticitySolver<DIM>::ComputeResidual(Vec currentGuess, Vec residualVector)
02264 {
02265
02266
02267
02268
02269 ReplicatableVector guess_repl(currentGuess);
02270 for(unsigned i=0; i<guess_repl.GetSize(); i++)
02271 {
02272 this->mCurrentSolution[i] = guess_repl[i];
02273 }
02274 AssembleSystem(true,false);
02275 VecCopy(this->mResidualVector, residualVector);
02276 }
02277
02278 template<unsigned DIM>
02279 void AbstractNonlinearElasticitySolver<DIM>::ComputeJacobian(Vec currentGuess, Mat* pJacobian, Mat* pPreconditioner)
02280 {
02281
02282
02283
02284
02285
02286
02287
02288 assert(mrJacobianMatrix==*pJacobian);
02289 assert(this->mPreconditionMatrix==*pPreconditioner);
02290
02291 MechanicsEventHandler::BeginEvent(MechanicsEventHandler::ASSEMBLE);
02292 ReplicatableVector guess_repl(currentGuess);
02293 for(unsigned i=0; i<guess_repl.GetSize(); i++)
02294 {
02295 this->mCurrentSolution[i] = guess_repl[i];
02296 }
02297
02298 AssembleSystem(false,true);
02299 MechanicsEventHandler::EndEvent(MechanicsEventHandler::ASSEMBLE);
02300 }
02301
02302
02303
02304 template<unsigned DIM>
02305 PetscErrorCode AbstractNonlinearElasticitySolver_ComputeResidual(SNES snes,
02306 Vec currentGuess,
02307 Vec residualVector,
02308 void* pContext)
02309 {
02310
02311 AbstractNonlinearElasticitySolver<DIM>* p_solver = (AbstractNonlinearElasticitySolver<DIM>*)pContext;
02312 p_solver->ComputeResidual(currentGuess, residualVector);
02313 return 0;
02314 }
02315
02316 template<unsigned DIM>
02317 #if ( PETSC_VERSION_MAJOR==3 && PETSC_VERSION_MINOR>=5 )
02318 PetscErrorCode AbstractNonlinearElasticitySolver_ComputeJacobian(SNES snes,
02319 Vec currentGuess,
02320 Mat globalJacobian,
02321 Mat preconditioner,
02322 void* pContext)
02323 {
02324
02325 AbstractNonlinearElasticitySolver<DIM>* p_solver = (AbstractNonlinearElasticitySolver<DIM>*) pContext;
02326 p_solver->ComputeJacobian(currentGuess, &globalJacobian, &preconditioner);
02327 return 0;
02328 }
02329 #else
02330 PetscErrorCode AbstractNonlinearElasticitySolver_ComputeJacobian(SNES snes,
02331 Vec currentGuess,
02332 Mat* pGlobalJacobian,
02333 Mat* pPreconditioner,
02334 MatStructure* pMatStructure,
02335 void* pContext)
02336 {
02337
02338 AbstractNonlinearElasticitySolver<DIM>* p_solver = (AbstractNonlinearElasticitySolver<DIM>*) pContext;
02339 p_solver->ComputeJacobian(currentGuess, pGlobalJacobian, pPreconditioner);
02340 return 0;
02341 }
02342 #endif
02343
02344
02345
02346
02347 template<unsigned DIM>
02348 double AbstractNonlinearElasticitySolver<DIM>::MAX_NEWTON_ABS_TOL = 1e-7;
02349
02350 template<unsigned DIM>
02351 double AbstractNonlinearElasticitySolver<DIM>::MIN_NEWTON_ABS_TOL = 1e-10;
02352
02353 template<unsigned DIM>
02354 double AbstractNonlinearElasticitySolver<DIM>::NEWTON_REL_TOL = 1e-4;
02355
02356 #endif