AbstractNonlinearElasticityAssembler.hpp

00001 /*
00002 
00003 Copyright (C) University of Oxford, 2005-2010
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 
00042 //#define MECH_VERBOSE
00043 //#define MECH_VERY_VERBOSE
00044 //#define MECH_USE_HYPRE
00045 
00046 #ifdef MECH_VERBOSE
00047 #include "Timer.hpp"
00048 #endif
00049 
00053 template<unsigned DIM>
00054 class AbstractNonlinearElasticityAssembler
00055 {
00056 protected:
00057 
00059     static const double MAX_NEWTON_ABS_TOL;
00060 
00062     static const double MIN_NEWTON_ABS_TOL;
00063 
00065     static const double NEWTON_REL_TOL;
00066 
00071     unsigned mNumDofs;
00072 
00078     std::vector<AbstractIncompressibleMaterialLaw<DIM>*> mMaterialLaws;
00079 
00085     LinearSystem* mpLinearSystem;
00086 
00105     LinearSystem* mpPreconditionMatrixLinearSystem;
00106 
00108     c_vector<double,DIM> mBodyForce;
00109 
00111     double mDensity;
00112 
00114     std::string mOutputDirectory;
00115 
00117     std::vector<unsigned> mFixedNodes;
00118 
00120     std::vector<c_vector<double,DIM> > mFixedNodeDisplacements;
00121 
00123     bool mWriteOutput;
00124 
00130     std::vector<double> mCurrentSolution;
00131 
00136     FourthOrderTensor<DIM> dTdE;
00137 
00139     unsigned mNumNewtonIterations;
00140 
00142     std::vector<c_vector<double,DIM> > mDeformedPosition;
00143 
00148     std::vector<double> mPressures;
00149 
00154     std::vector<c_vector<double,DIM> > mSurfaceTractions;
00155 
00157     c_vector<double,DIM> (*mpBodyForceFunction)(c_vector<double,DIM>&);
00158 
00163     c_vector<double,DIM> (*mpTractionBoundaryConditionFunction)(c_vector<double,DIM>&);
00164 
00166     bool mUsingBodyForceFunction;
00167 
00169     bool mUsingTractionBoundaryConditionFunction;
00170 
00175     virtual void FormInitialGuess()=0;
00176 
00187     virtual void AssembleSystem(bool assembleResidual, bool assembleJacobian)=0;
00188 
00193     virtual std::vector<c_vector<double,DIM> >& rGetDeformedPosition()=0;
00194 
00200     void ApplyBoundaryConditions(bool applyToMatrix);
00201 
00206     double ComputeResidualAndGetNorm();
00207 
00209     double CalculateResidualNorm();
00210 
00218     void VectorSum(std::vector<double>& rX, ReplicatableVector& rY, double a, std::vector<double>& rZ);
00219 
00226     void PrintLineSearchResult(double s, double residNorm);
00227 
00228 
00236     double TakeNewtonStep();
00237 
00245     double UpdateSolutionUsingLineSearch(Vec solution);
00246 
00247 
00255     virtual void PostNewtonStep(unsigned counter, double normResidual);
00256 
00257 public:
00258 
00269     AbstractNonlinearElasticityAssembler(unsigned numDofs,
00270                                          AbstractIncompressibleMaterialLaw<DIM>* pMaterialLaw,
00271                                          c_vector<double,DIM> bodyForce,
00272                                          double density,
00273                                          std::string outputDirectory,
00274                                          std::vector<unsigned>& fixedNodes);
00275 
00286     AbstractNonlinearElasticityAssembler(unsigned numDofs,
00287                                          std::vector<AbstractIncompressibleMaterialLaw<DIM>*>& rMaterialLaws,
00288                                          c_vector<double,DIM> bodyForce,
00289                                          double density,
00290                                          std::string outputDirectory,
00291                                          std::vector<unsigned>& fixedNodes);
00292 
00296     virtual ~AbstractNonlinearElasticityAssembler();
00297 
00306     void Solve(double tol=-1.0,
00307                unsigned offset=0,
00308                unsigned maxNumNewtonIterations=INT_MAX,
00309                bool quitIfNoConvergence=true);
00310 
00316     void WriteOutput(unsigned counter);
00317 
00321     unsigned GetNumNewtonIterations();
00322 
00329     void SetFunctionalBodyForce(c_vector<double,DIM> (*pFunction)(c_vector<double,DIM>&));
00330 
00336     void SetWriteOutput(bool writeOutput=true);
00337 
00338 };
00339 
00340 
00342 // Implementation
00344 
00345 
00346 template<unsigned DIM>
00347 void AbstractNonlinearElasticityAssembler<DIM>::ApplyBoundaryConditions(bool applyToMatrix)
00348 {
00349     assert(mFixedNodeDisplacements.size()==mFixedNodes.size());
00350 
00351     // The boundary conditions on the NONLINEAR SYSTEM are x=boundary_values
00352     // on the boundary nodes. However:
00353     // The boundary conditions on the LINEAR SYSTEM  Ju=f, where J is the
00354     // u the negative update vector and f is the residual is
00355     // u=current_soln-boundary_values on the boundary nodes
00356 
00357     std::vector<unsigned> rows;
00358     if(applyToMatrix)
00359     {
00360         rows.resize(DIM*mFixedNodes.size());
00361     }
00362 
00363     for (unsigned i=0; i<mFixedNodes.size(); i++)
00364     {
00365         unsigned node_index = mFixedNodes[i];
00366         for (unsigned j=0; j<DIM; j++)
00367         {
00368             unsigned dof_index = DIM*node_index+j;
00369 
00370             if(applyToMatrix)
00371             {
00372                 rows[DIM*i+j] = dof_index;
00373             }
00374 
00375             double value = mCurrentSolution[dof_index] - mFixedNodeDisplacements[i](j);
00376             mpLinearSystem->SetRhsVectorElement(dof_index, value);
00377         }
00378     }
00379 
00380     if(applyToMatrix)
00381     {
00382         mpLinearSystem->ZeroMatrixRowsWithValueOnDiagonal(rows, 1.0);
00383         mpPreconditionMatrixLinearSystem->ZeroMatrixRowsWithValueOnDiagonal(rows, 1.0);
00384     }
00385 }
00386 
00387 template<unsigned DIM>
00388 double AbstractNonlinearElasticityAssembler<DIM>::ComputeResidualAndGetNorm()
00389 {
00390       AssembleSystem(true, false);
00391 
00393 //    if(!allowException /* argument */)
00394 //    {
00395 //        // assemble the residual
00396 //        AssembleSystem(true, false);
00397 //    }
00398 //    else
00399 //    {
00400 //        try
00401 //        {
00402 //            // try to assemble the residual using this current solution
00403 //            AssembleSystem(true, false);
00404 //        }
00405 //        catch(Exception& e)
00406 //        {
00407 //            // if fail (because eg ODEs fail to solve), return infinity
00408 //            return DBL_MAX;
00409 //        }
00410 //    }
00411 
00412     // return the scaled norm of the residual
00413     return CalculateResidualNorm();
00414 }
00415 
00416 template<unsigned DIM>
00417 double AbstractNonlinearElasticityAssembler<DIM>::CalculateResidualNorm()
00418 {
00419     double norm;
00420     VecNorm(mpLinearSystem->rGetRhsVector(), NORM_2, &norm);
00421     return norm/mNumDofs;
00422 }
00423 
00424 template<unsigned DIM>
00425 void AbstractNonlinearElasticityAssembler<DIM>::VectorSum(std::vector<double>& rX,
00426                                                           ReplicatableVector& rY,
00427                                                           double a,
00428                                                           std::vector<double>& rZ)
00429 {
00430     assert(rX.size()==rY.GetSize());
00431     assert(rY.GetSize()==rZ.size());
00432     for(unsigned i=0; i<rX.size(); i++)
00433     {
00434         rZ[i] = rX[i] + a*rY[i];
00435     }
00436 }
00437 
00438 
00439 template<unsigned DIM>
00440 double AbstractNonlinearElasticityAssembler<DIM>::TakeNewtonStep()
00441 {
00442     #ifdef MECH_VERBOSE
00443     Timer::Reset();
00444     #endif
00445 
00447     // Assemble Jacobian (and preconditioner)
00449     MechanicsEventHandler::BeginEvent(MechanicsEventHandler::ASSEMBLE);
00450     AssembleSystem(true, true);
00451     MechanicsEventHandler::EndEvent(MechanicsEventHandler::ASSEMBLE);
00452     #ifdef MECH_VERBOSE
00453     Timer::PrintAndReset("AssembleSystem");
00454     #endif
00455 
00457     // Solve the linear system using Petsc GMRES and an LU
00458     // factorisation of the preconditioner. Note we
00459     // don't call Solve on the linear_system as we want to
00460     // set Petsc options..
00462     MechanicsEventHandler::BeginEvent(MechanicsEventHandler::SOLVE);
00463 
00464     KSP solver;
00465     Vec solution;
00466     VecDuplicate(mpLinearSystem->rGetRhsVector(),&solution);
00467 
00468     Mat& r_jac = mpLinearSystem->rGetLhsMatrix();
00469     Mat& r_precond_jac = mpPreconditionMatrixLinearSystem->rGetLhsMatrix();
00470 
00471     KSPCreate(PETSC_COMM_WORLD,&solver);
00472 
00473     KSPSetOperators(solver, r_jac, r_precond_jac, DIFFERENT_NONZERO_PATTERN /*in precond between successive solves*/);
00474 
00475     double ksp_relative_error_tol = 1e-6;
00476     unsigned num_restarts = 100;
00477     // set max iterations
00478     KSPSetTolerances(solver, ksp_relative_error_tol, PETSC_DEFAULT, PETSC_DEFAULT, 10000); //hopefully with the preconditioner this max is way too high
00479     KSPSetType(solver,KSPGMRES);
00480     KSPGMRESSetRestart(solver,num_restarts); // gmres num restarts
00481 
00482     KSPSetFromOptions(solver);
00483     KSPSetUp(solver);
00484     #ifdef MECH_VERBOSE
00485     Timer::PrintAndReset("KSP Setup");
00486     #endif
00487 
00488     PC pc;
00489     KSPGetPC(solver, &pc);
00490 
00491     #ifndef MECH_USE_HYPRE    
00492         PCSetType(pc, PCJACOBI);
00493     #else
00495         // Speed up linear solve time massively for larger simulations (in fact GMRES may stagnate without
00496         // this for larger problems), by using a AMG preconditioner -- needs HYPRE installed 
00498         PCSetType(pc, PCHYPRE);
00499         KSPSetPreconditionerSide(solver, PC_RIGHT);
00500     #endif     
00501 
00502     KSPSetFromOptions(solver);
00503     
00504     KSPSolve(solver,mpLinearSystem->rGetRhsVector(),solution);
00505     
00506     #ifdef MECH_VERBOSE
00507     Timer::PrintAndReset("KSP Solve");
00508     int num_iters;
00509     KSPGetIterationNumber(solver, &num_iters);
00510     std::cout << "[" << PetscTools::GetMyRank() << "]: Num iterations = " << num_iters << "\n" << std::flush;
00511     #endif
00512 
00513 
00514     MechanicsEventHandler::EndEvent(MechanicsEventHandler::SOLVE);
00515 
00517     // Update the solution
00518     //  Newton method:       sol = sol - update, where update=Jac^{-1}*residual
00519     //  Newton with damping: sol = sol - s*update, where s is chosen
00520     //   such that |residual(sol)| is minimised. Damping is important to
00521     //   avoid initial divergence.
00522     //
00523     // Normally, finding the best s from say 0.05,0.1,0.2,..,1.0 is cheap,
00524     // but this is not the case in cardiac electromechanics calculations.
00525     // Therefore, we initially check s=1 (expected to be the best most of the
00526     // time, then s=0.9. If the norm of the residual increases, we assume
00527     // s=1 is the best. Otherwise, check s=0.8 to see if s=0.9 is a local min.
00529     MechanicsEventHandler::BeginEvent(MechanicsEventHandler::UPDATE);
00530     double new_norm_resid = UpdateSolutionUsingLineSearch(solution);
00531     MechanicsEventHandler::EndEvent(MechanicsEventHandler::UPDATE);
00532 
00533     VecDestroy(solution);
00534     KSPDestroy(solver);
00535 
00536     return new_norm_resid;
00537 }
00538 
00539 template<unsigned DIM>
00540 void AbstractNonlinearElasticityAssembler<DIM>::PrintLineSearchResult(double s, double residNorm)
00541 {
00542     #ifdef MECH_VERBOSE
00543     std::cout << "\tTesting s = " << s << ", |f| = " << residNorm << "\n" << std::flush;
00544     #endif
00545 }
00546 
00547 template<unsigned DIM>
00548 double AbstractNonlinearElasticityAssembler<DIM>::UpdateSolutionUsingLineSearch(Vec solution)
00549 {
00550     double initial_norm_resid = CalculateResidualNorm();
00551     #ifdef MECH_VERBOSE
00552     std::cout << "\tInitial |f| [corresponding to s=0] is " << initial_norm_resid << "\n"  << std::flush;
00553     #endif
00554 
00555     ReplicatableVector update(solution);
00556 
00557     std::vector<double> old_solution = mCurrentSolution;
00558 
00559     std::vector<double> damping_values; // = {1.0, 0.9, .., 0.2, 0.1, 0.05} ie size 11
00560     for (unsigned i=10; i>=1; i--)
00561     {
00562         damping_values.push_back((double)i/10.0);
00563     }
00564     damping_values.push_back(0.05);
00565     assert(damping_values.size()==11);
00566 
00567 
00569     // let mCurrentSolution = old_solution - damping_val[0]*update; and compute residual
00570     unsigned index = 0;
00571     VectorSum(old_solution, update, -damping_values[index], mCurrentSolution);
00572     double current_resid_norm = ComputeResidualAndGetNorm();
00573     PrintLineSearchResult(damping_values[index], current_resid_norm);
00574 
00576     // let mCurrentSolution = old_solution - damping_val[1]*update; and compute residual
00577     index = 1;
00578     VectorSum(old_solution, update, -damping_values[index], mCurrentSolution);
00579     double next_resid_norm = ComputeResidualAndGetNorm();
00580     PrintLineSearchResult(damping_values[index], next_resid_norm);
00581 
00582     index = 2;
00583     // While f(s_next) < f(s_current), [f = residnorm], keep trying new damping values,
00584     // ie exit thus loop when next norm of the residual first increases
00585     while ( (next_resid_norm < current_resid_norm)  && index<damping_values.size())
00586     {
00587         current_resid_norm = next_resid_norm;
00588 
00589         // let mCurrentSolution = old_solution - damping_val*update; and compute residual
00590         VectorSum(old_solution, update, -damping_values[index], mCurrentSolution);
00591         next_resid_norm = ComputeResidualAndGetNorm();
00592         PrintLineSearchResult(damping_values[index], next_resid_norm);
00593 
00594         index++;
00595     }
00596 
00597     unsigned best_index;
00598 
00599     if(index==damping_values.size() && (next_resid_norm < current_resid_norm))
00600     {
00601         // Difficult to come up with large forces/tractions such that it had to
00602         // test right down to s=0.05, but overall doesn't fail.
00603         // The possible damping values have been manually temporarily altered to
00604         // get this code to be called, it appears to work correctly. Even if it
00605         // didn't tests wouldn't fail, they would just be v. slightly less efficient.
00606         #define COVERAGE_IGNORE
00607         // if we exited because we got to the end of the possible damping values, the
00608         // best one was the last one (excl the final index++ at the end)
00609         current_resid_norm = next_resid_norm;
00610         best_index = index-1;
00611         #undef COVERAGE_IGNORE
00612     }
00613     else
00614     {
00615         // else the best one must have been the second last one (excl the final index++ at the end)
00616         // (as we would have exited when the resid norm first increased)
00617         best_index = index-2;
00618     }
00619 
00620     // check out best was better than the original residual-norm
00621     if (initial_norm_resid < current_resid_norm)
00622     {
00623         #define COVERAGE_IGNORE
00624         // Have to use an assert/exit here as the following exception causes a seg fault (in cardiac mech problems?)
00625         // Don't know why
00626         std::cout << "CHASTE ERROR: (AbstractNonlinearElasticityAssembler.hpp): Residual does not appear to decrease in newton direction, quitting.\n" << std::flush;
00627         exit(0);
00628         //EXCEPTION("Residual does not appear to decrease in newton direction, quitting");
00629         #undef COVERAGE_IGNORE
00630     }
00631 
00632     #ifdef MECH_VERBOSE
00633     std::cout << "\tBest s = " << damping_values[best_index] << "\n"  << std::flush;
00634     #endif
00635     VectorSum(old_solution, update, -damping_values[best_index], mCurrentSolution);
00636 
00637     return current_resid_norm;
00638 
00640 //        double best_norm_resid = DBL_MAX;
00641 //        double best_damping_value = 0.0;
00642 //
00643 //        std::vector<double> damping_values;
00644 //        damping_values.reserve(12);
00645 //        damping_values.push_back(0.0);
00646 //        damping_values.push_back(0.05);
00647 //        for (unsigned i=1; i<=10; i++)
00648 //        {
00649 //            damping_values.push_back((double)i/10.0);
00650 //        }
00651 //
00652 //        for (unsigned i=0; i<damping_values.size(); i++)
00653 //        {
00654 //            for (unsigned j=0; j<mNumDofs; j++)
00655 //            {
00656 //                mCurrentSolution[j] = old_solution[j] - damping_values[i]*update[j];
00657 //            }
00658 //
00659 //            // compute residual
00660 //            double norm_resid = ComputeResidualAndGetNorm();
00661 //
00662 //            std::cout << "\tTesting s = " << damping_values[i] << ", |f| = " << norm_resid << "\n" << std::flush;
00663 //            if (norm_resid < best_norm_resid)
00664 //            {
00665 //                best_norm_resid = norm_resid;
00666 //                best_damping_value = damping_values[i];
00667 //            }
00668 //        }
00669 //
00670 //        if (best_damping_value == 0.0)
00671 //        {
00672 //            #define COVERAGE_IGNORE
00673 //            assert(0);
00674 //            EXCEPTION("Residual does not decrease in newton direction, quitting");
00675 //            #undef COVERAGE_IGNORE
00676 //        }
00677 //        else
00678 //        {
00679 //            std::cout << "\tBest s = " << best_damping_value << "\n"  << std::flush;
00680 //        }
00681 //        //Timer::PrintAndReset("Find best damping");
00682 //
00683 //        // implement best update and recalculate residual
00684 //        for (unsigned j=0; j<mNumDofs; j++)
00685 //        {
00686 //            mCurrentSolution[j] = old_solution[j] - best_damping_value*update[j];
00687 //        }
00688 }
00689 
00690 
00691 
00692 template<unsigned DIM>
00693 void AbstractNonlinearElasticityAssembler<DIM>::PostNewtonStep(unsigned counter, double normResidual)
00694 {
00695 }
00696 
00697 
00698 template<unsigned DIM>
00699 AbstractNonlinearElasticityAssembler<DIM>::AbstractNonlinearElasticityAssembler(unsigned numDofs,
00700                                          AbstractIncompressibleMaterialLaw<DIM>* pMaterialLaw,
00701                                          c_vector<double,DIM> bodyForce,
00702                                          double density,
00703                                          std::string outputDirectory,
00704                                          std::vector<unsigned>& fixedNodes)
00705     : mNumDofs(numDofs),
00706       mBodyForce(bodyForce),
00707       mDensity(density),
00708       mOutputDirectory(outputDirectory),
00709       mFixedNodes(fixedNodes),
00710       mNumNewtonIterations(0),
00711       mUsingBodyForceFunction(false),
00712       mUsingTractionBoundaryConditionFunction(false)
00713 {
00714     assert(pMaterialLaw != NULL);
00715     mMaterialLaws.push_back(pMaterialLaw);
00716 
00717     assert(DIM==2 || DIM==3);
00718     assert(density > 0);
00719     assert(fixedNodes.size() > 0);
00720     mWriteOutput = (mOutputDirectory != "");
00721 }
00722 
00723 
00724 template<unsigned DIM>
00725 AbstractNonlinearElasticityAssembler<DIM>::AbstractNonlinearElasticityAssembler(unsigned numDofs,
00726                                          std::vector<AbstractIncompressibleMaterialLaw<DIM>*>& rMaterialLaws,
00727                                          c_vector<double,DIM> bodyForce,
00728                                          double density,
00729                                          std::string outputDirectory,
00730                                          std::vector<unsigned>& fixedNodes)
00731     : mNumDofs(numDofs),
00732       mBodyForce(bodyForce),
00733       mDensity(density),
00734       mOutputDirectory(outputDirectory),
00735       mFixedNodes(fixedNodes),
00736       mUsingBodyForceFunction(false),
00737       mUsingTractionBoundaryConditionFunction(false)
00738 {
00739     mMaterialLaws.resize(rMaterialLaws.size(), NULL);
00740     for (unsigned i=0; i<mMaterialLaws.size(); i++)
00741     {
00742         assert(rMaterialLaws[i] != NULL);
00743         mMaterialLaws[i] = rMaterialLaws[i];
00744     }
00745 
00746     assert(DIM==2 || DIM==3);
00747     assert(density > 0);
00748     assert(fixedNodes.size() > 0);
00749     mWriteOutput = (mOutputDirectory != "");
00750 }
00751 
00752 
00753 template<unsigned DIM>
00754 AbstractNonlinearElasticityAssembler<DIM>::~AbstractNonlinearElasticityAssembler()
00755 {
00756     delete mpLinearSystem;
00757     delete mpPreconditionMatrixLinearSystem;
00758 }
00759 
00760 
00761 template<unsigned DIM>
00762 void AbstractNonlinearElasticityAssembler<DIM>::Solve(double tol,
00763                                                       unsigned offset,
00764                                                       unsigned maxNumNewtonIterations,
00765                                                       bool quitIfNoConvergence)
00766 {
00767     if (mWriteOutput)
00768     {
00769         WriteOutput(0+offset);
00770     }
00771 
00772     // compute residual
00773     double norm_resid = this->ComputeResidualAndGetNorm();
00774     #ifdef MECH_VERBOSE
00775     std::cout << "\nNorm of residual is " << norm_resid << "\n";
00776     #endif
00777 
00778     mNumNewtonIterations = 0;
00779     unsigned counter = 1;
00780 
00781     if (tol < 0) // ie if wasn't passed in as a parameter
00782     {
00783         tol = NEWTON_REL_TOL*norm_resid;
00784 
00785         #define COVERAGE_IGNORE // not going to have tests in cts for everything
00786         if (tol > MAX_NEWTON_ABS_TOL)
00787         {
00788             tol = MAX_NEWTON_ABS_TOL;
00789         }
00790         if (tol < MIN_NEWTON_ABS_TOL)
00791         {
00792             tol = MIN_NEWTON_ABS_TOL;
00793         }
00794         #undef COVERAGE_IGNORE
00795     }
00796 
00797     #ifdef MECH_VERBOSE
00798     std::cout << "Solving with tolerance " << tol << "\n";
00799     #endif
00800 
00801     while (norm_resid > tol && counter<=maxNumNewtonIterations)
00802     {
00803         #ifdef MECH_VERBOSE
00804         std::cout <<  "\n-------------------\n"
00805                   <<   "Newton iteration " << counter
00806                   << ":\n-------------------\n";
00807         #endif
00808 
00809         // take newton step (and get returned residual)
00810         norm_resid = TakeNewtonStep();
00811 
00812         #ifdef MECH_VERBOSE
00813         std::cout << "Norm of residual is " << norm_resid << "\n";
00814         #endif
00815         if (mWriteOutput)
00816         {
00817             WriteOutput(counter+offset);
00818         }
00819 
00820         mNumNewtonIterations = counter;
00821 
00822         PostNewtonStep(counter,norm_resid);
00823 
00824         counter++;
00825         if (counter==20)
00826         {
00827             #define COVERAGE_IGNORE
00828             EXCEPTION("Not converged after 20 newton iterations, quitting");
00829             #undef COVERAGE_IGNORE
00830         }
00831     }
00832 
00833     if ((norm_resid > tol) && quitIfNoConvergence)
00834     {
00835         #define COVERAGE_IGNORE
00836         EXCEPTION("Failed to converge");
00837         #undef COVERAGE_IGNORE
00838     }
00839 }
00840 
00841 
00842 template<unsigned DIM>
00843 void AbstractNonlinearElasticityAssembler<DIM>::WriteOutput(unsigned counter)
00844 {
00845     // only write output if the flag mWriteOutput has been set
00846     if (!mWriteOutput)
00847     {
00848         return;
00849     }
00850 
00851     OutputFileHandler output_file_handler(mOutputDirectory, (counter==0));
00852     std::stringstream file_name;
00853     file_name << "/solution_" << counter << ".nodes";
00854 
00855     out_stream p_file = output_file_handler.OpenOutputFile(file_name.str());
00856 
00857     std::vector<c_vector<double,DIM> >& r_deformed_position = rGetDeformedPosition();
00858     for (unsigned i=0; i<r_deformed_position.size(); i++)
00859     {
00860         for (unsigned j=0; j<DIM; j++)
00861         {
00862            * p_file << r_deformed_position[i](j) << " ";
00863         }
00864        * p_file << "\n";
00865     }
00866     p_file->close();
00867 }
00868 
00869 
00870 template<unsigned DIM>
00871 unsigned AbstractNonlinearElasticityAssembler<DIM>::GetNumNewtonIterations()
00872 {
00873     return mNumNewtonIterations;
00874 }
00875 
00876 
00877 template<unsigned DIM>
00878 void AbstractNonlinearElasticityAssembler<DIM>::SetFunctionalBodyForce(c_vector<double,DIM> (*pFunction)(c_vector<double,DIM>&))
00879 {
00880     mUsingBodyForceFunction = true;
00881     mpBodyForceFunction = pFunction;
00882 }
00883 
00884 
00885 template<unsigned DIM>
00886 void AbstractNonlinearElasticityAssembler<DIM>::SetWriteOutput(bool writeOutput)
00887 {
00888     if (writeOutput && (mOutputDirectory==""))
00889     {
00890         EXCEPTION("Can't write output if no output directory was given in constructor");
00891     }
00892     mWriteOutput = writeOutput;
00893 }
00894 
00895 //
00896 // Constant setting definitions
00897 //
00898 template<unsigned DIM>
00899 const double AbstractNonlinearElasticityAssembler<DIM>::MAX_NEWTON_ABS_TOL = 1e-7;
00900 
00901 template<unsigned DIM>
00902 const double AbstractNonlinearElasticityAssembler<DIM>::MIN_NEWTON_ABS_TOL = 1e-10;
00903 
00904 template<unsigned DIM>
00905 const double AbstractNonlinearElasticityAssembler<DIM>::NEWTON_REL_TOL = 1e-4;
00906 
00907 
00908 #endif /*ABSTRACTNONLINEARELASTICITYASSEMBLER_HPP_*/

Generated by  doxygen 1.6.2