https://mooseframework.inl.gov
NestedSolve.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
10 #include "NestedSolve.h"
11 
12 #include "libmesh/utility.h"
13 
14 template <bool is_ad>
17 {
19 
20  // Newton iteration control parameters
21  params.addParam<Real>("relative_tolerance",
22  relativeToleranceDefault(),
23  "Relative convergence tolerance for Newton iteration");
24  params.addParam<Real>("absolute_tolerance",
25  absoluteToleranceDefault(),
26  "Absolute convergence tolerance for Newton iteration");
27  params.addParam<Real>("step_size_tolerance",
28  xToleranceDefault(),
29  "Minimum step size of linear iterations relative to value of the solution");
30  params.addParam<unsigned int>(
31  "min_iterations",
32  minIterationsDefault(),
33  "Minimum number of nonlinear iterations to execute before accepting convergence");
34  params.addParam<unsigned int>(
35  "max_iterations", maxIterationsDefault(), "Maximum number of nonlinear iterations");
36  params.addParam<Real>("acceptable_multiplier",
37  acceptableMultiplierDefault(),
38  "Factor applied to relative and absolute "
39  "tolerance for acceptable nonlinear convergence if "
40  "iterations are no longer making progress");
41  params.addParam<Real>("damping_factor",
42  dampingFactorDefault(),
43  "Factor applied to step size if guess does not satisfy damping criteria");
44  params.addParam<unsigned int>(
45  "max_damping_iterations",
46  maxDampingIterationsDefault(),
47  "Maximum number of damping steps per linear iteration of nested solve");
48  return params;
49 }
50 
51 template <bool is_ad>
53  : _relative_tolerance_square(Utility::pow<2>(relativeToleranceDefault())),
54  _absolute_tolerance_square(Utility::pow<2>(absoluteToleranceDefault())),
55  _delta_thresh(xToleranceDefault()),
56  _damping_factor(dampingFactorDefault()),
57  _max_damping_iterations(maxDampingIterationsDefault()),
58  _min_iterations(minIterationsDefault()),
59  _max_iterations(maxIterationsDefault()),
60  _acceptable_multiplier(acceptableMultiplierDefault()),
61  _state(State::NONE),
62  _n_iterations(0)
63 {
64 }
65 
66 template <bool is_ad>
68  : _relative_tolerance_square(Utility::pow<2>(params.get<Real>("relative_tolerance"))),
69  _absolute_tolerance_square(Utility::pow<2>(params.get<Real>("absolute_tolerance"))),
70  _delta_thresh(params.get<Real>("step_size_tolerance")),
71  _damping_factor(params.get<Real>("damping_factor")),
72  _max_damping_iterations(params.get<unsigned int>("max_damping_iterations")),
73  _min_iterations(params.get<unsigned int>("min_iterations")),
74  _max_iterations(params.get<unsigned int>("max_iterations")),
75  _acceptable_multiplier(params.get<Real>("acceptable_multiplier")),
76  _state(State::NONE),
77  _n_iterations(0)
78 {
79 }
80 
81 template <bool is_ad>
82 void
86 {
87  const auto N = guess.size();
88  residual.resize(N, 1);
89  jacobian.resize(N, N);
90 }
91 
92 template <bool is_ad>
93 void
96  const NSRealVectorValue & b) const
97 {
98  x = A.inverse() * b;
99 }
100 
101 template <bool is_ad>
102 Real
104 {
106 }
107 
108 template <bool is_ad>
109 Real
111 {
113 }
114 
115 template <bool is_ad>
116 bool
117 NestedSolveTempl<is_ad>::isRelSmall(const NSReal & a, const NSReal & b, const NSReal & c)
118 {
119  return (abs(MetaPhysicL::raw_value(a)) <
121 }
122 
123 template <bool is_ad>
124 bool
126  const NSRealVectorValue & b,
127  const NSReal & c)
128 {
129  for (const auto i : make_range(LIBMESH_DIM))
130  {
131  if (abs(MetaPhysicL::raw_value(a)(i)) >=
133  return false;
134  }
135  return true;
136 }
137 
138 template <bool is_ad>
139 bool
141  const DynamicVector & b,
142  const NSReal & c)
143 {
144  return (a.cwiseAbs().array() < b.cwiseAbs().array() * MetaPhysicL::raw_value(c)).all();
145 }
146 
147 template class NestedSolveTempl<false>;
148 template class NestedSolveTempl<true>;
MetaPhysicL::DualNumber< V, D, asd > abs(const MetaPhysicL::DualNumber< V, D, asd > &a)
Definition: EigenADReal.h:42
Moose::GenericType< RankTwoTensor, is_ad > NSRankTwoTensor
Definition: NestedSolve.h:53
State
possible solver states
Definition: NestedSolve.h:147
static InputParameters validParams()
Definition: NestedSolve.C:16
T * get(const std::unique_ptr< T > &u)
The MooseUtils::get() specializations are used to support making forwards-compatible code changes fro...
Definition: MooseUtils.h:1155
static bool isRelSmall(const V &a, const V &b, const V &c)
Check if |a| < |b * c| for all elements in a and b. This checks if &#39;a&#39; is small relative to...
auto raw_value(const Eigen::Map< T > &in)
Definition: EigenADReal.h:73
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
Eigen::Matrix< NSReal, Eigen::Dynamic, Eigen::Dynamic > DynamicMatrix
Definition: NestedSolve.h:58
InputParameters emptyInputParameters()
void sizeItems(const NestedSolveTempl< is_ad >::DynamicVector &guess, NestedSolveTempl< is_ad >::DynamicVector &residual, NestedSolveTempl< is_ad >::DynamicMatrix &jacobian) const
Size a dynamic Jacobian matrix correctly.
Definition: NestedSolve.C:83
T pow(const T &x)
void linear(const J &A, V &x, const V &b) const
Solve A*x=b for x.
Definition: NestedSolve.h:546
Eigen::Matrix< NSReal, Eigen::Dynamic, 1 > DynamicVector
Eigen type shortcuts.
Definition: NestedSolve.h:57
static Real normSquare(const V &v)
Compute squared norm of v (dropping derivatives as this is only for convergence checking) ...
Definition: NestedSolve.h:555
NestedSolveTempl<is_ad> and its instantiations NestedSolve and ADNestedSolve are utility classes that...
Definition: NestedSolve.h:42
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Moose::GenericType< Real, is_ad > NSReal
AD/non-AD switched type shortcuts.
Definition: NestedSolve.h:51
IntRange< T > make_range(T beg, T end)
void addParam(const std::string &name, const S &value, const std::string &doc_string)
These methods add an optional parameter and a documentation string to the InputParameters object...
Moose::GenericType< RealVectorValue, is_ad > NSRealVectorValue
Definition: NestedSolve.h:52
void ErrorVector unsigned int
CTSub CT_OPERATOR_BINARY CTMul CTCompareLess CTCompareGreater CTCompareEqual _arg template pow< 2 >(tan(_arg))+1.0) *_arg.template D< dtag >()) CT_SIMPLE_UNARY_FUNCTION(sqrt