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 #include "AbstractMaterialLaw.hpp"
00045 #include "Warnings.hpp"
00046 
00047 
00048 
00049 
00050 
00051 //#define MECH_VERBOSE      // Print output on how nonlinear solve is progressing
00052 //#define MECH_VERY_VERBOSE // See number of elements done whilst assembling vectors or matrices
00053 //#define MECH_USE_HYPRE    // uses HYPRE to solve linear systems, requires Petsc to be installed with HYPRE
00054 //#define MECH_KSP_MONITOR  // Print residual norm each iteration in linear solve (ie -ksp_monitor).
00055 
00056 
00057 
00058 //EMTODO: allow change of max_iter, better preconditioner
00059 
00060 
00061 #ifdef MECH_VERBOSE
00062 #include "Timer.hpp"
00063 #endif
00064 
00065 typedef enum CompressibilityType_
00066 {
00067     COMPRESSIBLE,
00068     INCOMPRESSIBLE
00069 } CompressibilityType;
00070 
00071 
00075 template<unsigned DIM>
00076 class AbstractNonlinearElasticitySolver
00077 {
00078 protected:
00080     static const size_t NUM_VERTICES_PER_ELEMENT = DIM+1;
00082     static const size_t NUM_NODES_PER_ELEMENT = (DIM+1)*(DIM+2)/2;      // assuming quadratic
00084     static const size_t NUM_NODES_PER_BOUNDARY_ELEMENT = DIM*(DIM+1)/2; // assuming quadratic
00085 
00086 
00092     static double MAX_NEWTON_ABS_TOL;
00093 
00095     static double MIN_NEWTON_ABS_TOL;
00096 
00098     static double NEWTON_REL_TOL;
00099 
00100 
00105     QuadraticMesh<DIM>* mpQuadMesh;
00106 
00108     std::vector<BoundaryElement<DIM-1,DIM>*> mBoundaryElements;
00109 
00111     GaussianQuadratureRule<DIM>* mpQuadratureRule;
00112 
00114     GaussianQuadratureRule<DIM-1>* mpBoundaryQuadratureRule;
00115     
00116 
00120     double mKspAbsoluteTol;
00121 
00127     unsigned mNumDofs;
00128 
00129 
00134     Vec mResidualVector;
00135 
00152     Vec mLinearSystemRhsVector;
00153 
00157     Mat mJacobianMatrix;
00158 
00162     Vec mDirichletBoundaryConditionsVector;
00163 
00164 
00183     Mat mPreconditionMatrix;
00184 
00186     c_vector<double,DIM> mBodyForce;
00187 
00189     double mDensity;
00190 
00192     std::vector<unsigned> mFixedNodes;
00193 
00195     std::vector<c_vector<double,DIM> > mFixedNodeDisplacements;
00196 
00198     bool mWriteOutput;
00199 
00201     std::string mOutputDirectory;
00202 
00204     OutputFileHandler* mpOutputFileHandler;
00205 
00209     bool mWriteOutputEachNewtonIteration;
00210 
00211 
00218     std::vector<double> mCurrentSolution;
00219 
00223     FourthOrderTensor<DIM,DIM,DIM,DIM> dTdE;
00224 
00226     unsigned mNumNewtonIterations;
00227 
00229     std::vector<c_vector<double,DIM> > mDeformedPosition;
00230 
00231 
00235     std::vector<c_vector<double,DIM> > mSurfaceTractions;
00236 
00238     c_vector<double,DIM> (*mpBodyForceFunction)(c_vector<double,DIM>& X, double t);
00239 
00244     c_vector<double,DIM> (*mpTractionBoundaryConditionFunction)(c_vector<double,DIM>& X, double t);
00245 
00247     bool mUsingBodyForceFunction;
00248 
00250     bool mUsingTractionBoundaryConditionFunction;
00251     
00256     double mCurrentTime;
00257 
00261     CompressibilityType mCompressibilityType;
00262 
00275     virtual void AssembleSystem(bool assembleResidual, bool assembleLinearSystem)=0;
00276 
00277 
00278 
00284     void Initialise(std::vector<c_vector<double,DIM> >* pFixedNodeLocations);
00285 
00286 
00290     void AllocateMatrixMemory();
00291 
00305     void ApplyBoundaryConditions(bool applyToMatrix);
00306 
00314     void FinishAssembleSystem(bool assembleResidual, bool assembleLinearSystem);
00315 
00326     double ComputeResidualAndGetNorm(bool allowException);
00327 
00329     double CalculateResidualNorm();
00330 
00338     void VectorSum(std::vector<double>& rX, ReplicatableVector& rY, double a, std::vector<double>& rZ);
00339 
00346     void PrintLineSearchResult(double s, double residNorm);
00347 
00348 
00356     double TakeNewtonStep();
00357 
00365     double UpdateSolutionUsingLineSearch(Vec solution);
00366 
00367 
00375     virtual void PostNewtonStep(unsigned counter, double normResidual);
00376 
00377 
00396     virtual void ComputeStressAndStressDerivative(AbstractMaterialLaw<DIM>* pMaterialLaw,
00397                                                   c_matrix<double,DIM,DIM>& rC,
00398                                                   c_matrix<double,DIM,DIM>& rInvC,
00399                                                   double pressure,
00400                                                   unsigned elementIndex,
00401                                                   unsigned currentQuadPointGlobalIndex,
00402                                                   c_matrix<double,DIM,DIM>& rT,
00403                                                   FourthOrderTensor<DIM,DIM,DIM,DIM>& rDTdE,
00404                                                   bool computeDTdE)
00405     {
00406         // just call the method on the material law
00407         pMaterialLaw->ComputeStressAndStressDerivative(rC,rInvC,pressure,rT,rDTdE,computeDTdE);
00408     }
00409 
00410 public:
00411 
00424     AbstractNonlinearElasticitySolver(QuadraticMesh<DIM>* pQuadMesh,
00425                                       c_vector<double,DIM> bodyForce,
00426                                       double density,
00427                                       std::string outputDirectory,
00428                                       std::vector<unsigned>& fixedNodes,
00429                                       CompressibilityType compressibilityType);
00430 
00431 
00435     virtual ~AbstractNonlinearElasticitySolver();
00436 
00444     void Solve(double tol=-1.0,
00445                unsigned maxNumNewtonIterations=INT_MAX,
00446                bool quitIfNoConvergence=true);
00447 
00458     void WriteCurrentDeformation(std::string fileName, int counterToAppend=-1);
00459 
00463     unsigned GetNumNewtonIterations();
00464 
00472     void SetFunctionalBodyForce(c_vector<double,DIM> (*pFunction)(c_vector<double,DIM>& X, double t));
00473 
00479     void SetWriteOutput(bool writeOutput=true);
00480 
00487     void SetWriteOutputEachNewtonIteration(bool writeOutputEachNewtonIteration=true)
00488     {
00489         mWriteOutputEachNewtonIteration = writeOutputEachNewtonIteration;
00490     }
00491 
00497     void SetKspAbsoluteTolerance(double kspAbsoluteTolerance)
00498     {
00499         assert(kspAbsoluteTolerance>0);
00500         mKspAbsoluteTol = kspAbsoluteTolerance;
00501     }
00502 
00507     std::vector<double>& rGetCurrentSolution()
00508     {
00509         return mCurrentSolution;
00510     }
00511 
00520     void SetSurfaceTractionBoundaryConditions(std::vector<BoundaryElement<DIM-1,DIM>*>& rBoundaryElements,
00521                                               std::vector<c_vector<double,DIM> >& rSurfaceTractions);
00522 
00531     void SetFunctionalTractionBoundaryCondition(std::vector<BoundaryElement<DIM-1,DIM>*> rBoundaryElements,
00532                                                 c_vector<double,DIM> (*pFunction)(c_vector<double,DIM>& X, double t));
00533 
00537     std::vector<c_vector<double,DIM> >& rGetDeformedPosition();
00538 
00543     void SetCurrentTime(double time)
00544     {
00545         mCurrentTime = time;
00546     }
00547 
00551     void CreateCmguiOutput();
00552 };
00553 
00554 
00556 // Implementation
00558 
00559 
00560 template<unsigned DIM>
00561 void AbstractNonlinearElasticitySolver<DIM>::Initialise(std::vector<c_vector<double,DIM> >* pFixedNodeLocations)
00562 {
00563     assert(mpQuadMesh);
00564 
00565     AllocateMatrixMemory();
00566 
00567     for (unsigned i=0; i<mFixedNodes.size(); i++)
00568     {
00569         assert(mFixedNodes[i] < mpQuadMesh->GetNumNodes());
00570     }
00571 
00572     mpQuadratureRule = new GaussianQuadratureRule<DIM>(3);
00573     mpBoundaryQuadratureRule = new GaussianQuadratureRule<DIM-1>(3);
00574 
00575     mCurrentSolution.resize(mNumDofs, 0.0);
00576 
00577     // compute the displacements at each of the fixed nodes, given the
00578     // fixed nodes locations.
00579     if (pFixedNodeLocations == NULL)
00580     {
00581         mFixedNodeDisplacements.clear();
00582         for (unsigned i=0; i<mFixedNodes.size(); i++)
00583         {
00584             mFixedNodeDisplacements.push_back(zero_vector<double>(DIM));
00585         }
00586     }
00587     else
00588     {
00589         assert(pFixedNodeLocations->size()==mFixedNodes.size());
00590         for (unsigned i=0; i<mFixedNodes.size(); i++)
00591         {
00592             unsigned index = mFixedNodes[i];
00593             c_vector<double,DIM> displacement = (*pFixedNodeLocations)[i] - mpQuadMesh->GetNode(index)->rGetLocation();
00594             mFixedNodeDisplacements.push_back(displacement);
00595         }
00596     }
00597     assert(mFixedNodeDisplacements.size()==mFixedNodes.size());
00598 }
00599 
00600 
00601 
00602 template<unsigned DIM>
00603 void AbstractNonlinearElasticitySolver<DIM>::AllocateMatrixMemory()
00604 {
00605     mResidualVector = PetscTools::CreateVec(mNumDofs);
00606     VecDuplicate(mResidualVector, &mLinearSystemRhsVector);
00607     if(mCompressibilityType==COMPRESSIBLE)
00608     {
00609         VecDuplicate(mResidualVector, &mDirichletBoundaryConditionsVector);
00610     }
00611 
00612     if(DIM==2)
00613     {
00614         // 2D: N elements around a point => 7N+3 non-zeros in that row? Assume N<=10 (structured mesh would have N_max=6) => 73.
00615         unsigned num_non_zeros = std::min(75u, mNumDofs);
00616 
00617         PetscTools::SetupMat(mJacobianMatrix, mNumDofs, mNumDofs, num_non_zeros, PETSC_DECIDE, PETSC_DECIDE);
00618         PetscTools::SetupMat(mPreconditionMatrix, mNumDofs, mNumDofs, num_non_zeros, PETSC_DECIDE, PETSC_DECIDE);
00619 
00620     }
00621     else
00622     {
00623         assert(DIM==3);
00624 
00625         // 3D: N elements around a point. nz < (3*10+6)N (lazy estimate). Better estimate is 23N+4?. Assume N<20 => 500ish
00626 
00627         // in 3d we get the number of containing elements for each node and use that to obtain an upper bound
00628         // for the number of non-zeros for each DOF associated with that node.
00629 
00630         int* num_non_zeros_each_row = new int[mNumDofs];
00631         for(unsigned i=0; i<mNumDofs; i++)
00632         {
00633             num_non_zeros_each_row[i] = 0;
00634         }
00635 
00636         for(unsigned i=0; i<mpQuadMesh->GetNumNodes(); i++)
00637         {
00638             // this upper bound neglects the fact that two containing elements will share the same nodes..
00639             // 4 = max num dofs associated with this node
00640             // 30 = 3*9+3 = 3 dimensions x 9 other nodes on this element   +  3 vertices with a pressure unknown
00641             unsigned num_non_zeros_upper_bound = 4 + 30*mpQuadMesh->GetNode(i)->GetNumContainingElements();
00642 
00643             num_non_zeros_upper_bound = std::min(num_non_zeros_upper_bound, mNumDofs);
00644 
00645             num_non_zeros_each_row[DIM*i + 0] = num_non_zeros_upper_bound;
00646             num_non_zeros_each_row[DIM*i + 1] = num_non_zeros_upper_bound;
00647             num_non_zeros_each_row[DIM*i + 2] = num_non_zeros_upper_bound;
00648 
00649             if(mCompressibilityType==INCOMPRESSIBLE)
00650             {
00651                 //Could do !mpQuadMesh->GetNode(i)->IsInternal()
00652                 if(i<mpQuadMesh->GetNumVertices()) // then this is a vertex
00653                 {
00654                     num_non_zeros_each_row[DIM*mpQuadMesh->GetNumNodes() + i] = num_non_zeros_upper_bound;
00655                 }
00656             }
00657         }
00658 
00659         // NOTE: PetscTools::SetupMat() or the below creates a MATAIJ matrix, which means the matrix will
00660         // be of type MATSEQAIJ if num_procs=1 and MATMPIAIJ otherwise. In the former case
00661         // MatSeqAIJSetPreallocation MUST be called [MatMPIAIJSetPreallocation will have
00662         // no effect (silently)], and vice versa in the latter case
00663 
00666         //PetscTools::SetupMat(mJacobianMatrix, mNumDofs, mNumDofs, 0, PETSC_DECIDE, PETSC_DECIDE);
00667         //PetscTools::SetupMat(mPreconditionMatrix, mNumDofs, mNumDofs, 0, PETSC_DECIDE, PETSC_DECIDE);
00669 
00670         // possible todo: create a PetscTools::SetupMatNoAllocation()
00671 
00672         // In the future, when parallelising, remember to think about MAT_IGNORE_OFF_PROC_ENTRIES (see #1682)
00673 
00674 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00675         MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,mNumDofs,mNumDofs,&mJacobianMatrix);
00676         MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,mNumDofs,mNumDofs,&mPreconditionMatrix);
00677 #else //New API
00678         MatCreate(PETSC_COMM_WORLD,&mJacobianMatrix);
00679         MatCreate(PETSC_COMM_WORLD,&mPreconditionMatrix);
00680         MatSetSizes(mJacobianMatrix,PETSC_DECIDE,PETSC_DECIDE,mNumDofs,mNumDofs);
00681         MatSetSizes(mPreconditionMatrix,PETSC_DECIDE,PETSC_DECIDE,mNumDofs,mNumDofs);
00682 #endif
00683 
00684         if(PetscTools::IsSequential())
00685         {
00686             MatSetType(mJacobianMatrix, MATSEQAIJ);
00687             MatSetType(mPreconditionMatrix, MATSEQAIJ);
00688             MatSeqAIJSetPreallocation(mJacobianMatrix,     PETSC_NULL, num_non_zeros_each_row);
00689             MatSeqAIJSetPreallocation(mPreconditionMatrix, PETSC_NULL, num_non_zeros_each_row);
00690         }
00691         else
00692         {
00693             PetscInt lo, hi;
00694             VecGetOwnershipRange(mResidualVector, &lo, &hi);
00695 
00696             int* num_non_zeros_each_row_this_proc = new int[hi-lo];
00697             int* zero = new int[hi-lo];
00698             for(unsigned i=0; i<unsigned(hi-lo); i++)
00699             {
00700                 num_non_zeros_each_row_this_proc[i] = num_non_zeros_each_row[lo+i];
00701                 zero[i] = 0;
00702             }
00703 
00704             MatSetType(mJacobianMatrix, MATMPIAIJ);
00705             MatSetType(mPreconditionMatrix, MATMPIAIJ);
00706             MatMPIAIJSetPreallocation(mJacobianMatrix,     PETSC_NULL, num_non_zeros_each_row_this_proc, PETSC_NULL, num_non_zeros_each_row_this_proc);
00707             MatMPIAIJSetPreallocation(mPreconditionMatrix, PETSC_NULL, num_non_zeros_each_row_this_proc, PETSC_NULL, num_non_zeros_each_row_this_proc);
00708         }
00709 
00710         MatSetFromOptions(mJacobianMatrix);
00711         MatSetFromOptions(mPreconditionMatrix);
00712 
00713 
00714         //unsigned total_non_zeros = 0;
00715         //for(unsigned i=0; i<mNumDofs; i++)
00716         //{
00717         //   total_non_zeros += num_non_zeros_each_row[i];
00718         //}
00719         //std::cout << total_non_zeros << " versus " << 500*mNumDofs << "\n" << std::flush;
00720 
00721         delete [] num_non_zeros_each_row;
00722     }
00723 }
00724 
00725 /*
00726  * This method applies the appropriate BCs to the residual vector, and possible also the linear system.
00727  *
00728  * For the latter, in the compressible case, the BCs are imposed in such a way as to ensure that a
00729  * symmetric linear system remains symmetric. For each row with boundary condition applied, both the
00730  * row and column are zero'd and the RHS vector modified to take into account the zero'd column.
00731  *
00732  * Suppose we have a matrix
00733  * [a b c] [x] = [ b1 ]
00734  * [d e f] [y]   [ b2 ]
00735  * [g h i] [z]   [ b3 ]
00736  * and we want to apply the boundary condition x=v without losing symmetry if the matrix is
00737  * symmetric. We apply the boundary condition
00738  * [1 0 0] [x] = [ v  ]
00739  * [d e f] [y]   [ b2 ]
00740  * [g h i] [z]   [ b3 ]
00741  * and then zero the column as well, adding a term to the RHS to take account for the
00742  * zero-matrix components
00743  * [1 0 0] [x] = [ v  ] - v[ 0 ]
00744  * [0 e f] [y]   [ b2 ]    [ d ]
00745  * [0 h i] [z]   [ b3 ]    [ g ]
00746  * Note the last term is the first column of the matrix, with one component zeroed, and
00747  * multiplied by the boundary condition value. This last term is then stored in
00748  * rLinearSystem.rGetDirichletBoundaryConditionsVector(), and in general form is the
00749  * SUM_{d=1..D} v_d a'_d
00750  * where v_d is the boundary value of boundary condition d (d an index into the matrix),
00751  * and a'_d is the dth-column of the matrix but with the d-th component zeroed, and where
00752  * there are D boundary conditions
00753  */
00754 template<unsigned DIM>
00755 void AbstractNonlinearElasticitySolver<DIM>::ApplyBoundaryConditions(bool applyToLinearSystem)
00756 {
00757     assert(mFixedNodeDisplacements.size()==mFixedNodes.size());
00758 
00759     assert(mResidualVector); // BCs will be added to all the time
00760     if(applyToLinearSystem)
00761     {
00762         assert(mJacobianMatrix);
00763         assert(mLinearSystemRhsVector);
00764     }
00765 
00766     // The boundary conditions on the NONLINEAR SYSTEM are x=boundary_values
00767     // on the boundary nodes. However:
00768     // The boundary conditions on the LINEAR SYSTEM  Ju=f, where J is the
00769     // u the negative update vector and f is the residual is
00770     // u=current_soln-boundary_values on the boundary nodes
00771 
00772     std::vector<unsigned> rows;
00773     std::vector<double> values;
00774 
00775     rows.resize(DIM*mFixedNodes.size());
00776     values.resize(DIM*mFixedNodes.size());
00777 
00778     // whether to apply symmetrically, ie alter columns as well as rows (see comment above)
00779     bool applySymmetrically = (applyToLinearSystem) && (mCompressibilityType==COMPRESSIBLE);
00780 
00781     if(applySymmetrically)
00782     {
00783         assert(applyToLinearSystem);
00784         PetscVecTools::Zero(mDirichletBoundaryConditionsVector);
00785         PetscMatTools::AssembleFinal(mJacobianMatrix);
00786     }
00787 
00788 
00789     for (unsigned i=0; i<mFixedNodes.size(); i++)
00790     {
00791         unsigned node_index = mFixedNodes[i];
00792         for (unsigned j=0; j<DIM; j++)
00793         {
00794             unsigned dof_index = DIM*node_index+j;
00795             rows[DIM*i+j] = dof_index;
00796             values[DIM*i+j] = mCurrentSolution[dof_index] - mFixedNodeDisplacements[i](j);;
00797         }
00798     }
00799 
00800 
00801     if(applySymmetrically)
00802     {
00803         // Modify the matrix columns
00804         for (unsigned i=0; i<rows.size(); i++)
00805         {
00806             unsigned col = rows[i];
00807             double minus_value = -values[i];
00808 
00809             // Get a vector which will store the column of the matrix (column d, where d is
00810             // the index of the row (and column) to be altered for the boundary condition.
00811             // Since the matrix is symmetric when get row number "col" and treat it as a column.
00812             // PETSc uses compressed row format and therefore getting rows is far more efficient
00813             // than getting columns.
00814             Vec matrix_col = PetscMatTools::GetMatrixRowDistributed(mJacobianMatrix,col);
00815 
00816             // Zero the correct entry of the column
00817             PetscVecTools::SetElement(matrix_col, col, 0.0);
00818 
00819             // Set up the RHS Dirichlet boundary conditions vector
00820             // Assuming one boundary at the zeroth node (x_0 = value), this is equal to
00821             //   -value*[0 a_21 a_31 .. a_N1]
00822             // and will be added to the RHS.
00823             PetscVecTools::AddScaledVector(mDirichletBoundaryConditionsVector, matrix_col, minus_value);
00824             VecDestroy(matrix_col);
00825         }
00826     }
00827 
00828     if(applyToLinearSystem)
00829     {
00830         // Now zero the appropriate rows and columns of the matrix. If the matrix is symmetric we apply the
00831         // boundary conditions in a way the symmetry isn't lost (rows and columns). If not only the row is
00832         // zeroed.
00833         if (applySymmetrically)
00834         {
00835             PetscMatTools::ZeroRowsAndColumnsWithValueOnDiagonal(mJacobianMatrix, rows, 1.0);
00836             PetscMatTools::ZeroRowsAndColumnsWithValueOnDiagonal(mPreconditionMatrix, rows, 1.0);
00837 
00838             // Apply the RHS boundary conditions modification if required.
00839             PetscVecTools::AddScaledVector(mLinearSystemRhsVector, mDirichletBoundaryConditionsVector, 1.0);
00840         }
00841         else
00842         {
00843             PetscMatTools::ZeroRowsWithValueOnDiagonal(mJacobianMatrix, rows, 1.0);
00844             PetscMatTools::ZeroRowsWithValueOnDiagonal(mPreconditionMatrix, rows, 1.0);
00845         }
00846 
00847     }
00848 
00849     for (unsigned i=0; i<rows.size(); i++)
00850     {
00851         PetscVecTools::SetElement(mResidualVector, rows[i], values[i]);
00852     }
00853 
00854     if(applyToLinearSystem)
00855     {
00856         for (unsigned i=0; i<rows.size(); i++)
00857         {
00858             PetscVecTools::SetElement(mLinearSystemRhsVector, rows[i], values[i]);
00859         }
00860     }
00861 }
00862 
00863 
00864 
00865 template<unsigned DIM>
00866 void AbstractNonlinearElasticitySolver<DIM>::FinishAssembleSystem(bool assembleResidual, bool assembleJacobian)
00867 {
00868     if (assembleResidual)
00869     {
00870         PetscVecTools::Assemble(mResidualVector);
00871     }
00872     if (assembleJacobian)
00873     {
00874         PetscMatTools::AssembleIntermediate(mJacobianMatrix);
00875         PetscMatTools::AssembleIntermediate(mPreconditionMatrix);
00876 
00877         VecCopy(mResidualVector, mLinearSystemRhsVector);
00878     }
00879 
00880     // Apply Dirichlet boundary conditions
00881     ApplyBoundaryConditions(assembleJacobian);
00882 
00883     if (assembleResidual)
00884     {
00885         PetscVecTools::Assemble(mResidualVector);
00886     }
00887     if (assembleJacobian)
00888     {
00889         PetscMatTools::AssembleFinal(mJacobianMatrix);
00890         PetscMatTools::AssembleFinal(mPreconditionMatrix);
00891         PetscVecTools::Assemble(mLinearSystemRhsVector);
00892     }
00893 }
00894 
00895 template<unsigned DIM>
00896 double AbstractNonlinearElasticitySolver<DIM>::ComputeResidualAndGetNorm(bool allowException)
00897 {
00898     if(!allowException)
00899     {
00900         // assemble the residual
00901         AssembleSystem(true, false);
00902     }
00903     else
00904     {
00905         try
00906         {
00907             // try to assemble the residual using this current solution
00908             AssembleSystem(true, false);
00909         }
00910         catch(Exception& e)
00911         {
00912             // if fail (because eg ODEs fail to solve, or strains are too large for material law), 
00913             // return infinity
00914             return DBL_MAX;
00915         }
00916     }
00917 
00918     // return the scaled norm of the residual
00919     return CalculateResidualNorm();
00920 }
00921 
00922 template<unsigned DIM>
00923 double AbstractNonlinearElasticitySolver<DIM>::CalculateResidualNorm()
00924 {
00925     double norm;
00926     VecNorm(mResidualVector, NORM_2, &norm);
00927     return norm/mNumDofs;
00928 }
00929 
00930 template<unsigned DIM>
00931 void AbstractNonlinearElasticitySolver<DIM>::VectorSum(std::vector<double>& rX,
00932                                                        ReplicatableVector& rY,
00933                                                        double a,
00934                                                        std::vector<double>& rZ)
00935 {
00936     assert(rX.size()==rY.GetSize());
00937     assert(rY.GetSize()==rZ.size());
00938     for(unsigned i=0; i<rX.size(); i++)
00939     {
00940         rZ[i] = rX[i] + a*rY[i];
00941     }
00942 }
00943 
00944 
00945 template<unsigned DIM>
00946 double AbstractNonlinearElasticitySolver<DIM>::TakeNewtonStep()
00947 {
00948     #ifdef MECH_VERBOSE
00949     Timer::Reset();
00950     #endif
00951 
00953     // Assemble Jacobian (and preconditioner)
00955     MechanicsEventHandler::BeginEvent(MechanicsEventHandler::ASSEMBLE);
00956     AssembleSystem(true, true);
00957     MechanicsEventHandler::EndEvent(MechanicsEventHandler::ASSEMBLE);
00958     #ifdef MECH_VERBOSE
00959     Timer::PrintAndReset("AssembleSystem");
00960     #endif
00961 
00962 
00964     //
00965     // Solve the linear system.
00966     // Three alternatives
00967     //   (a) GMRES with ILU preconditioner (or bjacobi=ILU on each process) [default]. Very poor on large problems.
00968     //   (b) GMRES with AMG preconditioner. Uncomment #define MECH_USE_HYPRE above. Requires Petsc3 with HYPRE installed.
00969     //
00971     MechanicsEventHandler::BeginEvent(MechanicsEventHandler::SOLVE);
00972 
00973     Vec solution;
00974     VecDuplicate(mResidualVector,&solution);
00975 
00976     KSP solver;
00977     KSPCreate(PETSC_COMM_WORLD,&solver);
00978     PC pc;
00979     KSPGetPC(solver, &pc);
00980     
00982     // // no need for the precondition matrix
00983     // KSPSetOperators(solver, mJacobianMatrix, mJacobianMatrix, DIFFERENT_NONZERO_PATTERN /*in precond between successive solves*/);
00984     // KSPSetType(solver, KSPPREONLY);
00985     // PCSetType(pc, PCLU);
00986     // PCFactorSetMatSolverPackage(pc, MAT_SOLVER_MUMPS);
00987     // // Allow matrices with zero diagonals to be solved
00988     // PCFactorSetShiftNonzero(pc, PETSC_DECIDE);
00989     // // "Do LU factorization in-place (saves memory)"
00990     // PCASMSetUseInPlace(pc);
00992 
00993     KSPSetOperators(solver, mJacobianMatrix, mPreconditionMatrix, DIFFERENT_NONZERO_PATTERN /*in precond between successive solves*/);
00994 
00995     if(mCompressibilityType==COMPRESSIBLE)
00996     {
00997         KSPSetType(solver,KSPCG);
00998     }
00999     else
01000     {
01001         unsigned num_restarts = 100;
01002         KSPSetType(solver,KSPGMRES);
01003         KSPGMRESSetRestart(solver,num_restarts);
01004     }
01005 
01006     if(mKspAbsoluteTol < 0)
01007     {
01008         double ksp_rel_tol = 1e-6;
01009         KSPSetTolerances(solver, ksp_rel_tol, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT /* max iters */); // Note, max iters seems to be 1000 whatever we give here
01010     }
01011     else
01012     {
01013         KSPSetTolerances(solver, 1e-16, mKspAbsoluteTol, PETSC_DEFAULT, PETSC_DEFAULT /* max iters */); // Note, max iters seems to be 1000 whatever we give here
01014     }
01015 
01016     #ifndef MECH_USE_HYPRE
01017         PCSetType(pc, PCBJACOBI); // BJACOBI = ILU on each block (block = part of matrix on each process)
01018     #else
01020         // Speed up linear solve time massively for larger simulations (in fact GMRES may stagnate without
01021         // this for larger problems), by using a AMG preconditioner -- needs HYPRE installed
01023         PetscOptionsSetValue("-pc_hypre_type", "boomeramg");
01024         // PetscOptionsSetValue("-pc_hypre_boomeramg_max_iter", "1");
01025         // PetscOptionsSetValue("-pc_hypre_boomeramg_strong_threshold", "0.0");
01026     
01027         PCSetType(pc, PCHYPRE);
01028         KSPSetPreconditionerSide(solver, PC_RIGHT);
01029         
01030         // other possible preconditioners..
01031         //PCBlockDiagonalMechanics* p_custom_pc = new PCBlockDiagonalMechanics(solver, mPreconditionMatrix, mBlock1Size, mBlock2Size);
01032         //PCLDUFactorisationMechanics* p_custom_pc = new PCLDUFactorisationMechanics(solver, mPreconditionMatrix, mBlock1Size, mBlock2Size);
01033         //remember to delete memory..
01034         //KSPSetPreconditionerSide(solver, PC_RIGHT);
01035     #endif
01036 
01037     #ifdef MECH_KSP_MONITOR
01038     PetscOptionsSetValue("-ksp_monitor","");
01039     #endif
01040 
01041     KSPSetFromOptions(solver);
01042     KSPSetUp(solver);
01043 
01044 
01045     #ifdef MECH_VERBOSE
01046     Timer::PrintAndReset("KSP Setup");
01047     #endif
01048     
01049 //todo
01050     KSPSolve(solver,mLinearSystemRhsVector,solution);
01051     
01052     int num_iters;
01053     KSPGetIterationNumber(solver, &num_iters);
01054     if(num_iters==0)
01055     {
01056         VecDestroy(solution);
01057         KSPDestroy(solver);
01058         EXCEPTION("KSP Absolute tolerance was too high, linear system wasn't solved - there will be no decrease in Newton residual. Decrease KspAbsoluteTolerance");
01059     }
01060 
01061     // See comment on max_iters above
01062     if(num_iters==1000)
01063     {
01064         #define COVERAGE_IGNORE
01065         WARNING("Linear solver in mechanics solve may not have converged");
01066         #undef COVERAGE_IGNORE
01067     }
01068 
01069     
01070     #ifdef MECH_VERBOSE
01071     Timer::PrintAndReset("KSP Solve");
01072     std::cout << "[" << PetscTools::GetMyRank() << "]: Num iterations = " << num_iters << "\n" << std::flush;
01073     #endif
01074 
01075     MechanicsEventHandler::EndEvent(MechanicsEventHandler::SOLVE);
01076 
01077 
01079     // Update the solution
01080     //  Newton method:       sol = sol - update, where update=Jac^{-1}*residual
01081     //  Newton with damping: sol = sol - s*update, where s is chosen
01082     //   such that |residual(sol)| is minimised. Damping is important to
01083     //   avoid initial divergence.
01084     //
01085     // Normally, finding the best s from say 0.05,0.1,0.2,..,1.0 is cheap,
01086     // but this is not the case in cardiac electromechanics calculations.
01087     // Therefore, we initially check s=1 (expected to be the best most of the
01088     // time, then s=0.9. If the norm of the residual increases, we assume
01089     // s=1 is the best. Otherwise, check s=0.8 to see if s=0.9 is a local min.
01091     MechanicsEventHandler::BeginEvent(MechanicsEventHandler::UPDATE);
01092     double new_norm_resid = UpdateSolutionUsingLineSearch(solution);
01093     MechanicsEventHandler::EndEvent(MechanicsEventHandler::UPDATE);
01094 
01095     VecDestroy(solution);
01096     KSPDestroy(solver);
01097 
01098     return new_norm_resid;
01099 }
01100 
01101 template<unsigned DIM>
01102 void AbstractNonlinearElasticitySolver<DIM>::PrintLineSearchResult(double s, double residNorm)
01103 {
01104     #ifdef MECH_VERBOSE
01105     std::cout << "\tTesting s = " << s << ", |f| = " << residNorm << "\n" << std::flush;
01106     #endif
01107 }
01108 
01109 template<unsigned DIM>
01110 double AbstractNonlinearElasticitySolver<DIM>::UpdateSolutionUsingLineSearch(Vec solution)
01111 {
01112     double initial_norm_resid = CalculateResidualNorm();
01113     #ifdef MECH_VERBOSE
01114     std::cout << "\tInitial |f| [corresponding to s=0] is " << initial_norm_resid << "\n"  << std::flush;
01115     #endif
01116 
01117     ReplicatableVector update(solution);
01118 
01119     std::vector<double> old_solution = mCurrentSolution;
01120 
01121     std::vector<double> damping_values; // = {1.0, 0.9, .., 0.2, 0.1, 0.05} ie size 11
01122     for (unsigned i=10; i>=1; i--)
01123     {
01124         damping_values.push_back((double)i/10.0);
01125     }
01126     damping_values.push_back(0.05);
01127     assert(damping_values.size()==11);
01128 
01129 
01131     // let mCurrentSolution = old_solution - damping_val[0]*update; and compute residual
01132     unsigned index = 0;
01133     VectorSum(old_solution, update, -damping_values[index], mCurrentSolution);
01134     double current_resid_norm = ComputeResidualAndGetNorm(true);
01135     PrintLineSearchResult(damping_values[index], current_resid_norm);
01136 
01138     // let mCurrentSolution = old_solution - damping_val[1]*update; and compute residual
01139     index = 1;
01140     VectorSum(old_solution, update, -damping_values[index], mCurrentSolution);
01141     double next_resid_norm = ComputeResidualAndGetNorm(true);
01142     PrintLineSearchResult(damping_values[index], next_resid_norm);
01143 
01144     index = 2;
01145     // While f(s_next) < f(s_current), [f = residnorm], keep trying new damping values,
01146     // ie exit thus loop when next norm of the residual first increases
01147     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
01148             || ( (next_resid_norm < current_resid_norm) && (index<damping_values.size()) ) )
01149     {
01150         current_resid_norm = next_resid_norm;
01151 
01152         // let mCurrentSolution = old_solution - damping_val*update; and compute residual
01153         VectorSum(old_solution, update, -damping_values[index], mCurrentSolution);
01154         next_resid_norm = ComputeResidualAndGetNorm(true);
01155         PrintLineSearchResult(damping_values[index], next_resid_norm);
01156 
01157         index++;
01158     }
01159 
01160     unsigned best_index;
01161 
01162     if(index==damping_values.size() && (next_resid_norm < current_resid_norm))
01163     {
01164         // Difficult to come up with large forces/tractions such that it had to
01165         // test right down to s=0.05, but overall doesn't fail.
01166         // The possible damping values have been manually temporarily altered to
01167         // get this code to be called, it appears to work correctly. Even if it
01168         // didn't tests wouldn't fail, they would just be v. slightly less efficient.
01169         #define COVERAGE_IGNORE
01170         // if we exited because we got to the end of the possible damping values, the
01171         // best one was the last one (excl the final index++ at the end)
01172         current_resid_norm = next_resid_norm;
01173         best_index = index-1;
01174         #undef COVERAGE_IGNORE
01175     }
01176     else
01177     {
01178         // else the best one must have been the second last one (excl the final index++ at the end)
01179         // (as we would have exited when the resid norm first increased)
01180         best_index = index-2;
01181     }
01182 
01183     // check out best was better than the original residual-norm
01184     if (initial_norm_resid < current_resid_norm)
01185     {
01186         #define COVERAGE_IGNORE
01187         // Have to use an assert/exit here as the following exception causes a seg fault (in cardiac mech problems?)
01188         // Don't know why
01189         std::cout << "CHASTE ERROR: (AbstractNonlinearElasticitySolver.hpp): Residual does not appear to decrease in newton direction, quitting.\n" << std::flush;
01190         exit(0);
01191         //EXCEPTION("Residual does not appear to decrease in newton direction, quitting");
01192         #undef COVERAGE_IGNORE
01193     }
01194 
01195     #ifdef MECH_VERBOSE
01196     std::cout << "\tBest s = " << damping_values[best_index] << "\n"  << std::flush;
01197     #endif
01198     VectorSum(old_solution, update, -damping_values[best_index], mCurrentSolution);
01199 
01200     return current_resid_norm;
01201 
01203 //        double best_norm_resid = DBL_MAX;
01204 //        double best_damping_value = 0.0;
01205 //
01206 //        std::vector<double> damping_values;
01207 //        damping_values.reserve(12);
01208 //        damping_values.push_back(0.0);
01209 //        damping_values.push_back(0.05);
01210 //        for (unsigned i=1; i<=10; i++)
01211 //        {
01212 //            damping_values.push_back((double)i/10.0);
01213 //        }
01214 //
01215 //        for (unsigned i=0; i<damping_values.size(); i++)
01216 //        {
01217 //            for (unsigned j=0; j<mNumDofs; j++)
01218 //            {
01219 //                mCurrentSolution[j] = old_solution[j] - damping_values[i]*update[j];
01220 //            }
01221 //
01222 //            // compute residual
01223 //            double norm_resid = ComputeResidualAndGetNorm();
01224 //
01225 //            std::cout << "\tTesting s = " << damping_values[i] << ", |f| = " << norm_resid << "\n" << std::flush;
01226 //            if (norm_resid < best_norm_resid)
01227 //            {
01228 //                best_norm_resid = norm_resid;
01229 //                best_damping_value = damping_values[i];
01230 //            }
01231 //        }
01232 //
01233 //        if (best_damping_value == 0.0)
01234 //        {
01235 //            #define COVERAGE_IGNORE
01236 //            assert(0);
01237 //            EXCEPTION("Residual does not decrease in newton direction, quitting");
01238 //            #undef COVERAGE_IGNORE
01239 //        }
01240 //        else
01241 //        {
01242 //            std::cout << "\tBest s = " << best_damping_value << "\n"  << std::flush;
01243 //        }
01244 //        //Timer::PrintAndReset("Find best damping");
01245 //
01246 //        // implement best update and recalculate residual
01247 //        for (unsigned j=0; j<mNumDofs; j++)
01248 //        {
01249 //            mCurrentSolution[j] = old_solution[j] - best_damping_value*update[j];
01250 //        }
01251 }
01252 
01253 
01254 
01255 template<unsigned DIM>
01256 void AbstractNonlinearElasticitySolver<DIM>::PostNewtonStep(unsigned counter, double normResidual)
01257 {
01258 }
01259 
01260 
01261 template<unsigned DIM>
01262 AbstractNonlinearElasticitySolver<DIM>::AbstractNonlinearElasticitySolver(QuadraticMesh<DIM>* pQuadMesh,
01263                                                                           c_vector<double,DIM> bodyForce,
01264                                                                           double density,
01265                                                                           std::string outputDirectory,
01266                                                                           std::vector<unsigned>& fixedNodes,
01267                                                                           CompressibilityType compressibilityType)
01268     : mpQuadMesh(pQuadMesh),
01269       mpQuadratureRule(NULL),
01270       mpBoundaryQuadratureRule(NULL),
01271       mKspAbsoluteTol(-1),
01272       mResidualVector(NULL),
01273       mJacobianMatrix(NULL),
01274       mPreconditionMatrix(NULL),
01275       mBodyForce(bodyForce),
01276       mDensity(density),
01277       mFixedNodes(fixedNodes),
01278       mOutputDirectory(outputDirectory),
01279       mpOutputFileHandler(NULL),
01280       mWriteOutputEachNewtonIteration(false),
01281       mNumNewtonIterations(0),
01282       mUsingBodyForceFunction(false),
01283       mUsingTractionBoundaryConditionFunction(false),
01284       mCurrentTime(0.0),
01285       mCompressibilityType(compressibilityType)
01286 {
01287     assert(DIM==2 || DIM==3);
01288     assert(density > 0);
01289     assert(fixedNodes.size() > 0);
01290     assert(pQuadMesh != NULL);
01291 
01292     if(mCompressibilityType==COMPRESSIBLE)
01293     {
01294         mNumDofs = DIM*mpQuadMesh->GetNumNodes();
01295     }
01296     else
01297     {
01298         mNumDofs = DIM*mpQuadMesh->GetNumNodes() + mpQuadMesh->GetNumVertices();
01299     }
01300 
01301     mWriteOutput = (mOutputDirectory != "");
01302     if(mWriteOutput)
01303     {
01304         mpOutputFileHandler = new OutputFileHandler(mOutputDirectory);
01305     }
01306 }
01307 
01308 
01309 
01310 template<unsigned DIM>
01311 AbstractNonlinearElasticitySolver<DIM>::~AbstractNonlinearElasticitySolver()
01312 {
01313     if(mResidualVector)
01314     {
01315         VecDestroy(mResidualVector);
01316         VecDestroy(mLinearSystemRhsVector);
01317         MatDestroy(mJacobianMatrix);
01318         MatDestroy(mPreconditionMatrix);
01319         if(mCompressibilityType==COMPRESSIBLE)
01320         {
01321             VecDestroy(mDirichletBoundaryConditionsVector);
01322         }
01323     }
01324 
01325     if(mpQuadratureRule)
01326     {
01327         delete mpQuadratureRule;
01328         delete mpBoundaryQuadratureRule;
01329     }
01330     if(mpOutputFileHandler)
01331     {
01332         delete mpOutputFileHandler;
01333     }
01334 }
01335 
01336 
01337 template<unsigned DIM>
01338 void AbstractNonlinearElasticitySolver<DIM>::Solve(double tol,
01339                                                    unsigned maxNumNewtonIterations,
01340                                                    bool quitIfNoConvergence)
01341 {
01342     WriteCurrentDeformation("initial");
01343 
01344 
01345     if(mWriteOutputEachNewtonIteration)
01346     {
01347         WriteCurrentDeformation("newton_iteration", 0);
01348     }
01349 
01350     // compute residual
01351     double norm_resid = ComputeResidualAndGetNorm(false);
01352     #ifdef MECH_VERBOSE
01353     std::cout << "\nNorm of residual is " << norm_resid << "\n";
01354     #endif
01355 
01356     mNumNewtonIterations = 0;
01357     unsigned iteration_number = 1;
01358 
01359     if (tol < 0) // ie if wasn't passed in as a parameter
01360     {
01361         tol = NEWTON_REL_TOL*norm_resid;
01362 
01363         #define COVERAGE_IGNORE // not going to have tests in cts for everything
01364         if (tol > MAX_NEWTON_ABS_TOL)
01365         {
01366             tol = MAX_NEWTON_ABS_TOL;
01367         }
01368         if (tol < MIN_NEWTON_ABS_TOL)
01369         {
01370             tol = MIN_NEWTON_ABS_TOL;
01371         }
01372         #undef COVERAGE_IGNORE
01373     }
01374 
01375     #ifdef MECH_VERBOSE
01376     std::cout << "Solving with tolerance " << tol << "\n";
01377     #endif
01378 
01379     while (norm_resid > tol && iteration_number<=maxNumNewtonIterations)
01380     {
01381         #ifdef MECH_VERBOSE
01382         std::cout <<  "\n-------------------\n"
01383                   <<   "Newton iteration " << iteration_number
01384                   << ":\n-------------------\n";
01385         #endif
01386 
01387         // take newton step (and get returned residual)
01388         norm_resid = TakeNewtonStep();
01389 
01390         #ifdef MECH_VERBOSE
01391         std::cout << "Norm of residual is " << norm_resid << "\n";
01392         #endif
01393         if (mWriteOutputEachNewtonIteration)
01394         {
01395             WriteCurrentDeformation("newton_iteration", iteration_number);
01396         }
01397 
01398         mNumNewtonIterations = iteration_number;
01399 
01400         PostNewtonStep(iteration_number,norm_resid);
01401 
01402         iteration_number++;
01403         if (iteration_number==20)
01404         {
01405             #define COVERAGE_IGNORE
01406             EXCEPTION("Not converged after 20 newton iterations, quitting");
01407             #undef COVERAGE_IGNORE
01408         }
01409     }
01410 
01411     if ((norm_resid > tol) && quitIfNoConvergence)
01412     {
01413         #define COVERAGE_IGNORE
01414         EXCEPTION("Failed to converge");
01415         #undef COVERAGE_IGNORE
01416     }
01417 
01418     // write the final solution
01419     WriteCurrentDeformation("solution");
01420 }
01421 
01422 
01423 template<unsigned DIM>
01424 void AbstractNonlinearElasticitySolver<DIM>::WriteCurrentDeformation(std::string fileName, int counterToAppend)
01425 {
01426     // only write output if the flag mWriteOutput has been set
01427     if (!mWriteOutput)
01428     {
01429         return;
01430     }
01431 
01432     std::stringstream file_name;
01433     file_name << fileName;
01434     if(counterToAppend >= 0)
01435     {
01436         file_name << "_" << counterToAppend;
01437     }
01438     file_name << ".nodes";
01439 
01440     out_stream p_file = mpOutputFileHandler->OpenOutputFile(file_name.str());
01441 
01442     std::vector<c_vector<double,DIM> >& r_deformed_position = rGetDeformedPosition();
01443     for (unsigned i=0; i<r_deformed_position.size(); i++)
01444     {
01445         for (unsigned j=0; j<DIM; j++)
01446         {
01447            * p_file << r_deformed_position[i](j) << " ";
01448         }
01449        * p_file << "\n";
01450     }
01451     p_file->close();
01452 }
01453 
01454 
01455 template<unsigned DIM>
01456 unsigned AbstractNonlinearElasticitySolver<DIM>::GetNumNewtonIterations()
01457 {
01458     return mNumNewtonIterations;
01459 }
01460 
01461 
01462 template<unsigned DIM>
01463 void AbstractNonlinearElasticitySolver<DIM>::SetFunctionalBodyForce(c_vector<double,DIM> (*pFunction)(c_vector<double,DIM>& X, double t))
01464 {
01465     mUsingBodyForceFunction = true;
01466     mpBodyForceFunction = pFunction;
01467 }
01468 
01469 
01470 template<unsigned DIM>
01471 void AbstractNonlinearElasticitySolver<DIM>::SetWriteOutput(bool writeOutput)
01472 {
01473     if (writeOutput && (mOutputDirectory==""))
01474     {
01475         EXCEPTION("Can't write output if no output directory was given in constructor");
01476     }
01477     mWriteOutput = writeOutput;
01478 }
01479 
01480 template<unsigned DIM>
01481 void AbstractNonlinearElasticitySolver<DIM>::SetSurfaceTractionBoundaryConditions(
01482             std::vector<BoundaryElement<DIM-1,DIM>*>& rBoundaryElements,
01483             std::vector<c_vector<double,DIM> >& rSurfaceTractions)
01484 {
01485     assert(rBoundaryElements.size()==rSurfaceTractions.size());
01486     mBoundaryElements = rBoundaryElements;
01487     mSurfaceTractions = rSurfaceTractions;
01488 }
01489 
01490 template<unsigned DIM>
01491 void AbstractNonlinearElasticitySolver<DIM>::SetFunctionalTractionBoundaryCondition(
01492             std::vector<BoundaryElement<DIM-1,DIM>*> rBoundaryElements,
01493             c_vector<double,DIM> (*pFunction)(c_vector<double,DIM>& X, double t))
01494 {
01495     mBoundaryElements = rBoundaryElements;
01496     mUsingTractionBoundaryConditionFunction = true;
01497     mpTractionBoundaryConditionFunction = pFunction;
01498 }
01499 
01500 
01501 template<unsigned DIM>
01502 std::vector<c_vector<double,DIM> >& AbstractNonlinearElasticitySolver<DIM>::rGetDeformedPosition()
01503 {
01504     mDeformedPosition.resize(mpQuadMesh->GetNumNodes(), zero_vector<double>(DIM));
01505     for (unsigned i=0; i<mpQuadMesh->GetNumNodes(); i++)
01506     {
01507         for (unsigned j=0; j<DIM; j++)
01508         {
01509             mDeformedPosition[i](j) = mpQuadMesh->GetNode(i)->rGetLocation()[j] + mCurrentSolution[DIM*i+j];
01510         }
01511     }
01512     return mDeformedPosition;
01513 }
01514 
01515 
01516 template<unsigned DIM>
01517 void AbstractNonlinearElasticitySolver<DIM>::CreateCmguiOutput()
01518 {
01519     if(mOutputDirectory=="")
01520     {
01521         EXCEPTION("No output directory was given so no output was written, cannot convert to cmgui format");
01522     }
01523 
01524     CmguiDeformedSolutionsWriter<DIM> writer(mOutputDirectory + "/cmgui",
01525                                              "solution",
01526                                              *mpQuadMesh,
01527                                              WRITE_QUADRATIC_MESH);
01528 
01529     std::vector<c_vector<double,DIM> >& r_deformed_positions = rGetDeformedPosition();
01530     writer.WriteInitialMesh(); // this writes solution_0.exnode and .exelem
01531     writer.WriteDeformationPositions(r_deformed_positions, 1); // this writes the final solution as solution_1.exnode
01532     writer.WriteCmguiScript(); // writes LoadSolutions.com
01533 }
01534 
01535 //
01536 // Constant setting definitions
01537 //
01538 template<unsigned DIM>
01539 double AbstractNonlinearElasticitySolver<DIM>::MAX_NEWTON_ABS_TOL = 1e-7;
01540 
01541 template<unsigned DIM>
01542 double AbstractNonlinearElasticitySolver<DIM>::MIN_NEWTON_ABS_TOL = 1e-10;
01543 
01544 template<unsigned DIM>
01545 double AbstractNonlinearElasticitySolver<DIM>::NEWTON_REL_TOL = 1e-4;
01546 
01547 
01548 #endif /*ABSTRACTNONLINEARELASTICITYSOLVER_HPP_*/

Generated on Tue May 31 14:31:49 2011 for Chaste by  doxygen 1.5.5