Line data Source code
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>
15 : InputParameters
16 14565 : NestedSolveTempl<is_ad>::validParams()
17 : {
18 14565 : InputParameters params = emptyInputParameters();
19 :
20 : // Newton iteration control parameters
21 43695 : params.addParam<Real>("relative_tolerance",
22 29130 : relativeToleranceDefault(),
23 : "Relative convergence tolerance for Newton iteration");
24 43695 : params.addParam<Real>("absolute_tolerance",
25 29130 : absoluteToleranceDefault(),
26 : "Absolute convergence tolerance for Newton iteration");
27 43695 : params.addParam<Real>("step_size_tolerance",
28 29130 : xToleranceDefault(),
29 : "Minimum step size of linear iterations relative to value of the solution");
30 43695 : params.addParam<unsigned int>(
31 : "min_iterations",
32 29130 : minIterationsDefault(),
33 : "Minimum number of nonlinear iterations to execute before accepting convergence");
34 43695 : params.addParam<unsigned int>(
35 29130 : "max_iterations", maxIterationsDefault(), "Maximum number of nonlinear iterations");
36 43695 : params.addParam<Real>("acceptable_multiplier",
37 29130 : acceptableMultiplierDefault(),
38 : "Factor applied to relative and absolute "
39 : "tolerance for acceptable nonlinear convergence if "
40 : "iterations are no longer making progress");
41 43695 : params.addParam<Real>("damping_factor",
42 29130 : dampingFactorDefault(),
43 : "Factor applied to step size if guess does not satisfy damping criteria");
44 43695 : params.addParam<unsigned int>(
45 : "max_damping_iterations",
46 29130 : maxDampingIterationsDefault(),
47 : "Maximum number of damping steps per linear iteration of nested solve");
48 14565 : return params;
49 0 : }
50 :
51 : template <bool is_ad>
52 10 : NestedSolveTempl<is_ad>::NestedSolveTempl()
53 10 : : _relative_tolerance_square(Utility::pow<2>(relativeToleranceDefault())),
54 10 : _absolute_tolerance_square(Utility::pow<2>(absoluteToleranceDefault())),
55 10 : _delta_thresh(xToleranceDefault()),
56 10 : _damping_factor(dampingFactorDefault()),
57 10 : _max_damping_iterations(maxDampingIterationsDefault()),
58 10 : _min_iterations(minIterationsDefault()),
59 10 : _max_iterations(maxIterationsDefault()),
60 10 : _acceptable_multiplier(acceptableMultiplierDefault()),
61 10 : _state(State::NONE),
62 10 : _n_iterations(0)
63 : {
64 10 : }
65 :
66 : template <bool is_ad>
67 230 : NestedSolveTempl<is_ad>::NestedSolveTempl(const InputParameters & params)
68 230 : : _relative_tolerance_square(Utility::pow<2>(params.get<Real>("relative_tolerance"))),
69 230 : _absolute_tolerance_square(Utility::pow<2>(params.get<Real>("absolute_tolerance"))),
70 230 : _delta_thresh(params.get<Real>("step_size_tolerance")),
71 230 : _damping_factor(params.get<Real>("damping_factor")),
72 230 : _max_damping_iterations(params.get<unsigned int>("max_damping_iterations")),
73 230 : _min_iterations(params.get<unsigned int>("min_iterations")),
74 230 : _max_iterations(params.get<unsigned int>("max_iterations")),
75 230 : _acceptable_multiplier(params.get<Real>("acceptable_multiplier")),
76 230 : _state(State::NONE),
77 230 : _n_iterations(0)
78 : {
79 230 : }
80 :
81 : template <bool is_ad>
82 : void
83 1 : NestedSolveTempl<is_ad>::sizeItems(const NestedSolveTempl<is_ad>::DynamicVector & guess,
84 : NestedSolveTempl<is_ad>::DynamicVector & residual,
85 : NestedSolveTempl<is_ad>::DynamicMatrix & jacobian) const
86 : {
87 1 : const auto N = guess.size();
88 1 : residual.resize(N, 1);
89 1 : jacobian.resize(N, N);
90 1 : }
91 :
92 : template <bool is_ad>
93 : void
94 2 : NestedSolveTempl<is_ad>::linear(const NSRankTwoTensor & A,
95 : NSRealVectorValue & x,
96 : const NSRealVectorValue & b) const
97 : {
98 2 : x = A.inverse() * b;
99 2 : }
100 :
101 : template <bool is_ad>
102 : Real
103 121259 : NestedSolveTempl<is_ad>::normSquare(const NSReal & v)
104 : {
105 121259 : return Utility::pow<2>(MetaPhysicL::raw_value(v));
106 : }
107 :
108 : template <bool is_ad>
109 : Real
110 2 : NestedSolveTempl<is_ad>::normSquare(const NSRealVectorValue & v)
111 : {
112 2 : return (MetaPhysicL::raw_value(v) * MetaPhysicL::raw_value(v));
113 : }
114 :
115 : template <bool is_ad>
116 : bool
117 113626 : NestedSolveTempl<is_ad>::isRelSmall(const NSReal & a, const NSReal & b, const NSReal & c)
118 : {
119 113626 : return (abs(MetaPhysicL::raw_value(a)) <
120 113626 : abs(MetaPhysicL::raw_value(c) * MetaPhysicL::raw_value(b)));
121 : }
122 :
123 : template <bool is_ad>
124 : bool
125 2 : NestedSolveTempl<is_ad>::isRelSmall(const NSRealVectorValue & a,
126 : const NSRealVectorValue & b,
127 : const NSReal & c)
128 : {
129 5 : for (const auto i : make_range(LIBMESH_DIM))
130 : {
131 4 : if (abs(MetaPhysicL::raw_value(a)(i)) >=
132 4 : abs(MetaPhysicL::raw_value(b)(i) * MetaPhysicL::raw_value(c)))
133 1 : return false;
134 : }
135 1 : return true;
136 : }
137 :
138 : template <bool is_ad>
139 : bool
140 48 : NestedSolveTempl<is_ad>::isRelSmall(const DynamicVector & a,
141 : const DynamicVector & b,
142 : const NSReal & c)
143 : {
144 48 : return (a.cwiseAbs().array() < b.cwiseAbs().array() * MetaPhysicL::raw_value(c)).all();
145 : }
146 :
147 : template class NestedSolveTempl<false>;
148 : template class NestedSolveTempl<true>;
|