Chaste  Release::2017.1
SimpleNewtonNonlinearSolver.cpp
1 /*
2 
3 Copyright (c) 2005-2017, University of Oxford.
4 All rights reserved.
5 
6 University of Oxford means the Chancellor, Masters and Scholars of the
7 University of Oxford, having an administrative office at Wellington
8 Square, Oxford OX1 2JD, UK.
9 
10 This file is part of Chaste.
11 
12 Redistribution and use in source and binary forms, with or without
13 modification, are permitted provided that the following conditions are met:
14  * Redistributions of source code must retain the above copyright notice,
15  this list of conditions and the following disclaimer.
16  * Redistributions in binary form must reproduce the above copyright notice,
17  this list of conditions and the following disclaimer in the documentation
18  and/or other materials provided with the distribution.
19  * Neither the name of the University of Oxford nor the names of its
20  contributors may be used to endorse or promote products derived from this
21  software without specific prior written permission.
22 
23 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
24 AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
25 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
26 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
27 LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
28 CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
29 GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
30 HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
31 LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
32 OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
33 
34 */
35 
36 #include "SimpleNewtonNonlinearSolver.hpp"
37 #include <iostream>
38 #include <cassert>
39 #include "Exception.hpp"
40 
42  : mTolerance(1e-5),
43  mWriteStats(false)
44 {
45  mTestDampingValues.push_back(-0.1);
46  mTestDampingValues.push_back(0.05);
47  for (unsigned i=1; i<=12; i++)
48  {
49  double val = double(i)/10;
50  mTestDampingValues.push_back(val);
51  }
52 }
53 
55 {}
56 
57 Vec SimpleNewtonNonlinearSolver::Solve(PetscErrorCode (*pComputeResidual)(SNES,Vec,Vec,void*),
58 #if (PETSC_VERSION_MAJOR==3 && PETSC_VERSION_MINOR>=5)
59  PetscErrorCode (*pComputeJacobian)(SNES,Vec,Mat,Mat,void*),
60 #else
61  PetscErrorCode (*pComputeJacobian)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),
62 #endif
63  Vec initialGuess,
64  unsigned fill,
65  void* pContext)
66 {
67  PetscInt size;
68  VecGetSize(initialGuess, &size);
69 
70  Vec current_solution;
71  VecDuplicate(initialGuess, &current_solution);
72  VecCopy(initialGuess, current_solution);
73 
74  // The "false" says that we are allowed to do new mallocs without PETSc 3.3 causing an error
75  LinearSystem linear_system(current_solution, fill, false);
76 
77  (*pComputeResidual)(nullptr, current_solution, linear_system.rGetRhsVector(), pContext);
78 
79 
80  double residual_norm;
81  VecNorm(linear_system.rGetRhsVector(), NORM_2, &residual_norm);
82  double scaled_residual_norm = residual_norm/size;
83 
84  if (mWriteStats)
85  {
86  std::cout << "Newton's method:\n Initial ||residual||/N = " << scaled_residual_norm
87  << "\n Attempting to solve to tolerance " << mTolerance << "..\n";
88  }
89 
90  double old_scaled_residual_norm;
91  unsigned counter = 0;
92  while (scaled_residual_norm > mTolerance)
93  {
94  counter++;
95 
96  // Store the old norm to check with the new later
97  old_scaled_residual_norm = scaled_residual_norm;
98 
99  // Compute Jacobian and solve J dx = f for the (negative) update dx, (J the jacobian, f the residual)
100 #if (PETSC_VERSION_MAJOR==3 && PETSC_VERSION_MINOR>=5)
101  (*pComputeJacobian)(nullptr, current_solution, (linear_system.rGetLhsMatrix()), nullptr, pContext);
102 #else
103  (*pComputeJacobian)(nullptr, current_solution, &(linear_system.rGetLhsMatrix()), nullptr, nullptr, pContext);
104 #endif
105 
106  Vec negative_update = linear_system.Solve();
107 
108 
109  Vec test_vec;
110  VecDuplicate(initialGuess, &test_vec);
111 
112  double best_damping_factor = 1.0;
113  double best_scaled_residual = 1e20; // large
114 
115  // Loop over all the possible damping value and determine which gives smallest residual
116  for (unsigned i=0; i<mTestDampingValues.size(); i++)
117  {
118  // Note: WAXPY calls VecWAXPY(w,a,x,y) which computes w = ax+y
119  PetscVecTools::WAXPY(test_vec,-mTestDampingValues[i],negative_update,current_solution);
120 
121  // Compute new residual
122  linear_system.ZeroLinearSystem();
123  (*pComputeResidual)(nullptr, test_vec, linear_system.rGetRhsVector(), pContext);
124  VecNorm(linear_system.rGetRhsVector(), NORM_2, &residual_norm);
125  scaled_residual_norm = residual_norm/size;
126 
127  if (scaled_residual_norm < best_scaled_residual)
128  {
129  best_scaled_residual = scaled_residual_norm;
130  best_damping_factor = mTestDampingValues[i];
131  }
132  }
133  PetscTools::Destroy(test_vec);
134 
135  // Check the smallest residual was actually smaller than the previous; if not, quit
136  if (best_scaled_residual > old_scaled_residual_norm)
137  {
138  // Free memory
139  PetscTools::Destroy(current_solution);
140  PetscTools::Destroy(negative_update);
141 
142  // Raise error
143  EXCEPTION("Iteration " << counter << ", unable to find damping factor such that residual decreases in update direction");
144  }
145 
146  if (mWriteStats)
147  {
148  std::cout << " Best damping factor = " << best_damping_factor << "\n";
149  }
150 
151  // Update solution: current_guess = current_solution - best_damping_factor*negative_update
152  PetscVecTools::AddScaledVector(current_solution, negative_update, -best_damping_factor);
153  scaled_residual_norm = best_scaled_residual;
154  PetscTools::Destroy(negative_update);
155 
156  // Compute best residual vector again and store in linear_system for next Solve()
157  linear_system.ZeroLinearSystem();
158  (*pComputeResidual)(nullptr, current_solution, linear_system.rGetRhsVector(), pContext);
159 
160  if (mWriteStats)
161  {
162  std::cout << " Iteration " << counter << ": ||residual||/N = " << scaled_residual_norm << "\n";
163  }
164  }
165 
166  if (mWriteStats)
167  {
168  std::cout << " ..solved!\n\n";
169  }
170 
171  return current_solution;
172 }
173 
175 {
176  assert(tolerance > 0);
177  mTolerance = tolerance;
178 }
void ZeroLinearSystem()
Vec Solve(Vec lhsGuess=nullptr)
#define EXCEPTION(message)
Definition: Exception.hpp:143
virtual Vec Solve(PetscErrorCode(*pComputeResidual)(SNES, Vec, Vec, void *), PetscErrorCode(*pComputeJacobian)(SNES, Vec, Mat *, Mat *, MatStructure *, void *), Vec initialGuess, unsigned fill, void *pContext)
Mat & rGetLhsMatrix()
Vec & rGetRhsVector()
static void AddScaledVector(Vec y, Vec x, double scaleFactor)
static void Destroy(Vec &rVec)
Definition: PetscTools.hpp:352
static void WAXPY(Vec w, double a, Vec x, Vec y)