19 template <
typename QRSolver,
typename VectorType>
24 typename VectorType::Scalar m_delta,
25 typename VectorType::Scalar &par,
31 typedef typename QRSolver::MatrixType
MatrixType;
32 typedef typename QRSolver::Scalar
Scalar;
50 const Scalar dwarf = (std::numeric_limits<Scalar>::min)();
63 wa1.tail(
n-rank).setZero();
67 x =
qr.colsPermutation()*wa1;
73 wa2 = diag.cwiseProduct(
x);
86 wa1 =
qr.colsPermutation().inverse() * diag.cwiseProduct(wa2)/
dxnorm;
88 temp = wa1.blueNorm();
89 parl =
fp / m_delta / temp / temp;
93 for (
j = 0;
j <
n; ++
j)
94 wa1[
j] =
s.col(
j).head(
j+1).dot(
qtb.head(
j+1)) / diag[
qr.colsPermutation().indices()(
j)];
96 gnorm = wa1.stableNorm();
97 paru = gnorm / m_delta;
103 par = (std::max)(par,
parl);
104 par = (std::min)(par,
paru);
115 wa1 = sqrt(par)* diag;
120 wa2 = diag.cwiseProduct(
x);
128 if (
abs(
fp) <=
Scalar(0.1) * m_delta || (
parl == 0. &&
fp <= temp && temp < 0.) || iter == 10)
132 wa1 =
qr.colsPermutation().inverse() * diag.cwiseProduct(wa2/
dxnorm);
134 for (
j = 0;
j <
n; ++
j) {
138 wa1[
i] -=
s.coeff(
i,
j) * temp;
140 temp = wa1.blueNorm();
141 parc =
fp / m_delta / temp / temp;
int n
Definition BiCGSTAB_simple.cpp:1
int i
Definition BiCGSTAB_step_by_step.cpp:9
#define eigen_assert(x)
Definition Macros.h:1037
SCALAR Scalar
Definition bench_gemm.cpp:46
MatrixXf MatrixType
Definition benchmark-blocking-sizes.cpp:52
#define abs(x)
Definition datatypes.h:17
set noclip points set clip one set noclip two set bar set border lt lw set xdata set ydata set zdata set x2data set y2data set boxwidth set dummy x
Definition gnuplot_common_settings.hh:12
A triangularView< Lower >().adjoint().solveInPlace(B)
RealScalar s
Definition level1_cplx_impl.h:126
void lmqrsolv(Matrix< Scalar, Rows, Cols > &s, const PermutationMatrix< Dynamic, Dynamic, PermIndex > &iPerm, const Matrix< Scalar, Dynamic, 1 > &diag, const Matrix< Scalar, Dynamic, 1 > &qtb, Matrix< Scalar, Dynamic, 1 > &x, Matrix< Scalar, Dynamic, 1 > &sdiag)
Definition LMqrsolv.h:23
void lmpar2(const QRSolver &qr, const VectorType &diag, const VectorType &qtb, typename VectorType::Scalar m_delta, typename VectorType::Scalar &par, VectorType &x)
Definition LMpar.h:20
Namespace containing all symbols from the Eigen library.
Definition bench_norm.cpp:85
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:74
Definition BandTriangularSolver.h:13
void qr()
Definition qr_colpivoting.cpp:99
Definition ForwardDeclarations.h:17
std::ptrdiff_t j
Definition tut_arithmetic_redux_minmax.cpp:2