VTK-m  2.2
NewtonsMethod.h
Go to the documentation of this file.
1 //============================================================================
2 // Copyright (c) Kitware, Inc.
3 // All rights reserved.
4 // See LICENSE.txt for details.
5 //
6 // This software is distributed WITHOUT ANY WARRANTY; without even
7 // the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
8 // PURPOSE. See the above copyright notice for more information.
9 //============================================================================
10 #ifndef vtk_m_NewtonsMethod_h
11 #define vtk_m_NewtonsMethod_h
12 
13 #include <vtkm/Math.h>
14 #include <vtkm/Matrix.h>
15 
16 namespace vtkm
17 {
18 
21 template <typename ScalarType, vtkm::IdComponent Size>
23 {
25  bool Valid;
27  bool Converged;
32 };
33 
62 template <typename ScalarType,
63  vtkm::IdComponent Size,
64  typename JacobianFunctor,
65  typename FunctionFunctor>
67  JacobianFunctor jacobianEvaluator,
68  FunctionFunctor functionEvaluator,
69  vtkm::Vec<ScalarType, Size> desiredFunctionOutput,
70  vtkm::Vec<ScalarType, Size> initialGuess = vtkm::Vec<ScalarType, Size>(ScalarType(0)),
71  ScalarType convergeDifference = ScalarType(1e-3),
72  vtkm::IdComponent maxIterations = 10)
73 {
74  using VectorType = vtkm::Vec<ScalarType, Size>;
75  using MatrixType = vtkm::Matrix<ScalarType, Size, Size>;
76 
77  VectorType x = initialGuess;
78 
79  bool valid = false;
80  bool converged = false;
81  for (vtkm::IdComponent iteration = 0; !converged && (iteration < maxIterations); iteration++)
82  {
83  // For Newton's method, we solve the linear system
84  //
85  // Jacobian x deltaX = currentFunctionOutput - desiredFunctionOutput
86  //
87  // The subtraction on the right side simply makes the target of the solve
88  // at zero, which is what Newton's method solves for. The deltaX tells us
89  // where to move to to solve for a linear system, which we assume will be
90  // closer for our nonlinear system.
91 
92  MatrixType jacobian = jacobianEvaluator(x);
93  VectorType currentFunctionOutput = functionEvaluator(x);
94 
95  VectorType deltaX =
96  vtkm::SolveLinearSystem(jacobian, currentFunctionOutput - desiredFunctionOutput, valid);
97  if (!valid)
98  {
99  break;
100  }
101 
102  x = x - deltaX;
103 
104  converged = true;
105  for (vtkm::IdComponent index = 0; index < Size; index++)
106  {
107  converged &= (vtkm::Abs(deltaX[index]) < convergeDifference);
108  }
109  }
110 
111  // Not checking whether converged.
112  return { valid, converged, x };
113 }
114 
115 } // namespace vtkm
116 
117 #endif //vtk_m_NewtonsMethod_h
vtkm::NewtonsMethodResult::Valid
bool Valid
True if Newton's method ran into a singularity.
Definition: NewtonsMethod.h:25
vtkm
Groups connected points that have the same field value.
Definition: Atomic.h:19
VTKM_EXEC_CONT
#define VTKM_EXEC_CONT
Definition: ExportMacros.h:52
vtkm::IdComponent
vtkm::Int32 IdComponent
Base type to use to index small lists.
Definition: Types.h:194
vtkm::NewtonsMethodResult
An object returned from NewtonsMethod() that contains the result and other information about the fina...
Definition: NewtonsMethod.h:22
Matrix.h
vtkm::NewtonsMethod
NewtonsMethodResult< ScalarType, Size > NewtonsMethod(JacobianFunctor jacobianEvaluator, FunctionFunctor functionEvaluator, vtkm::Vec< ScalarType, Size > desiredFunctionOutput, vtkm::Vec< ScalarType, Size > initialGuess=vtkm::Vec< ScalarType, Size >(ScalarType(0)), ScalarType convergeDifference=ScalarType(1e-3), vtkm::IdComponent maxIterations=10)
Uses Newton's method (a.k.a.
Definition: NewtonsMethod.h:66
Math.h
vtkm::NewtonsMethodResult::Solution
vtkm::Vec< ScalarType, Size > Solution
The solution found by Newton's method.
Definition: NewtonsMethod.h:31
vtkm::SolveLinearSystem
vtkm::Vec< T, Size > SolveLinearSystem(const vtkm::Matrix< T, Size, Size > &A, const vtkm::Vec< T, Size > &b, bool &valid)
Solve the linear system Ax = b for x.
Definition: Matrix.h:448
vtkm::Vec< ScalarType, Size >
vtkm::Matrix
Basic Matrix type.
Definition: Matrix.h:33
vtkm::NewtonsMethodResult::Converged
bool Converged
True if Newton's method converted to below the convergence value.
Definition: NewtonsMethod.h:27
VTKM_SUPPRESS_EXEC_WARNINGS
#define VTKM_SUPPRESS_EXEC_WARNINGS
Definition: ExportMacros.h:53