SimpleNewtonNonlinearSolver.cpp

00001 /*
00002 
00003 Copyright (c) 2005-2015, University of Oxford.
00004 All rights reserved.
00005 
00006 University of Oxford means the Chancellor, Masters and Scholars of the
00007 University of Oxford, having an administrative office at Wellington
00008 Square, Oxford OX1 2JD, UK.
00009 
00010 This file is part of Chaste.
00011 
00012 Redistribution and use in source and binary forms, with or without
00013 modification, are permitted provided that the following conditions are met:
00014  * Redistributions of source code must retain the above copyright notice,
00015    this list of conditions and the following disclaimer.
00016  * Redistributions in binary form must reproduce the above copyright notice,
00017    this list of conditions and the following disclaimer in the documentation
00018    and/or other materials provided with the distribution.
00019  * Neither the name of the University of Oxford nor the names of its
00020    contributors may be used to endorse or promote products derived from this
00021    software without specific prior written permission.
00022 
00023 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
00024 AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00025 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
00026 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
00027 LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
00028 CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
00029 GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
00030 HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
00031 LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
00032 OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
00033 
00034 */
00035 
00036 #include "SimpleNewtonNonlinearSolver.hpp"
00037 #include <iostream>
00038 #include <cassert>
00039 #include "Exception.hpp"
00040 
00041 SimpleNewtonNonlinearSolver::SimpleNewtonNonlinearSolver()
00042         : mTolerance(1e-5),
00043         mWriteStats(false)
00044 {
00045     mTestDampingValues.push_back(-0.1);
00046     mTestDampingValues.push_back(0.05);
00047     for (unsigned i=1; i<=12; i++)
00048     {
00049         double val = double(i)/10;
00050         mTestDampingValues.push_back(val);
00051     }
00052 }
00053 
00054 SimpleNewtonNonlinearSolver::~SimpleNewtonNonlinearSolver()
00055 {}
00056 
00057 Vec SimpleNewtonNonlinearSolver::Solve(PetscErrorCode (*pComputeResidual)(SNES,Vec,Vec,void*),
00058 #if ( PETSC_VERSION_MAJOR==3 && PETSC_VERSION_MINOR>=5 )
00059                                        PetscErrorCode (*pComputeJacobian)(SNES,Vec,Mat,Mat,void*),
00060 #else
00061                                        PetscErrorCode (*pComputeJacobian)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),
00062 #endif
00063                                        Vec initialGuess,
00064                                        unsigned fill,
00065                                        void* pContext)
00066 {
00067     PetscInt size;
00068     VecGetSize(initialGuess, &size);
00069 
00070     Vec current_solution;
00071     VecDuplicate(initialGuess, &current_solution);
00072     VecCopy(initialGuess, current_solution);
00073 
00074     // The "false" says that we are allowed to do new mallocs without PETSc 3.3 causing an error
00075     LinearSystem linear_system(current_solution, fill, false);
00076 
00077     (*pComputeResidual)(NULL, current_solution, linear_system.rGetRhsVector(), pContext);
00078 
00079 
00080     double residual_norm;
00081     VecNorm(linear_system.rGetRhsVector(), NORM_2, &residual_norm);
00082     double scaled_residual_norm = residual_norm/size;
00083 
00084     if (mWriteStats)
00085     {
00086         std::cout << "Newton's method:\n  Initial ||residual||/N = " << scaled_residual_norm
00087                   << "\n  Attempting to solve to tolerance " << mTolerance << "..\n";
00088     }
00089 
00090     double old_scaled_residual_norm;
00091     unsigned counter = 0;
00092     while (scaled_residual_norm > mTolerance)
00093     {
00094         counter++;
00095 
00096         // Store the old norm to check with the new later
00097         old_scaled_residual_norm = scaled_residual_norm;
00098 
00099         // Compute Jacobian and solve J dx = f for the (negative) update dx, (J the jacobian, f the residual)
00100 #if ( PETSC_VERSION_MAJOR==3 && PETSC_VERSION_MINOR>=5 )
00101         (*pComputeJacobian)(NULL, current_solution, (linear_system.rGetLhsMatrix()), NULL, pContext);
00102 #else
00103         (*pComputeJacobian)(NULL, current_solution, &(linear_system.rGetLhsMatrix()), NULL, NULL, pContext);
00104 #endif
00105 
00106         Vec negative_update = linear_system.Solve();
00107 
00108 
00109         Vec test_vec;
00110         VecDuplicate(initialGuess, &test_vec);
00111 
00112         double best_damping_factor = 1.0;
00113         double best_scaled_residual = 1e20; // large
00114 
00115         // Loop over all the possible damping value and determine which gives smallest residual
00116         for (unsigned i=0; i<mTestDampingValues.size(); i++)
00117         {
00118             // Note: WAXPY calls VecWAXPY(w,a,x,y) which computes w = ax+y
00119             PetscVecTools::WAXPY(test_vec,-mTestDampingValues[i],negative_update,current_solution);
00120 
00121             // Compute new residual
00122             linear_system.ZeroLinearSystem();
00123             (*pComputeResidual)(NULL, test_vec, linear_system.rGetRhsVector(), pContext);
00124             VecNorm(linear_system.rGetRhsVector(), NORM_2, &residual_norm);
00125             scaled_residual_norm = residual_norm/size;
00126 
00127             if (scaled_residual_norm < best_scaled_residual)
00128             {
00129                 best_scaled_residual = scaled_residual_norm;
00130                 best_damping_factor = mTestDampingValues[i];
00131             }
00132         }
00133         PetscTools::Destroy(test_vec);
00134 
00135         // Check the smallest residual was actually smaller than the previous; if not, quit
00136         if (best_scaled_residual > old_scaled_residual_norm)
00137         {
00138             // Free memory
00139             PetscTools::Destroy(current_solution);
00140             PetscTools::Destroy(negative_update);
00141 
00142             // Raise error
00143             EXCEPTION("Iteration " << counter << ", unable to find damping factor such that residual decreases in update direction");
00144         }
00145 
00146         if (mWriteStats)
00147         {
00148             std::cout << "    Best damping factor = " << best_damping_factor << "\n";
00149         }
00150 
00151         // Update solution: current_guess = current_solution - best_damping_factor*negative_update
00152         PetscVecTools::AddScaledVector(current_solution, negative_update, -best_damping_factor);
00153         scaled_residual_norm = best_scaled_residual;
00154         PetscTools::Destroy(negative_update);
00155 
00156         // Compute best residual vector again and store in linear_system for next Solve()
00157         linear_system.ZeroLinearSystem();
00158         (*pComputeResidual)(NULL, current_solution, linear_system.rGetRhsVector(), pContext);
00159 
00160         if (mWriteStats)
00161         {
00162             std::cout << "    Iteration " << counter << ": ||residual||/N = " << scaled_residual_norm << "\n";
00163         }
00164     }
00165 
00166     if (mWriteStats)
00167     {
00168         std::cout << "  ..solved!\n\n";
00169     }
00170 
00171     return current_solution;
00172 }
00173 
00174 void SimpleNewtonNonlinearSolver::SetTolerance(double tolerance)
00175 {
00176     assert(tolerance > 0);
00177     mTolerance = tolerance;
00178 }

Generated by  doxygen 1.6.2