NewtonsMethod.h 3.59 KB
Newer Older
1 2 3 4 5 6 7 8
//============================================================================
//  Copyright (c) Kitware, Inc.
//  All rights reserved.
//  See LICENSE.txt for details.
//  This software is distributed WITHOUT ANY WARRANTY; without even
//  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
//  PURPOSE.  See the above copyright notice for more information.
//
9
//  Copyright 2015 National Technology & Engineering Solutions of Sandia, LLC (NTESS).
10 11 12
//  Copyright 2015 UT-Battelle, LLC.
//  Copyright 2015 Los Alamos National Security.
//
13
//  Under the terms of Contract DE-NA0003525 with NTESS,
14 15 16 17 18 19
//  the U.S. Government retains certain rights in this software.
//
//  Under the terms of Contract DE-AC52-06NA25396 with Los Alamos National
//  Laboratory (LANL), the U.S. Government retains certain rights in
//  this software.
//============================================================================
20 21
#ifndef vtk_m_NewtonsMethod_h
#define vtk_m_NewtonsMethod_h
22 23 24 25

#include <vtkm/Math.h>
#include <vtkm/Matrix.h>

26 27
namespace vtkm
{
28

29 30 31 32 33 34 35 36
template <typename ScalarType, vtkm::IdComponent Size>
struct NewtonsMethodResult
{
  bool Valid;
  bool Converged;
  vtkm::Vec<ScalarType, Size> Solution;
};

37 38 39 40 41 42 43 44 45 46
/// Uses Newton's method (a.k.a. Newton-Raphson method) to solve a nonlinear
/// system of equations. This function assumes that the number of variables
/// equals the number of equations. Newton's method operates on an iterative
/// evaluate and search. Evaluations are performed using the functors passed
/// into the NewtonsMethod. The first functor returns the NxN matrix of the
/// Jacobian at a given input point. The second functor returns the N tuple
/// that is the function evaluation at the given input point. The input point
/// that evaluates to the desired output, or the closest point found, is
/// returned.
///
47
VTKM_SUPPRESS_EXEC_WARNINGS
48 49 50
template <typename ScalarType,
          vtkm::IdComponent Size,
          typename JacobianFunctor,
51
          typename FunctionFunctor>
52
VTKM_EXEC_CONT NewtonsMethodResult<ScalarType, Size> NewtonsMethod(
53 54
  JacobianFunctor jacobianEvaluator,
  FunctionFunctor functionEvaluator,
55 56
  vtkm::Vec<ScalarType, Size> desiredFunctionOutput,
  vtkm::Vec<ScalarType, Size> initialGuess = vtkm::Vec<ScalarType, Size>(ScalarType(0)),
57 58
  ScalarType convergeDifference = ScalarType(1e-3),
  vtkm::IdComponent maxIterations = 10)
59
{
60 61
  using VectorType = vtkm::Vec<ScalarType, Size>;
  using MatrixType = vtkm::Matrix<ScalarType, Size, Size>;
62 63 64

  VectorType x = initialGuess;

65
  bool valid = false;
66
  bool converged = false;
67 68
  for (vtkm::IdComponent iteration = 0; !converged && (iteration < maxIterations); iteration++)
  {
69 70 71 72 73 74 75 76 77 78 79 80 81
    // For Newton's method, we solve the linear system
    //
    // Jacobian x deltaX = currentFunctionOutput - desiredFunctionOutput
    //
    // The subtraction on the right side simply makes the target of the solve
    // at zero, which is what Newton's method solves for. The deltaX tells us
    // where to move to to solve for a linear system, which we assume will be
    // closer for our nonlinear system.

    MatrixType jacobian = jacobianEvaluator(x);
    VectorType currentFunctionOutput = functionEvaluator(x);

    VectorType deltaX =
82
      vtkm::SolveLinearSystem(jacobian, currentFunctionOutput - desiredFunctionOutput, valid);
83 84 85 86
    if (!valid)
    {
      break;
    }
87 88 89 90 91

    x = x - deltaX;

    converged = true;
    for (vtkm::IdComponent index = 0; index < Size; index++)
92
    {
93 94
      converged &= (vtkm::Abs(deltaX[index]) < convergeDifference);
    }
95
  }
96 97

  // Not checking whether converged.
98
  return { valid, converged, x };
99 100
}

101
} // namespace vtkm
102

103
#endif //vtk_m_NewtonsMethod_h