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