12#include <Eigen/IterativeLinearSolvers>
13#include <Eigen/Sparse>
31template<
typename MatrixType,
34 typename FilterMatrixType,
39 const FilterMatrixType &
filter,
42 typename Dest::RealScalar &tol_error)
46 using RealScalar =
typename Dest::RealScalar;
47 using Scalar =
typename Dest::Scalar;
48 using VectorType = Matrix<Scalar, Dynamic, 1>;
50 RealScalar tol = tol_error;
55 VectorType residual =
filter * (
rhs - mat *
x);
57 RealScalar rhsNorm2 = (
filter *
rhs).squaredNorm();
65 RealScalar threshold = tol * tol * rhsNorm2;
66 RealScalar residualNorm2 = residual.squaredNorm();
67 if (residualNorm2 < threshold) {
69 tol_error =
sqrt(residualNorm2 / rhsNorm2);
74 p =
filter * precond.solve(residual);
76 VectorType
z(n), tmp(n);
77 RealScalar absNew = numext::real(
80 while (
i < maxIters) {
81 tmp.noalias() =
filter * (mat * p);
83 Scalar alpha = absNew / p.dot(tmp);
85 residual -= alpha * tmp;
87 residualNorm2 = residual.squaredNorm();
88 if (residualNorm2 < threshold) {
92 z = precond.solve(residual);
94 RealScalar absOld = absNew;
95 absNew = numext::real(residual.dot(
z));
98 RealScalar
beta = absNew / absOld;
103 tol_error =
sqrt(residualNorm2 / rhsNorm2);
110template<
typename MatrixType>
struct MatrixFilter {
111 MatrixFilter() : m_cmat(NULL) {}
113 MatrixFilter(
const MatrixType &cmat) : m_cmat(&cmat) {}
115 void setMatrix(
const MatrixType &cmat)
120 template<
typename VectorType>
void apply(VectorType
v)
const
126 const MatrixType *m_cmat;
130template<
typename _MatrixType,
132 typename _FilterMatrixType = _MatrixType,
133 typename _Preconditioner = DiagonalPreconditioner<typename _MatrixType::Scalar>>
138template<
typename _MatrixType,
int _UpLo,
typename _FilterMatrixType,
typename _Preconditioner>
197template<
typename _MatrixType,
int _UpLo,
typename _FilterMatrixType,
typename _Preconditioner>
199 :
public IterativeSolverBase<
200 ConstrainedConjugateGradient<_MatrixType, _UpLo, _FilterMatrixType, _Preconditioner>> {
201 using Base = IterativeSolverBase<ConstrainedConjugateGradient>;
204 using Base::m_isInitialized;
205 using Base::m_iterations;
206 using Base::mp_matrix;
210 using Scalar =
typename MatrixType::Scalar;
211 using Index =
typename MatrixType::Index;
251 template<
typename Rhs,
typename Guess>
252 internal::solve_retval_with_guess<ConstrainedConjugateGradient, Rhs, Guess>
solveWithGuess(
253 const MatrixBase<Rhs> &
b,
const Guess &x0)
const
255 eigen_assert(m_isInitialized &&
"ConjugateGradient is not initialized.");
257 Base::rows() ==
b.rows() &&
258 "ConjugateGradient::solve(): invalid number of rows of the right hand side matrix b");
259 return internal::solve_retval_with_guess<ConstrainedConjugateGradient, Rhs, Guess>(
260 *
this,
b.derived(), x0);
266 m_iterations = Base::maxIterations();
267 m_error = Base::m_tolerance;
269 for (
int j = 0; j <
b.cols(); j++) {
270 m_iterations = Base::maxIterations();
271 m_error = Base::m_tolerance;
273 typename Dest::ColXpr xj(
x, j);
278 Base::m_preconditioner,
283 m_isInitialized =
true;
284 m_info = m_error <= Base::m_tolerance ?
Success : NoConvergence;
288 template<
typename Rhs,
typename Dest>
void _solve(
const Rhs &
b, Dest &
x)
const
300template<
typename _MatrixType,
int _UpLo,
typename _Filter,
typename _Preconditioner,
typename Rhs>
303 : solve_retval_base<ConstrainedConjugateGradient<_MatrixType, _UpLo, _Filter, _Preconditioner>,
306 EIGEN_MAKE_SOLVE_HELPERS(
Dec, Rhs)
308 template<typename Dest>
void evalTo(Dest &dst)
const
310 dec()._solve(
rhs(), dst);
ATTR_WARN_UNUSED_RESULT const BMVert * v
SIMD_FORCE_INLINE const btScalar & z() const
Return the z value.
A conjugate gradient solver for sparse self-adjoint problems with additional constraints.
_FilterMatrixType FilterMatrixType
ConstrainedConjugateGradient()
void _solve(const Rhs &b, Dest &x) const
typename MatrixType::RealScalar RealScalar
FilterMatrixType & filter()
_Preconditioner Preconditioner
typename MatrixType::Index Index
internal::solve_retval_with_guess< ConstrainedConjugateGradient, Rhs, Guess > solveWithGuess(const MatrixBase< Rhs > &b, const Guess &x0) const
void _solveWithGuess(const Rhs &b, Dest &x) const
ConstrainedConjugateGradient(const MatrixType &A)
FilterMatrixType m_filter
~ConstrainedConjugateGradient()=default
typename MatrixType::Scalar Scalar
const FilterMatrixType & filter() const
ccl_device_inline float beta(const float x, const float y)
EIGEN_DONT_INLINE void constrained_conjugate_gradient(const MatrixType &mat, const Rhs &rhs, Dest &x, const FilterMatrixType &filter, const Preconditioner &precond, int &iters, typename Dest::RealScalar &tol_error)
void apply(bContext &C, GestureData &gesture_data, wmOperator &op)
void evalTo(Dest &dst) const
ConstrainedConjugateGradient< _MatrixType, _UpLo, _Filter, _Preconditioner > Dec
_Preconditioner Preconditioner
_FilterMatrixType FilterMatrixType