Eigen::HybridNonLinearSolver< FunctorType, Scalar > Class Template Reference

Finds a zero of a system of n nonlinear functions in n variables by a modification of the Powell hybrid method ("dogleg"). More...

Classes

struct  Parameters
 

Public Types

typedef Matrix< Scalar, Dynamic, 1 > FVectorType
 
typedef DenseIndex Index
 
typedef Matrix< Scalar, Dynamic, DynamicJacobianType
 
typedef Matrix< Scalar, Dynamic, DynamicUpperTriangularType
 

Public Member Functions

HybridNonLinearSolverSpace::Status hybrd1 (FVectorType &x, const Scalar tol=numext::sqrt(NumTraits< Scalar >::epsilon()))
 
 HybridNonLinearSolver (FunctorType &_functor)
 
HybridNonLinearSolverSpace::Status hybrj1 (FVectorType &x, const Scalar tol=numext::sqrt(NumTraits< Scalar >::epsilon()))
 
void resetParameters (void)
 
HybridNonLinearSolverSpace::Status solve (FVectorType &x)
 
HybridNonLinearSolverSpace::Status solveInit (FVectorType &x)
 
HybridNonLinearSolverSpace::Status solveNumericalDiff (FVectorType &x)
 
HybridNonLinearSolverSpace::Status solveNumericalDiffInit (FVectorType &x)
 
HybridNonLinearSolverSpace::Status solveNumericalDiffOneStep (FVectorType &x)
 
HybridNonLinearSolverSpace::Status solveOneStep (FVectorType &x)
 

Public Attributes

FVectorType diag
 
JacobianType fjac
 
Scalar fnorm
 
FVectorType fvec
 
Index iter
 
Index nfev
 
Index njev
 
Parameters parameters
 
FVectorType qtf
 
UpperTriangularType R
 
bool useExternalScaling
 

Private Member Functions

HybridNonLinearSolveroperator= (const HybridNonLinearSolver &)
 

Private Attributes

Scalar actred
 
Scalar delta
 
Scalar fnorm1
 
FunctorType & functor
 
bool jeval
 
Index n
 
Index ncfail
 
Index ncsuc
 
Index nslow1
 
Index nslow2
 
Scalar pnorm
 
Scalar prered
 
Scalar ratio
 
bool sing
 
Scalar sum
 
Scalar temp
 
FVectorType wa1
 
FVectorType wa2
 
FVectorType wa3
 
FVectorType wa4
 
Scalar xnorm
 

Detailed Description

template<typename FunctorType, typename Scalar = double>
class Eigen::HybridNonLinearSolver< FunctorType, Scalar >

Finds a zero of a system of n nonlinear functions in n variables by a modification of the Powell hybrid method ("dogleg").

The user must provide a subroutine which calculates the functions. The Jacobian is either provided by the user, or approximated using a forward-difference method.

Definition at line 45 of file HybridNonLinearSolver.h.

Member Typedef Documentation

◆ FVectorType

template<typename FunctorType , typename Scalar = double>
typedef Matrix< Scalar, Dynamic, 1 > Eigen::HybridNonLinearSolver< FunctorType, Scalar >::FVectorType

Definition at line 68 of file HybridNonLinearSolver.h.

◆ Index

template<typename FunctorType , typename Scalar = double>
typedef DenseIndex Eigen::HybridNonLinearSolver< FunctorType, Scalar >::Index

Definition at line 48 of file HybridNonLinearSolver.h.

◆ JacobianType

template<typename FunctorType , typename Scalar = double>
typedef Matrix< Scalar, Dynamic, Dynamic > Eigen::HybridNonLinearSolver< FunctorType, Scalar >::JacobianType

Definition at line 69 of file HybridNonLinearSolver.h.

◆ UpperTriangularType

template<typename FunctorType , typename Scalar = double>
typedef Matrix< Scalar, Dynamic, Dynamic > Eigen::HybridNonLinearSolver< FunctorType, Scalar >::UpperTriangularType

Definition at line 71 of file HybridNonLinearSolver.h.

Constructor & Destructor Documentation

◆ HybridNonLinearSolver()

template<typename FunctorType , typename Scalar = double>
Eigen::HybridNonLinearSolver< FunctorType, Scalar >::HybridNonLinearSolver ( FunctorType &  _functor)
inline

Member Function Documentation

◆ hybrd1()

template<typename FunctorType , typename Scalar >
HybridNonLinearSolverSpace::Status Eigen::HybridNonLinearSolver< FunctorType, Scalar >::hybrd1 ( FVectorType x,
const Scalar  tol = numext::sqrt(NumTraits<Scalar>::epsilon()) 
)

Definition at line 360 of file HybridNonLinearSolver.h.

364 {
365  n = x.size();
366 
367  /* check the input parameters for errors. */
368  if (n <= 0 || tol < 0.)
370 
371  resetParameters();
372  parameters.maxfev = 200*(n+1);
373  parameters.xtol = tol;
374 
375  diag.setConstant(n, 1.);
376  useExternalScaling = true;
377  return solveNumericalDiff(x);
378 }
HybridNonLinearSolverSpace::Status solveNumericalDiff(FVectorType &x)
Derived & setConstant(Index rows, Index cols, const Scalar &val)

◆ hybrj1()

template<typename FunctorType , typename Scalar >
HybridNonLinearSolverSpace::Status Eigen::HybridNonLinearSolver< FunctorType, Scalar >::hybrj1 ( FVectorType x,
const Scalar  tol = numext::sqrt(NumTraits<Scalar>::epsilon()) 
)

Definition at line 124 of file HybridNonLinearSolver.h.

128 {
129  n = x.size();
130 
131  /* check the input parameters for errors. */
132  if (n <= 0 || tol < 0.)
134 
135  resetParameters();
136  parameters.maxfev = 100*(n+1);
137  parameters.xtol = tol;
138  diag.setConstant(n, 1.);
139  useExternalScaling = true;
140  return solve(x);
141 }
HybridNonLinearSolverSpace::Status solve(FVectorType &x)

◆ operator=()

template<typename FunctorType , typename Scalar = double>
HybridNonLinearSolver& Eigen::HybridNonLinearSolver< FunctorType, Scalar >::operator= ( const HybridNonLinearSolver< FunctorType, Scalar > &  )
private

◆ resetParameters()

template<typename FunctorType , typename Scalar = double>
void Eigen::HybridNonLinearSolver< FunctorType, Scalar >::resetParameters ( void  )
inline

Definition at line 91 of file HybridNonLinearSolver.h.

91 { parameters = Parameters(); }

◆ solve()

template<typename FunctorType , typename Scalar >
HybridNonLinearSolverSpace::Status Eigen::HybridNonLinearSolver< FunctorType, Scalar >::solve ( FVectorType x)

Definition at line 346 of file HybridNonLinearSolver.h.

347 {
350  return status;
352  status = solveOneStep(x);
353  return status;
354 }
HybridNonLinearSolverSpace::Status solveInit(FVectorType &x)
HybridNonLinearSolverSpace::Status solveOneStep(FVectorType &x)

◆ solveInit()

template<typename FunctorType , typename Scalar >
HybridNonLinearSolverSpace::Status Eigen::HybridNonLinearSolver< FunctorType, Scalar >::solveInit ( FVectorType x)

Definition at line 145 of file HybridNonLinearSolver.h.

146 {
147  n = x.size();
148 
150  fvec.resize(n);
151  qtf.resize(n);
152  fjac.resize(n, n);
153  if (!useExternalScaling)
154  diag.resize(n);
155  eigen_assert( (!useExternalScaling || diag.size()==n) && "When useExternalScaling is set, the caller must provide a valid 'diag'");
156 
157  /* Function Body */
158  nfev = 0;
159  njev = 0;
160 
161  /* check the input parameters for errors. */
162  if (n <= 0 || parameters.xtol < 0. || parameters.maxfev <= 0 || parameters.factor <= 0. )
164  if (useExternalScaling)
165  for (Index j = 0; j < n; ++j)
166  if (diag[j] <= 0.)
168 
169  /* evaluate the function at the starting point */
170  /* and calculate its norm. */
171  nfev = 1;
172  if ( functor(x, fvec) < 0)
174  fnorm = fvec.stableNorm();
175 
176  /* initialize iteration counter and monitors. */
177  iter = 1;
178  ncsuc = 0;
179  ncfail = 0;
180  nslow1 = 0;
181  nslow2 = 0;
182 
184 }
#define eigen_assert(x)
constexpr void resize(Index rows, Index cols)
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
std::ptrdiff_t j

◆ solveNumericalDiff()

template<typename FunctorType , typename Scalar >
HybridNonLinearSolverSpace::Status Eigen::HybridNonLinearSolver< FunctorType, Scalar >::solveNumericalDiff ( FVectorType x)

Definition at line 589 of file HybridNonLinearSolver.h.

590 {
593  return status;
595  status = solveNumericalDiffOneStep(x);
596  return status;
597 }
HybridNonLinearSolverSpace::Status solveNumericalDiffInit(FVectorType &x)
HybridNonLinearSolverSpace::Status solveNumericalDiffOneStep(FVectorType &x)

◆ solveNumericalDiffInit()

template<typename FunctorType , typename Scalar >
HybridNonLinearSolverSpace::Status Eigen::HybridNonLinearSolver< FunctorType, Scalar >::solveNumericalDiffInit ( FVectorType x)

Definition at line 382 of file HybridNonLinearSolver.h.

383 {
384  n = x.size();
385 
388 
390  qtf.resize(n);
391  fjac.resize(n, n);
392  fvec.resize(n);
393  if (!useExternalScaling)
394  diag.resize(n);
395  eigen_assert( (!useExternalScaling || diag.size()==n) && "When useExternalScaling is set, the caller must provide a valid 'diag'");
396 
397  /* Function Body */
398  nfev = 0;
399  njev = 0;
400 
401  /* check the input parameters for errors. */
404  if (useExternalScaling)
405  for (Index j = 0; j < n; ++j)
406  if (diag[j] <= 0.)
408 
409  /* evaluate the function at the starting point */
410  /* and calculate its norm. */
411  nfev = 1;
412  if ( functor(x, fvec) < 0)
414  fnorm = fvec.stableNorm();
415 
416  /* initialize iteration counter and monitors. */
417  iter = 1;
418  ncsuc = 0;
419  ncfail = 0;
420  nslow1 = 0;
421  nslow2 = 0;
422 
424 }

◆ solveNumericalDiffOneStep()

template<typename FunctorType , typename Scalar >
HybridNonLinearSolverSpace::Status Eigen::HybridNonLinearSolver< FunctorType, Scalar >::solveNumericalDiffOneStep ( FVectorType x)

Definition at line 428 of file HybridNonLinearSolver.h.

429 {
430  using std::sqrt;
431  using std::abs;
432 
433  eigen_assert(x.size()==n); // check the caller is not cheating us
434 
435  Index j;
436  std::vector<JacobiRotation<Scalar> > v_givens(n), w_givens(n);
437 
438  jeval = true;
441 
442  /* calculate the jacobian matrix. */
446 
447  wa2 = fjac.colwise().blueNorm();
448 
449  /* on the first iteration and if external scaling is not used, scale according */
450  /* to the norms of the columns of the initial jacobian. */
451  if (iter == 1) {
452  if (!useExternalScaling)
453  for (j = 0; j < n; ++j)
454  diag[j] = (wa2[j]==0.) ? 1. : wa2[j];
455 
456  /* on the first iteration, calculate the norm of the scaled x */
457  /* and initialize the step bound delta. */
458  xnorm = diag.cwiseProduct(x).stableNorm();
460  if (delta == 0.)
462  }
463 
464  /* compute the qr factorization of the jacobian. */
465  HouseholderQR<JacobianType> qrfac(fjac); // no pivoting:
466 
467  /* copy the triangular factor of the qr factorization into r. */
468  R = qrfac.matrixQR();
469 
470  /* accumulate the orthogonal factor in fjac. */
471  fjac = qrfac.householderQ();
472 
473  /* form (q transpose)*fvec and store in qtf. */
474  qtf = fjac.transpose() * fvec;
475 
476  /* rescale if necessary. */
477  if (!useExternalScaling)
478  diag = diag.cwiseMax(wa2);
479 
480  while (true) {
481  /* determine the direction p. */
482  internal::dogleg<Scalar>(R, diag, qtf, delta, wa1);
483 
484  /* store the direction p and x + p. calculate the norm of p. */
485  wa1 = -wa1;
486  wa2 = x + wa1;
487  pnorm = diag.cwiseProduct(wa1).stableNorm();
488 
489  /* on the first iteration, adjust the initial step bound. */
490  if (iter == 1)
491  delta = (std::min)(delta,pnorm);
492 
493  /* evaluate the function at x + p and calculate its norm. */
494  if ( functor(wa2, wa4) < 0)
496  ++nfev;
497  fnorm1 = wa4.stableNorm();
498 
499  /* compute the scaled actual reduction. */
500  actred = -1.;
501  if (fnorm1 < fnorm) /* Computing 2nd power */
502  actred = 1. - numext::abs2(fnorm1 / fnorm);
503 
504  /* compute the scaled predicted reduction. */
505  wa3 = R.template triangularView<Upper>()*wa1 + qtf;
506  temp = wa3.stableNorm();
507  prered = 0.;
508  if (temp < fnorm) /* Computing 2nd power */
509  prered = 1. - numext::abs2(temp / fnorm);
510 
511  /* compute the ratio of the actual to the predicted reduction. */
512  ratio = 0.;
513  if (prered > 0.)
514  ratio = actred / prered;
515 
516  /* update the step bound. */
517  if (ratio < Scalar(.1)) {
518  ncsuc = 0;
519  ++ncfail;
520  delta = Scalar(.5) * delta;
521  } else {
522  ncfail = 0;
523  ++ncsuc;
524  if (ratio >= Scalar(.5) || ncsuc > 1)
525  delta = (std::max)(delta, pnorm / Scalar(.5));
526  if (abs(ratio - 1.) <= Scalar(.1)) {
527  delta = pnorm / Scalar(.5);
528  }
529  }
530 
531  /* test for successful iteration. */
532  if (ratio >= Scalar(1e-4)) {
533  /* successful iteration. update x, fvec, and their norms. */
534  x = wa2;
535  wa2 = diag.cwiseProduct(x);
536  fvec = wa4;
537  xnorm = wa2.stableNorm();
538  fnorm = fnorm1;
539  ++iter;
540  }
541 
542  /* determine the progress of the iteration. */
543  ++nslow1;
544  if (actred >= Scalar(.001))
545  nslow1 = 0;
546  if (jeval)
547  ++nslow2;
548  if (actred >= Scalar(.1))
549  nslow2 = 0;
550 
551  /* test for convergence. */
552  if (delta <= parameters.xtol * xnorm || fnorm == 0.)
554 
555  /* tests for termination and stringent tolerances. */
556  if (nfev >= parameters.maxfev)
558  if (Scalar(.1) * (std::max)(Scalar(.1) * delta, pnorm) <= NumTraits<Scalar>::epsilon() * xnorm)
560  if (nslow2 == 5)
562  if (nslow1 == 10)
564 
565  /* criterion for recalculating jacobian. */
566  if (ncfail == 2)
567  break; // leave inner loop and go for the next outer loop iteration
568 
569  /* calculate the rank one modification to the jacobian */
570  /* and update qtf if necessary. */
571  wa1 = diag.cwiseProduct( diag.cwiseProduct(wa1)/pnorm );
572  wa2 = fjac.transpose() * wa4;
573  if (ratio >= Scalar(1e-4))
574  qtf = wa2;
575  wa2 = (wa2-wa3)/pnorm;
576 
577  /* compute the qr factorization of the updated jacobian. */
578  internal::r1updt<Scalar>(R, wa1, v_givens, w_givens, wa2, wa3, &sing);
579  internal::r1mpyq<Scalar>(n, n, fjac.data(), v_givens, w_givens);
580  internal::r1mpyq<Scalar>(1, n, qtf.data(), v_givens, w_givens);
581 
582  jeval = false;
583  }
585 }
DenseIndex fdjac1(const FunctorType &Functor, Matrix< Scalar, Dynamic, 1 > &x, Matrix< Scalar, Dynamic, 1 > &fvec, Matrix< Scalar, Dynamic, Dynamic > &fjac, DenseIndex ml, DenseIndex mu, Scalar epsfcn)
Definition: fdjac1.h:8
bool abs2(bool x)
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)
CleanedUpDerType< DerType >::type() min(const AutoDiffScalar< DerType > &x, const T &y)
CleanedUpDerType< DerType >::type() max(const AutoDiffScalar< DerType > &x, const T &y)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)
adouble abs(const adouble &x)
Definition: AdolcForward:74

◆ solveOneStep()

template<typename FunctorType , typename Scalar >
HybridNonLinearSolverSpace::Status Eigen::HybridNonLinearSolver< FunctorType, Scalar >::solveOneStep ( FVectorType x)

Definition at line 188 of file HybridNonLinearSolver.h.

189 {
190  using std::abs;
191 
192  eigen_assert(x.size()==n); // check the caller is not cheating us
193 
194  Index j;
195  std::vector<JacobiRotation<Scalar> > v_givens(n), w_givens(n);
196 
197  jeval = true;
198 
199  /* calculate the jacobian matrix. */
200  if ( functor.df(x, fjac) < 0)
202  ++njev;
203 
204  wa2 = fjac.colwise().blueNorm();
205 
206  /* on the first iteration and if external scaling is not used, scale according */
207  /* to the norms of the columns of the initial jacobian. */
208  if (iter == 1) {
209  if (!useExternalScaling)
210  for (j = 0; j < n; ++j)
211  diag[j] = (wa2[j]==0.) ? 1. : wa2[j];
212 
213  /* on the first iteration, calculate the norm of the scaled x */
214  /* and initialize the step bound delta. */
215  xnorm = diag.cwiseProduct(x).stableNorm();
217  if (delta == 0.)
219  }
220 
221  /* compute the qr factorization of the jacobian. */
222  HouseholderQR<JacobianType> qrfac(fjac); // no pivoting:
223 
224  /* copy the triangular factor of the qr factorization into r. */
225  R = qrfac.matrixQR();
226 
227  /* accumulate the orthogonal factor in fjac. */
228  fjac = qrfac.householderQ();
229 
230  /* form (q transpose)*fvec and store in qtf. */
231  qtf = fjac.transpose() * fvec;
232 
233  /* rescale if necessary. */
234  if (!useExternalScaling)
235  diag = diag.cwiseMax(wa2);
236 
237  while (true) {
238  /* determine the direction p. */
239  internal::dogleg<Scalar>(R, diag, qtf, delta, wa1);
240 
241  /* store the direction p and x + p. calculate the norm of p. */
242  wa1 = -wa1;
243  wa2 = x + wa1;
244  pnorm = diag.cwiseProduct(wa1).stableNorm();
245 
246  /* on the first iteration, adjust the initial step bound. */
247  if (iter == 1)
248  delta = (std::min)(delta,pnorm);
249 
250  /* evaluate the function at x + p and calculate its norm. */
251  if ( functor(wa2, wa4) < 0)
253  ++nfev;
254  fnorm1 = wa4.stableNorm();
255 
256  /* compute the scaled actual reduction. */
257  actred = -1.;
258  if (fnorm1 < fnorm) /* Computing 2nd power */
259  actred = 1. - numext::abs2(fnorm1 / fnorm);
260 
261  /* compute the scaled predicted reduction. */
262  wa3 = R.template triangularView<Upper>()*wa1 + qtf;
263  temp = wa3.stableNorm();
264  prered = 0.;
265  if (temp < fnorm) /* Computing 2nd power */
266  prered = 1. - numext::abs2(temp / fnorm);
267 
268  /* compute the ratio of the actual to the predicted reduction. */
269  ratio = 0.;
270  if (prered > 0.)
271  ratio = actred / prered;
272 
273  /* update the step bound. */
274  if (ratio < Scalar(.1)) {
275  ncsuc = 0;
276  ++ncfail;
277  delta = Scalar(.5) * delta;
278  } else {
279  ncfail = 0;
280  ++ncsuc;
281  if (ratio >= Scalar(.5) || ncsuc > 1)
282  delta = (std::max)(delta, pnorm / Scalar(.5));
283  if (abs(ratio - 1.) <= Scalar(.1)) {
284  delta = pnorm / Scalar(.5);
285  }
286  }
287 
288  /* test for successful iteration. */
289  if (ratio >= Scalar(1e-4)) {
290  /* successful iteration. update x, fvec, and their norms. */
291  x = wa2;
292  wa2 = diag.cwiseProduct(x);
293  fvec = wa4;
294  xnorm = wa2.stableNorm();
295  fnorm = fnorm1;
296  ++iter;
297  }
298 
299  /* determine the progress of the iteration. */
300  ++nslow1;
301  if (actred >= Scalar(.001))
302  nslow1 = 0;
303  if (jeval)
304  ++nslow2;
305  if (actred >= Scalar(.1))
306  nslow2 = 0;
307 
308  /* test for convergence. */
309  if (delta <= parameters.xtol * xnorm || fnorm == 0.)
311 
312  /* tests for termination and stringent tolerances. */
313  if (nfev >= parameters.maxfev)
315  if (Scalar(.1) * (std::max)(Scalar(.1) * delta, pnorm) <= NumTraits<Scalar>::epsilon() * xnorm)
317  if (nslow2 == 5)
319  if (nslow1 == 10)
321 
322  /* criterion for recalculating jacobian. */
323  if (ncfail == 2)
324  break; // leave inner loop and go for the next outer loop iteration
325 
326  /* calculate the rank one modification to the jacobian */
327  /* and update qtf if necessary. */
328  wa1 = diag.cwiseProduct( diag.cwiseProduct(wa1)/pnorm );
329  wa2 = fjac.transpose() * wa4;
330  if (ratio >= Scalar(1e-4))
331  qtf = wa2;
332  wa2 = (wa2-wa3)/pnorm;
333 
334  /* compute the qr factorization of the updated jacobian. */
335  internal::r1updt<Scalar>(R, wa1, v_givens, w_givens, wa2, wa3, &sing);
336  internal::r1mpyq<Scalar>(n, n, fjac.data(), v_givens, w_givens);
337  internal::r1mpyq<Scalar>(1, n, qtf.data(), v_givens, w_givens);
338 
339  jeval = false;
340  }
342 }

Member Data Documentation

◆ actred

template<typename FunctorType , typename Scalar = double>
Scalar Eigen::HybridNonLinearSolver< FunctorType, Scalar >::actred
private

Definition at line 114 of file HybridNonLinearSolver.h.

◆ delta

template<typename FunctorType , typename Scalar = double>
Scalar Eigen::HybridNonLinearSolver< FunctorType, Scalar >::delta
private

Definition at line 107 of file HybridNonLinearSolver.h.

◆ diag

template<typename FunctorType , typename Scalar = double>
FVectorType Eigen::HybridNonLinearSolver< FunctorType, Scalar >::diag

Definition at line 93 of file HybridNonLinearSolver.h.

◆ fjac

template<typename FunctorType , typename Scalar = double>
JacobianType Eigen::HybridNonLinearSolver< FunctorType, Scalar >::fjac

Definition at line 94 of file HybridNonLinearSolver.h.

◆ fnorm

template<typename FunctorType , typename Scalar = double>
Scalar Eigen::HybridNonLinearSolver< FunctorType, Scalar >::fnorm

Definition at line 99 of file HybridNonLinearSolver.h.

◆ fnorm1

template<typename FunctorType , typename Scalar = double>
Scalar Eigen::HybridNonLinearSolver< FunctorType, Scalar >::fnorm1
private

Definition at line 111 of file HybridNonLinearSolver.h.

◆ functor

template<typename FunctorType , typename Scalar = double>
FunctorType& Eigen::HybridNonLinearSolver< FunctorType, Scalar >::functor
private

Definition at line 102 of file HybridNonLinearSolver.h.

◆ fvec

template<typename FunctorType , typename Scalar = double>
FVectorType Eigen::HybridNonLinearSolver< FunctorType, Scalar >::fvec

Definition at line 93 of file HybridNonLinearSolver.h.

◆ iter

template<typename FunctorType , typename Scalar = double>
Index Eigen::HybridNonLinearSolver< FunctorType, Scalar >::iter

Definition at line 98 of file HybridNonLinearSolver.h.

◆ jeval

template<typename FunctorType , typename Scalar = double>
bool Eigen::HybridNonLinearSolver< FunctorType, Scalar >::jeval
private

Definition at line 108 of file HybridNonLinearSolver.h.

◆ n

template<typename FunctorType , typename Scalar = double>
Index Eigen::HybridNonLinearSolver< FunctorType, Scalar >::n
private

Definition at line 103 of file HybridNonLinearSolver.h.

◆ ncfail

template<typename FunctorType , typename Scalar = double>
Index Eigen::HybridNonLinearSolver< FunctorType, Scalar >::ncfail
private

Definition at line 113 of file HybridNonLinearSolver.h.

◆ ncsuc

template<typename FunctorType , typename Scalar = double>
Index Eigen::HybridNonLinearSolver< FunctorType, Scalar >::ncsuc
private

Definition at line 109 of file HybridNonLinearSolver.h.

◆ nfev

template<typename FunctorType , typename Scalar = double>
Index Eigen::HybridNonLinearSolver< FunctorType, Scalar >::nfev

Definition at line 96 of file HybridNonLinearSolver.h.

◆ njev

template<typename FunctorType , typename Scalar = double>
Index Eigen::HybridNonLinearSolver< FunctorType, Scalar >::njev

Definition at line 97 of file HybridNonLinearSolver.h.

◆ nslow1

template<typename FunctorType , typename Scalar = double>
Index Eigen::HybridNonLinearSolver< FunctorType, Scalar >::nslow1
private

Definition at line 112 of file HybridNonLinearSolver.h.

◆ nslow2

template<typename FunctorType , typename Scalar = double>
Index Eigen::HybridNonLinearSolver< FunctorType, Scalar >::nslow2
private

Definition at line 112 of file HybridNonLinearSolver.h.

◆ parameters

template<typename FunctorType , typename Scalar = double>
Parameters Eigen::HybridNonLinearSolver< FunctorType, Scalar >::parameters

Definition at line 92 of file HybridNonLinearSolver.h.

◆ pnorm

template<typename FunctorType , typename Scalar = double>
Scalar Eigen::HybridNonLinearSolver< FunctorType, Scalar >::pnorm
private

Definition at line 111 of file HybridNonLinearSolver.h.

◆ prered

template<typename FunctorType , typename Scalar = double>
Scalar Eigen::HybridNonLinearSolver< FunctorType, Scalar >::prered
private

Definition at line 114 of file HybridNonLinearSolver.h.

◆ qtf

template<typename FunctorType , typename Scalar = double>
FVectorType Eigen::HybridNonLinearSolver< FunctorType, Scalar >::qtf

Definition at line 93 of file HybridNonLinearSolver.h.

◆ R

template<typename FunctorType , typename Scalar = double>
UpperTriangularType Eigen::HybridNonLinearSolver< FunctorType, Scalar >::R

Definition at line 95 of file HybridNonLinearSolver.h.

◆ ratio

template<typename FunctorType , typename Scalar = double>
Scalar Eigen::HybridNonLinearSolver< FunctorType, Scalar >::ratio
private

Definition at line 110 of file HybridNonLinearSolver.h.

◆ sing

template<typename FunctorType , typename Scalar = double>
bool Eigen::HybridNonLinearSolver< FunctorType, Scalar >::sing
private

Definition at line 105 of file HybridNonLinearSolver.h.

◆ sum

template<typename FunctorType , typename Scalar = double>
Scalar Eigen::HybridNonLinearSolver< FunctorType, Scalar >::sum
private

Definition at line 104 of file HybridNonLinearSolver.h.

◆ temp

template<typename FunctorType , typename Scalar = double>
Scalar Eigen::HybridNonLinearSolver< FunctorType, Scalar >::temp
private

Definition at line 106 of file HybridNonLinearSolver.h.

◆ useExternalScaling

template<typename FunctorType , typename Scalar = double>
bool Eigen::HybridNonLinearSolver< FunctorType, Scalar >::useExternalScaling

Definition at line 100 of file HybridNonLinearSolver.h.

◆ wa1

template<typename FunctorType , typename Scalar = double>
FVectorType Eigen::HybridNonLinearSolver< FunctorType, Scalar >::wa1
private

Definition at line 115 of file HybridNonLinearSolver.h.

◆ wa2

template<typename FunctorType , typename Scalar = double>
FVectorType Eigen::HybridNonLinearSolver< FunctorType, Scalar >::wa2
private

Definition at line 115 of file HybridNonLinearSolver.h.

◆ wa3

template<typename FunctorType , typename Scalar = double>
FVectorType Eigen::HybridNonLinearSolver< FunctorType, Scalar >::wa3
private

Definition at line 115 of file HybridNonLinearSolver.h.

◆ wa4

template<typename FunctorType , typename Scalar = double>
FVectorType Eigen::HybridNonLinearSolver< FunctorType, Scalar >::wa4
private

Definition at line 115 of file HybridNonLinearSolver.h.

◆ xnorm

template<typename FunctorType , typename Scalar = double>
Scalar Eigen::HybridNonLinearSolver< FunctorType, Scalar >::xnorm
private

Definition at line 111 of file HybridNonLinearSolver.h.


The documentation for this class was generated from the following file: