15 #include "../MatOp/internal/ArnoldiOp.h" 16 #include "../Util/TypeTraits.h" 17 #include "../Util/SimpleRandom.h" 30 template <
typename Scalar,
typename ArnoldiOpType>
71 SimpleRandom<Scalar>
rng(seed + 123 *
iter);
72 f.noalias() = rng.random_vec(m_n);
74 m_op.trans_product(V, f, Vf);
75 f.noalias() -= V * Vf;
88 m_op(op), m_n(op.
rows()), m_m(m), m_k(0),
103 void init(MapConstVec&
v0, Index& op_counter)
111 const Scalar v0norm = m_op.norm(v0);
112 if (v0norm < m_near_0)
113 throw std::invalid_argument(
"initial residual vector cannot be zero");
119 v.noalias() = v0 / v0norm;
123 m_op.perform_op(
v.data(), w.
data());
126 m_fac_H(0, 0) = m_op.inner_product(
v, w);
127 m_fac_f.noalias() = w -
v *
m_fac_H(0, 0);
131 if (m_fac_f.cwiseAbs().maxCoeff() <
m_eps)
138 m_beta = m_op.norm(m_fac_f);
155 std::stringstream
msg;
156 msg <<
"Arnoldi: from_k (= " << from_k <<
") is larger than the current subspace dimension (= " << m_k <<
")";
157 throw std::invalid_argument(msg.str());
167 m_fac_H.rightCols(m_m - from_k).
setZero();
168 m_fac_H.block(from_k, 0, m_m - from_k, from_k).
setZero();
170 for (Index
i = from_k;
i <= to_m - 1;
i++)
172 bool restart =
false;
176 if (m_beta < m_near_0)
184 m_fac_V.col(
i).noalias() = m_fac_f /
m_beta;
193 const Index i1 =
i + 1;
195 MapConstMat Vs(m_fac_V.
data(),
m_n, i1);
199 m_op.trans_product(Vs, w, h);
202 m_fac_f.noalias() = w - Vs *
h;
203 m_beta = m_op.norm(m_fac_f);
205 if (m_beta >
Scalar(0.717) * m_op.norm(h))
210 m_op.trans_product(Vs, m_fac_f, Vf.head(i1));
211 Scalar ortho_err = Vf.head(i1).cwiseAbs().maxCoeff();
214 while (count < 5 && ortho_err > m_eps * m_beta)
221 if (m_beta < beta_thresh)
229 m_fac_f.noalias() -= Vs * Vf.head(i1);
231 h.noalias() += Vf.head(i1);
233 m_beta = m_op.norm(m_fac_f);
235 m_op.trans_product(Vs, m_fac_f, Vf.head(i1));
236 ortho_err = Vf.head(i1).cwiseAbs().maxCoeff();
266 Matrix Vs(m_n, m_k + 1);
267 for (Index
i = 0;
i <
m_k;
i++)
269 const Index nnz = m_m - m_k +
i + 1;
270 MapConstVec
q(&
Q(0, i), nnz);
271 Vs.col(i).noalias() = m_fac_V.leftCols(nnz) *
q;
273 Vs.col(m_k).noalias() = m_fac_V * Q.col(m_k);
274 m_fac_V.leftCols(m_k + 1).noalias() = Vs;
276 Vector fk = m_fac_f *
Q(m_m - 1, m_k - 1) + m_fac_V.col(m_k) *
m_fac_H(m_k, m_k - 1);
278 m_beta = m_op.norm(m_fac_f);
EIGEN_DEVICE_FUNC Derived & setZero(Index size)
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const Scalar * data() const
void expand_basis(MapConstMat &V, const Index seed, Vector &f, Scalar &fnorm)
A matrix or vector expression mapping an existing array of data.
Eigen::Matrix< Scalar, Eigen::Dynamic, 1 > Vector
void compress_H(const UpperHessenbergQR< Scalar > &decomp)
Namespace containing all symbols from the Eigen library.
iterator iter(handle obj)
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void swap(DenseBase< OtherDerived > &other)
const Matrix & matrix_V() const
void init(MapConstVec &v0, Index &op_counter)
void matrix_QtHQ(Matrix &dest) const
virtual void matrix_QtHQ(Matrix &dest) const
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Index subspace_dim() const
Eigen::Matrix< Scalar, Eigen::Dynamic, Eigen::Dynamic > Matrix
Array< int, Dynamic, 1 > v
Point2(* f)(const Point3 &, OptionalJacobian< 2, 3 >)
Eigen::Map< Matrix > MapMat
EIGEN_DEVICE_FUNC const Scalar & q
const Matrix & matrix_H() const
void compress_V(const Matrix &Q)
Eigen::Map< const Matrix > MapConstMat
void compress_H(const DoubleShiftQR< Scalar > &decomp)
Eigen::Map< const Vector > MapConstVec
const Vector & vector_f() const
Eigen::Map< Vector > MapVec
Jet< T, N > sqrt(const Jet< T, N > &f)
virtual void factorize_from(Index from_k, Index to_m, Index &op_counter)
Arnoldi(const ArnoldiOpType &op, Index m)