SimpleNewtonNonlinearSolver.cpp

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 #include "SimpleNewtonNonlinearSolver.hpp"
00030 #include <iostream>
00031 #include <cassert>
00032 #include "Exception.hpp"
00033 
00034 SimpleNewtonNonlinearSolver::SimpleNewtonNonlinearSolver(double linearSolverRelativeTolerance)
00035         : mLinearSolverRelativeTolerance(linearSolverRelativeTolerance),
00036         mTolerance(1e-5),
00037         mWriteStats(false)
00038 {
00039     mTestDampingValues.push_back(-0.1);
00040     mTestDampingValues.push_back(0.05);
00041     for (unsigned i=1; i<=12; i++)
00042     {
00043         double val = double(i)/10;
00044         mTestDampingValues.push_back(val);
00045     }
00046 }
00047 
00048 SimpleNewtonNonlinearSolver::~SimpleNewtonNonlinearSolver()
00049 {}
00050 
00051 Vec SimpleNewtonNonlinearSolver::Solve(PetscErrorCode (*pComputeResidual)(SNES,Vec,Vec,void*),
00052                                        PetscErrorCode (*pComputeJacobian)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),
00053                                        Vec initialGuess,
00054                                        unsigned fill,
00055                                        void* pContext)
00056 {
00057     PetscInt size;
00058     VecGetSize(initialGuess, &size);
00059 
00060     Vec current_solution;
00061     VecDuplicate(initialGuess, &current_solution);
00062     VecCopy(initialGuess, current_solution);
00063 
00064     LinearSystem linear_system(current_solution, fill);
00065 
00066     (*pComputeResidual)(NULL, current_solution, linear_system.rGetRhsVector(), pContext);
00067 
00068 
00069     double residual_norm;
00070     VecNorm(linear_system.rGetRhsVector(), NORM_2, &residual_norm);
00071     double scaled_residual_norm = residual_norm/size;
00072 
00073     if (mWriteStats)
00074     {
00075         std::cout << "Newton's method:\n  Initial ||residual||/N = " << scaled_residual_norm
00076                   << "\n  Attempting to solve to tolerance " << mTolerance << "..\n";
00077     }
00078 
00079     double old_scaled_residual_norm;
00080     unsigned counter = 0;
00081     while (scaled_residual_norm > mTolerance)
00082     {
00083         counter++;
00084 
00085         // Store the old norm to check with the new later
00086         old_scaled_residual_norm = scaled_residual_norm;
00087 
00088         // Compute Jacobian and solve J dx = f for the (negative) update dx, (J the jacobian, f the residual)
00089         (*pComputeJacobian)(NULL, current_solution, &(linear_system.rGetLhsMatrix()), NULL, NULL, pContext);
00090 
00091         Vec negative_update = linear_system.Solve();
00092 
00093 
00094         Vec test_vec;
00095         VecDuplicate(initialGuess, &test_vec);
00096 
00097         double best_damping_factor = 1.0;
00098         double best_scaled_residual = 1e20; // large
00099 
00100         // Loop over all the possible damping value and determine which gives smallest residual
00101         for (unsigned i=0; i<mTestDampingValues.size(); i++)
00102         {
00103             // Note: WAXPY calls VecWAXPY(w,a,x,y) which computes w = ax+y
00104             PetscVecTools::WAXPY(test_vec,-mTestDampingValues[i],negative_update,current_solution);
00105 
00106             // Compute new residual
00107             linear_system.ZeroLinearSystem();
00108             (*pComputeResidual)(NULL, test_vec, linear_system.rGetRhsVector(), pContext);
00109             VecNorm(linear_system.rGetRhsVector(), NORM_2, &residual_norm);
00110             scaled_residual_norm = residual_norm/size;
00111 
00112             if (scaled_residual_norm < best_scaled_residual)
00113             {
00114                 best_scaled_residual = scaled_residual_norm;
00115                 best_damping_factor = mTestDampingValues[i];
00116             }
00117         }
00118         VecDestroy(test_vec);
00119 
00120         // Check the smallest residual was actually smaller than the previous; if not, quit
00121         if (best_scaled_residual > old_scaled_residual_norm)
00122         {
00123             // Free memory
00124             VecDestroy(current_solution);
00125             VecDestroy(negative_update);
00126 
00127             // Raise error
00128             EXCEPTION("Iteration " << counter << ", unable to find damping factor such that residual decreases in update direction");
00129         }
00130 
00131         if (mWriteStats)
00132         {
00133             std::cout << "    Best damping factor = " << best_damping_factor << "\n";
00134         }
00135 
00136         // Update solution: current_guess = current_solution - best_damping_factor*negative_update
00137         PetscVecTools::AddScaledVector(current_solution, negative_update, -best_damping_factor);
00138         scaled_residual_norm = best_scaled_residual;
00139         VecDestroy(negative_update);
00140 
00141         // Compute best residual vector again and store in linear_system for next Solve()
00142         linear_system.ZeroLinearSystem();
00143         (*pComputeResidual)(NULL, current_solution, linear_system.rGetRhsVector(), pContext);
00144 
00145         if (mWriteStats)
00146         {
00147             std::cout << "    Iteration " << counter << ": ||residual||/N = " << scaled_residual_norm << "\n";
00148         }
00149     }
00150 
00151     if (mWriteStats)
00152     {
00153         std::cout << "  ..solved!\n\n";
00154     }
00155 
00156     return current_solution;
00157 }
00158 
00159 void SimpleNewtonNonlinearSolver::SetTolerance(double tolerance)
00160 {
00161     assert(tolerance > 0);
00162     mTolerance = tolerance;
00163 }
Generated on Thu Dec 22 13:00:07 2011 for Chaste by  doxygen 1.6.3