19 #ifndef EIGEN_LEVENBERGMARQUARDT_H
20 #define EIGEN_LEVENBERGMARQUARDT_H
26 namespace LevenbergMarquardtSpace {
43 template <
typename Scalar_,
int NX=Dynamic,
int NY=Dynamic>
70 template <
typename Scalar_,
typename Index_>
98 template <
typename QRSolver,
typename VectorType>
99 void lmpar2(
const QRSolver &
qr,
const VectorType &diag,
const VectorType &qtb,
100 typename VectorType::Scalar m_delta,
typename VectorType::Scalar &par,
111 template<
typename FunctorType_>
118 typedef typename JacobianType::Scalar
Scalar;
277 template<
typename FunctorType>
283 m_isInitialized =
true;
288 status = minimizeOneStep(
x);
290 m_isInitialized =
true;
294 template<
typename FunctorType>
299 m = m_functor.values();
301 m_wa1.
resize(
n); m_wa2.resize(
n); m_wa3.resize(
n);
307 if (!m_useExternalScaling)
309 eigen_assert( (!m_useExternalScaling || m_diag.size()==
n) &&
"When m_useExternalScaling is set, the caller must provide a valid 'm_diag'");
317 if (
n <= 0 ||
m <
n || m_ftol < 0. || m_xtol < 0. || m_gtol < 0. || m_maxfev <= 0 || m_factor <= 0.){
322 if (m_useExternalScaling)
333 if ( m_functor(
x, m_fvec) < 0)
335 m_fnorm = m_fvec.stableNorm();
344 template<
typename FunctorType>
352 m = m_functor.values();
355 if (
n <= 0 ||
m <
n || tol < 0.)
361 m_maxfev = 100*(
n+1);
367 template<
typename FunctorType>
377 Index m = functor.values();
380 if (
n <= 0 ||
m <
n || tol < 0.)
HouseholderQR< MatrixXf > qr(A)
Performs non linear optimization over a non-linear function, using a variant of the Levenberg Marquar...
RealScalar epsilon() const
static LevenbergMarquardtSpace::Status lmdif1(FunctorType &functor, FVectorType &x, Index *nfev, const Scalar tol=std::sqrt(NumTraits< Scalar >::epsilon()))
JacobianType::Scalar Scalar
FunctorType::QRSolver QRSolver
LevenbergMarquardtSpace::Status minimizeOneStep(FVectorType &x)
ComputationInfo info() const
Reports whether the minimization was successful.
JacobianType::RealScalar RealScalar
void setFtol(RealScalar ftol)
bool m_useExternalScaling
void setExternalScaling(bool value)
LevenbergMarquardtSpace::Status minimize(FVectorType &x)
JacobianType & jacobian()
LevenbergMarquardtSpace::Status lmder1(FVectorType &x, const Scalar tol=std::sqrt(NumTraits< Scalar >::epsilon()))
LevenbergMarquardt(FunctorType &functor)
void setGtol(RealScalar gtol)
RealScalar lm_param(void)
void setEpsilon(RealScalar epsfcn)
void setXtol(RealScalar xtol)
PermutationType m_permutation
void setFactor(RealScalar factor)
LevenbergMarquardtSpace::Status minimizeInit(FVectorType &x)
FunctorType::JacobianType JacobianType
PermutationMatrix< Dynamic, Dynamic, int > PermutationType
PermutationType permutation()
RealScalar factor() const
QRSolver::StorageIndex PermIndex
void setMaxfev(Index maxfev)
Matrix< Scalar, Dynamic, 1 > FVectorType
constexpr void resize(Index rows, Index cols)
@ RelativeReductionTooSmall
@ TooManyFunctionEvaluation
@ ImproperInputParameters
@ RelativeErrorAndReductionTooSmall
void lmpar2(const QRSolver &qr, const VectorType &diag, const VectorType &qtb, typename VectorType::Scalar m_delta, typename VectorType::Scalar &par, VectorType &x)
: TensorContractionSycl.h, provides various tensor contraction kernel for SYCL backend
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
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
Matrix< Scalar, ValuesAtCompileTime, 1 > ValueType
ColPivHouseholderQR< JacobianType > QRSolver
DenseFunctor(int inputs, int values)
Matrix< Scalar, ValuesAtCompileTime, InputsAtCompileTime > JacobianType
Matrix< Scalar, InputsAtCompileTime, 1 > InputType
Matrix< Scalar, Dynamic, 1 > ValueType
SparseMatrix< Scalar, ColMajor, Index > JacobianType
SparseFunctor(int inputs, int values)
Matrix< Scalar, Dynamic, 1 > InputType
SparseQR< JacobianType, COLAMDOrdering< int > > QRSolver