AbstractNonlinearElasticityAssembler.hpp

00001 /*
00002 
00003 Copyright (C) University of Oxford, 2005-2009
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 ABSTRACTNONLINEARELASTICITYASSEMBLER_HPP_
00030 #define ABSTRACTNONLINEARELASTICITYASSEMBLER_HPP_
00031 
00032 #include <vector>
00033 #include <cmath>
00034 #include "LinearSystem.hpp"
00035 #include "AbstractIncompressibleMaterialLaw.hpp"
00036 #include "OutputFileHandler.hpp"
00037 #include "LogFile.hpp"
00038 #include "PetscTools.hpp"
00039 #include "MechanicsEventHandler.hpp"
00040 #include "ReplicatableVector.hpp"
00041 //#include "Timer.hpp" 
00042 
00046 template<unsigned DIM>
00047 class AbstractNonlinearElasticityAssembler
00048 {
00049 protected:
00050 
00052     static const double MAX_NEWTON_ABS_TOL;
00053 
00055     static const double MIN_NEWTON_ABS_TOL;
00056 
00058     static const double NEWTON_REL_TOL;
00059 
00064     unsigned mNumDofs;
00065 
00071     std::vector<AbstractIncompressibleMaterialLaw<DIM>*> mMaterialLaws;
00072 
00078     LinearSystem *mpLinearSystem;
00079 
00098     LinearSystem *mpPreconditionMatrixLinearSystem;
00099 
00101     c_vector<double,DIM> mBodyForce;
00102 
00104     double mDensity;
00105 
00107     std::string mOutputDirectory;
00108 
00110     std::vector<unsigned> mFixedNodes;
00111 
00113     std::vector<c_vector<double,DIM> > mFixedNodeDisplacements;
00114 
00116     bool mWriteOutput;
00117 
00123     std::vector<double> mCurrentSolution;
00124 
00129     FourthOrderTensor<DIM> dTdE;
00130 
00132     unsigned mNumNewtonIterations;
00133 
00135     std::vector<c_vector<double,DIM> > mDeformedPosition;
00136 
00141     std::vector<double> mPressures;
00142 
00147     std::vector<c_vector<double,DIM> > mSurfaceTractions;
00148 
00150     c_vector<double,DIM> (*mpBodyForceFunction)(c_vector<double,DIM>&);
00151 
00156     c_vector<double,DIM> (*mpTractionBoundaryConditionFunction)(c_vector<double,DIM>&);
00157 
00159     bool mUsingBodyForceFunction;
00160 
00162     bool mUsingTractionBoundaryConditionFunction;
00163 
00168     virtual void FormInitialGuess()=0;
00169 
00180     virtual void AssembleSystem(bool assembleResidual, bool assembleJacobian)=0;
00181 
00186     virtual std::vector<c_vector<double,DIM> >& rGetDeformedPosition()=0;
00187 
00193     void ApplyBoundaryConditions(bool applyToMatrix);
00194 
00196     double CalculateResidualNorm();
00197 
00205     double TakeNewtonStep();
00206 
00214     virtual void PostNewtonStep(unsigned counter, double normResidual);
00215 
00220     void AllocateMatrixMemory();
00221 
00222 public:
00223 
00234     AbstractNonlinearElasticityAssembler(unsigned numDofs,
00235                                          AbstractIncompressibleMaterialLaw<DIM>* pMaterialLaw,
00236                                          c_vector<double,DIM> bodyForce,
00237                                          double density,
00238                                          std::string outputDirectory,
00239                                          std::vector<unsigned>& fixedNodes);
00240 
00251     AbstractNonlinearElasticityAssembler(unsigned numDofs,
00252                                          std::vector<AbstractIncompressibleMaterialLaw<DIM>*>& rMaterialLaws,
00253                                          c_vector<double,DIM> bodyForce,
00254                                          double density,
00255                                          std::string outputDirectory,
00256                                          std::vector<unsigned>& fixedNodes);
00257 
00261     virtual ~AbstractNonlinearElasticityAssembler();
00262 
00271     void Solve(double tol=-1.0,
00272                unsigned offset=0,
00273                unsigned maxNumNewtonIterations=INT_MAX,
00274                bool quitIfNoConvergence=true);
00275 
00281     void WriteOutput(unsigned counter);
00282 
00286     unsigned GetNumNewtonIterations();
00287 
00294     void SetFunctionalBodyForce(c_vector<double,DIM> (*pFunction)(c_vector<double,DIM>&));
00295 
00301     void SetWriteOutput(bool writeOutput=true);
00302 
00303 };
00304 
00305 
00307 // Implementation
00309 
00310 
00311 template<unsigned DIM>
00312 void AbstractNonlinearElasticityAssembler<DIM>::ApplyBoundaryConditions(bool applyToMatrix)
00313 {
00314     assert(mFixedNodeDisplacements.size()==mFixedNodes.size());
00315 
00316     // The boundary conditions on the NONLINEAR SYSTEM are x=boundary_values
00317     // on the boundary nodes. However:
00318     // The boundary conditions on the LINEAR SYSTEM  Ju=f, where J is the
00319     // u the negative update vector and f is the residual is
00320     // u=current_soln-boundary_values on the boundary nodes
00321     for (unsigned i=0; i<mFixedNodes.size(); i++)
00322     {
00323         unsigned node_index = mFixedNodes[i];
00324         for (unsigned j=0; j<DIM; j++)
00325         {
00326             unsigned dof_index = DIM*node_index+j;
00327             double value = mCurrentSolution[dof_index] - mFixedNodeDisplacements[i](j);
00328             if (applyToMatrix)
00329             {
00330                 mpLinearSystem->ZeroMatrixRow(dof_index);
00331                 mpLinearSystem->SetMatrixElement(dof_index,dof_index,1);
00332 
00333                 // apply same bcs to preconditioner matrix
00334                 mpPreconditionMatrixLinearSystem->ZeroMatrixRow(dof_index);
00335                 mpPreconditionMatrixLinearSystem->SetMatrixElement(dof_index,dof_index,1);
00336             }
00337             mpLinearSystem->SetRhsVectorElement(dof_index, value);
00338         }
00339     }
00340 }
00341 
00342 template<unsigned DIM>
00343 double AbstractNonlinearElasticityAssembler<DIM>::CalculateResidualNorm()
00344 {
00345     double norm;
00346     VecNorm(mpLinearSystem->rGetRhsVector(), NORM_2, &norm);
00347     return norm/mNumDofs;
00348 }
00349 
00350 template<unsigned DIM>
00351 void AbstractNonlinearElasticityAssembler<DIM>::AllocateMatrixMemory()
00352 {
00353     mpLinearSystem = new LinearSystem(mNumDofs, (MatType)MATAIJ); // default Mat tyype is MATMPIAIJ, see below
00354     mpPreconditionMatrixLinearSystem = new LinearSystem(mNumDofs, (MatType)MATAIJ); //MATAIJ is needed for precond to work
00355 
00356     // 2D: N elements around a point => 7N+3 non-zeros in that row? Assume N<=10 (structured mesh would have N_max=6) => 73.  
00357     // 3D: N elements around a point. nz < (3*10+6)N (lazy estimate). Better estimate is 23N+4?. Assume N<20 => 500ish
00358     unsigned num_non_zeros = DIM < 3 ? 75 : 500;
00359 
00361     // without leaking memory. This is the call to preallocate an MPI AIJ matrix: 
00362     // MatSeqAIJSetPreallocation(mpLinearSystem->rGetLhsMatrix(), num_non_zeros, PETSC_NULL, (PetscInt) (num_non_zeros*0.5), PETSC_NULL);
00363 
00364     MatSeqAIJSetPreallocation(mpLinearSystem->rGetLhsMatrix(),                   num_non_zeros, PETSC_NULL);
00365     MatSeqAIJSetPreallocation(mpPreconditionMatrixLinearSystem->rGetLhsMatrix(), num_non_zeros, PETSC_NULL);
00366 }
00367 
00368 template<unsigned DIM>
00369 double AbstractNonlinearElasticityAssembler<DIM>::TakeNewtonStep()
00370 {
00371     //Timer::Reset();
00372 
00374     // Assemble Jacobian (and preconditioner)
00376     MechanicsEventHandler::BeginEvent(MechanicsEventHandler::ASSEMBLE);
00377     AssembleSystem(true, true);
00378     MechanicsEventHandler::EndEvent(MechanicsEventHandler::ASSEMBLE);
00379     //Timer::PrintAndReset("AssembleSystem");
00380 
00382     // Solve the linear system using Petsc GMRES and an LU
00383     // factorisation of the preconditioner. Note we
00384     // don't call Solve on the linear_system as we want to
00385     // set Petsc options..
00387     MechanicsEventHandler::BeginEvent(MechanicsEventHandler::SOLVE);
00388 
00389 
00390     KSP solver;
00391     Vec solution;
00392     VecDuplicate(mpLinearSystem->rGetRhsVector(),&solution);
00393 
00394     Mat& r_jac = mpLinearSystem->rGetLhsMatrix();
00395     Mat& r_precond_jac = mpPreconditionMatrixLinearSystem->rGetLhsMatrix();
00396 
00397     KSPCreate(MPI_COMM_SELF,&solver);
00398 
00399     KSPSetOperators(solver, r_jac, r_precond_jac, DIFFERENT_NONZERO_PATTERN /*in precond between successive solves*/);
00400 
00401     // set max iterations
00402     KSPSetTolerances(solver, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10000); //hopefully with the preconditioner this max is way too high
00403     KSPSetType(solver,KSPGMRES);
00404     KSPGMRESSetRestart(solver,100); // gmres num restarts
00405 
00406     KSPSetFromOptions(solver);
00407     KSPSetUp(solver);
00408 
00409     PC pc;
00410     KSPGetPC(solver, &pc);
00411     PCSetType(pc, PCLU);         // Note: ILU factorisation doesn't have much effect, but LU works well.
00412 
00413     KSPSetFromOptions(solver);
00414     KSPSolve(solver,mpLinearSystem->rGetRhsVector(),solution);
00415 
00416     //Timer::PrintAndReset("KSP Solve");
00417     ReplicatableVector update(solution);
00418 
00419     MechanicsEventHandler::EndEvent(MechanicsEventHandler::SOLVE);
00420 
00422     // Update the solution
00423     //  Newton method:       sol = sol - update, where update=Jac^{-1}*residual
00424     //  Newton with damping: sol = sol - s*update, where s is chosen
00425     //   such that |residual(sol)| is minimised. Damping is important to
00426     //   avoid initial divergence.
00427     //
00428     // Normally, finding the best s from say 0.05,0.1,0.2,..,1.0 is cheap,
00429     // but this is not the case in cardiac electromechanics calculations.
00430     // Therefore, we initially check s=1 (expected to be the best most of the
00431     // time, then s=0.9. If the norm of the residual increases, we assume
00432     // s=1 is the best. Otherwise, check s=0.8 to see if s=0.9 is a local min.
00434     MechanicsEventHandler::BeginEvent(MechanicsEventHandler::UPDATE);
00435     std::vector<double> old_solution(mNumDofs);
00436     for (unsigned i=0; i<mNumDofs; i++)
00437     {
00438         old_solution[i] = mCurrentSolution[i];
00439     }
00440 
00441     std::vector<double> damping_values; // = {1.0, 0.9, .., 0.2, 0.1, 0.05} ie size 11
00442     for (unsigned i=10; i>=1; i--)
00443     {
00444         damping_values.push_back((double)i/10.0);
00445     }
00446     damping_values.push_back(0.05);
00447     assert(damping_values.size()==11);
00448 
00449     double initial_norm_resid = CalculateResidualNorm();
00450     unsigned index = 0;
00451     for (unsigned j=0; j<mNumDofs; j++)
00452     {
00453         mCurrentSolution[j] = old_solution[j] - damping_values[index]*update[j];
00454     }
00455 
00456     // compute residual
00457     AssembleSystem(true, false);
00458     double norm_resid = CalculateResidualNorm();
00459     std::cout << "\tTesting s = " << damping_values[index] << ", |f| = " << norm_resid << "\n" << std::flush;
00460 
00461     double next_norm_resid = -DBL_MAX;
00462     index = 1;
00463 
00464     // exit loop when next norm of the residual first increases
00465     while (next_norm_resid < norm_resid  && index<damping_values.size())
00466     {
00467         if (index!=1)
00468         {
00469             norm_resid = next_norm_resid;
00470         }
00471 
00472         for (unsigned j=0; j<mNumDofs; j++)
00473         {
00474            mCurrentSolution[j] = old_solution[j] - damping_values[index]*update[j];
00475         }
00476 
00477         // compute residual
00478         AssembleSystem(true, false);
00479         next_norm_resid = CalculateResidualNorm();
00480         std::cout << "\tTesting s = " << damping_values[index] << ", |f| = " << next_norm_resid << "\n" << std::flush;
00481         index++;
00482     }
00483 
00484     if (initial_norm_resid < norm_resid)
00485     {
00486         #define COVERAGE_IGNORE
00487         assert(0); // assert here as sometimes the following causes a seg fault - don't know why
00488         EXCEPTION("Residual does not appear to decrease in newton direction, quitting");
00489         #undef COVERAGE_IGNORE
00490     }
00491     else
00492     {
00493         index-=2;
00494         std::cout << "\tBest s = " << damping_values[index] << "\n"  << std::flush;
00495         for (unsigned j=0; j<mNumDofs; j++)
00496         {
00497             mCurrentSolution[j] = old_solution[j] - damping_values[index]*update[j];
00498         }
00499     }
00500 
00501 //        double best_norm_resid = DBL_MAX;
00502 //        double best_damping_value = 0.0;
00503 //
00504 //        std::vector<double> damping_values;
00505 //        damping_values.reserve(12);
00506 //        damping_values.push_back(0.0);
00507 //        damping_values.push_back(0.05);
00508 //        for (unsigned i=1; i<=10; i++)
00509 //        {
00510 //            damping_values.push_back((double)i/10.0);
00511 //        }
00512 //
00513 //        for (unsigned i=0; i<damping_values.size(); i++)
00514 //        {
00515 //            for (unsigned j=0; j<mNumDofs; j++)
00516 //            {
00517 //                mCurrentSolution[j] = old_solution[j] - damping_values[i]*update[j];
00518 //            }
00519 //
00520 //            // compute residual
00521 //            AssembleSystem(true, false);
00522 //            double norm_resid = CalculateResidualNorm();
00523 //
00524 //            std::cout << "\tTesting s = " << damping_values[i] << ", |f| = " << norm_resid << "\n" << std::flush;
00525 //            if (norm_resid < best_norm_resid)
00526 //            {
00527 //                best_norm_resid = norm_resid;
00528 //                best_damping_value = damping_values[i];
00529 //            }
00530 //        }
00531 //
00532 //        if (best_damping_value == 0.0)
00533 //        {
00534 //            #define COVERAGE_IGNORE
00535 //            assert(0);
00536 //            EXCEPTION("Residual does not decrease in newton direction, quitting");
00537 //            #undef COVERAGE_IGNORE
00538 //        }
00539 //        else
00540 //        {
00541 //            std::cout << "\tBest s = " << best_damping_value << "\n"  << std::flush;
00542 //        }
00543 //        //Timer::PrintAndReset("Find best damping");
00544 //
00545 //        // implement best update and recalculate residual
00546 //        for (unsigned j=0; j<mNumDofs; j++)
00547 //        {
00548 //            mCurrentSolution[j] = old_solution[j] - best_damping_value*update[j];
00549 //        }
00550     MechanicsEventHandler::EndEvent(MechanicsEventHandler::UPDATE);
00551 
00552     VecDestroy(solution);
00553     KSPDestroy(solver);
00554 
00555     return norm_resid;
00556 }
00557 
00558 
00559 template<unsigned DIM>
00560 void AbstractNonlinearElasticityAssembler<DIM>::PostNewtonStep(unsigned counter, double normResidual)
00561 {
00562 }
00563 
00564 
00565 template<unsigned DIM>
00566 AbstractNonlinearElasticityAssembler<DIM>::AbstractNonlinearElasticityAssembler(unsigned numDofs,
00567                                          AbstractIncompressibleMaterialLaw<DIM>* pMaterialLaw,
00568                                          c_vector<double,DIM> bodyForce,
00569                                          double density,
00570                                          std::string outputDirectory,
00571                                          std::vector<unsigned>& fixedNodes)
00572     : mNumDofs(numDofs),
00573       mBodyForce(bodyForce),
00574       mDensity(density),
00575       mOutputDirectory(outputDirectory),
00576       mFixedNodes(fixedNodes),
00577       mNumNewtonIterations(0),
00578       mUsingBodyForceFunction(false),
00579       mUsingTractionBoundaryConditionFunction(false)
00580 {
00581     assert(pMaterialLaw != NULL);
00582     mMaterialLaws.push_back(pMaterialLaw);
00583 
00584     assert(DIM==2 || DIM==3);
00585     assert(density > 0);
00586     assert(fixedNodes.size() > 0);
00587     mWriteOutput = (mOutputDirectory != "");
00588     
00589     AllocateMatrixMemory();
00590 }
00591 
00592 
00593 template<unsigned DIM>
00594 AbstractNonlinearElasticityAssembler<DIM>::AbstractNonlinearElasticityAssembler(unsigned numDofs,
00595                                          std::vector<AbstractIncompressibleMaterialLaw<DIM>*>& rMaterialLaws,
00596                                          c_vector<double,DIM> bodyForce,
00597                                          double density,
00598                                          std::string outputDirectory,
00599                                          std::vector<unsigned>& fixedNodes)
00600     : mNumDofs(numDofs),
00601       mBodyForce(bodyForce),
00602       mDensity(density),
00603       mOutputDirectory(outputDirectory),
00604       mFixedNodes(fixedNodes),
00605       mUsingBodyForceFunction(false),
00606       mUsingTractionBoundaryConditionFunction(false)
00607 {
00608     mMaterialLaws.resize(rMaterialLaws.size(), NULL);
00609     for (unsigned i=0; i<mMaterialLaws.size(); i++)
00610     {
00611         assert(rMaterialLaws[i] != NULL);
00612         mMaterialLaws[i] = rMaterialLaws[i];
00613     }
00614 
00615     assert(DIM==2 || DIM==3);
00616     assert(density > 0);
00617     assert(fixedNodes.size() > 0);
00618     mWriteOutput = (mOutputDirectory != "");
00619 
00620     AllocateMatrixMemory();
00621 }
00622 
00623 
00624 template<unsigned DIM>
00625 AbstractNonlinearElasticityAssembler<DIM>::~AbstractNonlinearElasticityAssembler()
00626 {
00627     delete mpLinearSystem;
00628     delete mpPreconditionMatrixLinearSystem;
00629 }
00630 
00631 
00632 template<unsigned DIM>
00633 void AbstractNonlinearElasticityAssembler<DIM>::Solve(double tol,
00634            unsigned offset,
00635            unsigned maxNumNewtonIterations,
00636            bool quitIfNoConvergence)
00637 {
00638     if (mWriteOutput)
00639     {
00640         WriteOutput(0+offset);
00641     }
00642 
00643     // compute residual
00644     AssembleSystem(true, false);
00645     double norm_resid = this->CalculateResidualNorm();
00646     std::cout << "\nNorm of residual is " << norm_resid << "\n";
00647 
00648     mNumNewtonIterations = 0;
00649     unsigned counter = 1;
00650 
00651     if (tol < 0) // ie if wasn't passed in as a parameter
00652     {
00653         tol = NEWTON_REL_TOL*norm_resid;
00654 
00655         #define COVERAGE_IGNORE // not going to have tests in cts for everything
00656         if (tol > MAX_NEWTON_ABS_TOL)
00657         {
00658             tol = MAX_NEWTON_ABS_TOL;
00659         }
00660         if (tol < MIN_NEWTON_ABS_TOL)
00661         {
00662             tol = MIN_NEWTON_ABS_TOL;
00663         }
00664         #undef COVERAGE_IGNORE
00665     }
00666 
00667     std::cout << "Solving with tolerance " << tol << "\n";
00668 
00669     while (norm_resid > tol && counter<=maxNumNewtonIterations)
00670     {
00671         std::cout <<  "\n-------------------\n"
00672                   <<   "Newton iteration " << counter
00673                   << ":\n-------------------\n";
00674 
00675         // take newton step (and get returned residual)
00676         norm_resid = TakeNewtonStep();
00677 
00678         std::cout << "Norm of residual is " << norm_resid << "\n";
00679         if (mWriteOutput)
00680         {
00681             WriteOutput(counter+offset);
00682         }
00683 
00684         mNumNewtonIterations = counter;
00685 
00686         PostNewtonStep(counter,norm_resid);
00687 
00688         counter++;
00689         if (counter==20)
00690         {
00691             #define COVERAGE_IGNORE
00692             EXCEPTION("Not converged after 20 newton iterations, quitting");
00693             #undef COVERAGE_IGNORE
00694         }
00695     }
00696 
00697     if ((norm_resid > tol) && quitIfNoConvergence)
00698     {
00699         #define COVERAGE_IGNORE
00700         EXCEPTION("Failed to converge");
00701         #undef COVERAGE_IGNORE
00702     }
00703 
00704     // we have solved for a deformation so note this
00705     //mADeformedHasBeenSolved = true;
00706 }
00707 
00708 
00709 template<unsigned DIM>
00710 void AbstractNonlinearElasticityAssembler<DIM>::WriteOutput(unsigned counter)
00711 {
00712     // only write output if the flag mWriteOutput has been set
00713     if (!mWriteOutput)
00714     {
00715         return;
00716     }
00717 
00718     OutputFileHandler output_file_handler(mOutputDirectory, (counter==0));
00719     std::stringstream file_name;
00720     file_name << "/solution_" << counter << ".nodes";
00721 
00722     out_stream p_file = output_file_handler.OpenOutputFile(file_name.str());
00723 
00724     std::vector<c_vector<double,DIM> >& r_deformed_position = rGetDeformedPosition();
00725     for (unsigned i=0; i<r_deformed_position.size(); i++)
00726     {
00727         for (unsigned j=0; j<DIM; j++)
00728         {
00729             *p_file << r_deformed_position[i](j) << " ";
00730         }
00731         *p_file << "\n";
00732     }
00733     p_file->close();
00734 }
00735 
00736 
00737 template<unsigned DIM>
00738 unsigned AbstractNonlinearElasticityAssembler<DIM>::GetNumNewtonIterations()
00739 {
00740     return mNumNewtonIterations;
00741 }
00742 
00743 
00744 template<unsigned DIM>
00745 void AbstractNonlinearElasticityAssembler<DIM>::SetFunctionalBodyForce(c_vector<double,DIM> (*pFunction)(c_vector<double,DIM>&))
00746 {
00747     mUsingBodyForceFunction = true;
00748     mpBodyForceFunction = pFunction;
00749 }
00750 
00751 
00752 template<unsigned DIM>
00753 void AbstractNonlinearElasticityAssembler<DIM>::SetWriteOutput(bool writeOutput)
00754 {
00755     if (writeOutput && (mOutputDirectory==""))
00756     {
00757         EXCEPTION("Can't write output if no output directory was given in constructor");
00758     }
00759     mWriteOutput = writeOutput;
00760 }
00761 
00762 //
00763 // Constant setting definitions
00764 //
00765 template<unsigned DIM>
00766 const double AbstractNonlinearElasticityAssembler<DIM>::MAX_NEWTON_ABS_TOL = 1e-8;
00767 
00768 template<unsigned DIM>
00769 const double AbstractNonlinearElasticityAssembler<DIM>::MIN_NEWTON_ABS_TOL = 1e-12;
00770 
00771 template<unsigned DIM>
00772 const double AbstractNonlinearElasticityAssembler<DIM>::NEWTON_REL_TOL = 1e-4;
00773 
00774 
00775 #endif /*ABSTRACTNONLINEARELASTICITYASSEMBLER_HPP_*/

Generated on Tue Aug 4 16:10:24 2009 for Chaste by  doxygen 1.5.5