dogleg.h
Go to the documentation of this file.
2 
3 namespace Eigen {
4 
5 namespace internal {
6 
7 template <typename Scalar>
8 void dogleg(
10  const Matrix< Scalar, Dynamic, 1 > &diag,
12  Scalar delta,
14 {
15  using std::abs;
16  using std::sqrt;
17 
18  typedef DenseIndex Index;
19 
20  /* Local variables */
21  Index i, j;
22  Scalar sum, temp, alpha, bnorm;
23  Scalar gnorm, qnorm;
24  Scalar sgnorm;
25 
26  /* Function Body */
27  const Scalar epsmch = NumTraits<Scalar>::epsilon();
28  const Index n = qrfac.cols();
29  eigen_assert(n==qtb.size());
30  eigen_assert(n==x.size());
31  eigen_assert(n==diag.size());
32  Matrix< Scalar, Dynamic, 1 > wa1(n), wa2(n);
33 
34  /* first, calculate the gauss-newton direction. */
35  for (j = n-1; j >=0; --j) {
36  temp = qrfac(j,j);
37  if (temp == 0.) {
38  temp = epsmch * qrfac.col(j).head(j+1).maxCoeff();
39  if (temp == 0.)
40  temp = epsmch;
41  }
42  if (j==n-1)
43  x[j] = qtb[j] / temp;
44  else
45  x[j] = (qtb[j] - qrfac.row(j).tail(n-j-1).dot(x.tail(n-j-1))) / temp;
46  }
47 
48  /* test whether the gauss-newton direction is acceptable. */
49  qnorm = diag.cwiseProduct(x).stableNorm();
50  if (qnorm <= delta)
51  return;
52 
53  // TODO : this path is not tested by Eigen unit tests
54 
55  /* the gauss-newton direction is not acceptable. */
56  /* next, calculate the scaled gradient direction. */
57 
58  wa1.fill(0.);
59  for (j = 0; j < n; ++j) {
60  wa1.tail(n-j) += qrfac.row(j).tail(n-j) * qtb[j];
61  wa1[j] /= diag[j];
62  }
63 
64  /* calculate the norm of the scaled gradient and test for */
65  /* the special case in which the scaled gradient is zero. */
66  gnorm = wa1.stableNorm();
67  sgnorm = 0.;
68  alpha = delta / qnorm;
69  if (gnorm == 0.)
70  goto algo_end;
71 
72  /* calculate the point along the scaled gradient */
73  /* at which the quadratic is minimized. */
74  wa1.array() /= (diag*gnorm).array();
75  // TODO : once unit tests cover this part,:
76  // wa2 = qrfac.template triangularView<Upper>() * wa1;
77  for (j = 0; j < n; ++j) {
78  sum = 0.;
79  for (i = j; i < n; ++i) {
80  sum += qrfac(j,i) * wa1[i];
81  }
82  wa2[j] = sum;
83  }
84  temp = wa2.stableNorm();
85  sgnorm = gnorm / temp / temp;
86 
87  /* test whether the scaled gradient direction is acceptable. */
88  alpha = 0.;
89  if (sgnorm >= delta)
90  goto algo_end;
91 
92  /* the scaled gradient direction is not acceptable. */
93  /* finally, calculate the point along the dogleg */
94  /* at which the quadratic is minimized. */
95  bnorm = qtb.stableNorm();
96  temp = bnorm / gnorm * (bnorm / qnorm) * (sgnorm / delta);
97  temp = temp - delta / qnorm * numext::abs2(sgnorm / delta) + sqrt(numext::abs2(temp - delta / qnorm) + (1.-numext::abs2(delta / qnorm)) * (1.-numext::abs2(sgnorm / delta)));
98  alpha = delta / qnorm * (1. - numext::abs2(sgnorm / delta)) / temp;
99 algo_end:
100 
101  /* form appropriate convex combination of the gauss-newton */
102  /* direction and the scaled gradient direction. */
103  temp = (1.-alpha) * (std::min)(sgnorm,delta);
104  x = temp * wa1 + alpha * x;
105 }
106 
107 } // end namespace internal
108 
109 } // end namespace Eigen
int n
int i
#define eigen_assert(x)
EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
void dogleg(const Matrix< Scalar, Dynamic, Dynamic > &qrfac, const Matrix< Scalar, Dynamic, 1 > &diag, const Matrix< Scalar, Dynamic, 1 > &qtb, Scalar delta, Matrix< Scalar, Dynamic, 1 > &x)
Definition: dogleg.h:8
bool abs2(bool x)
: TensorContractionSycl.h, provides various tensor contraction kernel for SYCL backend
std::array< T, N > array
Eigen::AutoDiffScalar< EIGEN_EXPR_BINARYOP_SCALAR_RETURN_TYPE(Eigen::internal::remove_all_t< DerType >, typename Eigen::internal::traits< Eigen::internal::remove_all_t< DerType >>::Scalar, product) > sqrt(const Eigen::AutoDiffScalar< DerType > &x)
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
CleanedUpDerType< DerType >::type() min(const AutoDiffScalar< DerType > &x, const T &y)
EIGEN_DEFAULT_DENSE_INDEX_TYPE DenseIndex
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
adouble abs(const adouble &x)
Definition: AdolcForward:74
std::ptrdiff_t j