Blender V4.3
ConstrainedConjugateGradient.h
Go to the documentation of this file.
1/* SPDX-FileCopyrightText: Blender Authors
2 *
3 * SPDX-License-Identifier: GPL-2.0-or-later */
4
5#pragma once
6
11#include <Eigen/Core>
12
13namespace Eigen {
14
15namespace internal {
16
29template<typename MatrixType,
30 typename Rhs,
31 typename Dest,
32 typename FilterMatrixType,
33 typename Preconditioner>
34EIGEN_DONT_INLINE void constrained_conjugate_gradient(const MatrixType &mat,
35 const Rhs &rhs,
36 Dest &x,
37 const FilterMatrixType &filter,
38 const Preconditioner &precond,
39 int &iters,
40 typename Dest::RealScalar &tol_error)
41{
42 using std::abs;
43 using std::sqrt;
44 typedef typename Dest::RealScalar RealScalar;
45 typedef typename Dest::Scalar Scalar;
46 typedef Matrix<Scalar, Dynamic, 1> VectorType;
47
48 RealScalar tol = tol_error;
49 int maxIters = iters;
50
51 int n = mat.cols();
52
53 VectorType residual = filter * (rhs - mat * x); /* initial residual */
54
55 RealScalar rhsNorm2 = (filter * rhs).squaredNorm();
56 if (rhsNorm2 == 0) {
57 /* XXX TODO: set constrained result here. */
58 x.setZero();
59 iters = 0;
60 tol_error = 0;
61 return;
62 }
63 RealScalar threshold = tol * tol * rhsNorm2;
64 RealScalar residualNorm2 = residual.squaredNorm();
65 if (residualNorm2 < threshold) {
66 iters = 0;
67 tol_error = sqrt(residualNorm2 / rhsNorm2);
68 return;
69 }
70
71 VectorType p(n);
72 p = filter * precond.solve(residual); /* initial search direction */
73
74 VectorType z(n), tmp(n);
75 RealScalar absNew = numext::real(
76 residual.dot(p)); /* The square of the absolute value of `r` scaled by `invM`. */
77 int i = 0;
78 while (i < maxIters) {
79 tmp.noalias() = filter * (mat * p); /* The bottleneck of the algorithm. */
80
81 Scalar alpha = absNew / p.dot(tmp); /* The amount we travel on direction. */
82 x += alpha * p; /* Update solution. */
83 residual -= alpha * tmp; /* Update residue. */
84
85 residualNorm2 = residual.squaredNorm();
86 if (residualNorm2 < threshold) {
87 break;
88 }
89
90 z = precond.solve(residual); /* Approximately solve for `A z = residual`. */
91
92 RealScalar absOld = absNew;
93 absNew = numext::real(residual.dot(z)); /* Update the absolute value of `r`. */
94
95 /* Calculate the Gram-Schmidt value used to create the new search direction. */
96 RealScalar beta = absNew / absOld;
97
98 p = filter * (z + beta * p); /* Update search direction. */
99 i++;
100 }
101 tol_error = sqrt(residualNorm2 / rhsNorm2);
102 iters = i;
103}
104
105} // namespace internal
106
107#if 0 /* unused */
108template<typename MatrixType> struct MatrixFilter {
109 MatrixFilter() : m_cmat(NULL) {}
110
111 MatrixFilter(const MatrixType &cmat) : m_cmat(&cmat) {}
112
113 void setMatrix(const MatrixType &cmat)
114 {
115 m_cmat = &cmat;
116 }
117
118 template<typename VectorType> void apply(VectorType v) const
119 {
120 v = (*m_cmat) * v;
121 }
122
123 protected:
124 const MatrixType *m_cmat;
125};
126#endif
127
128template<typename _MatrixType,
129 int _UpLo = Lower,
130 typename _FilterMatrixType = _MatrixType,
131 typename _Preconditioner = DiagonalPreconditioner<typename _MatrixType::Scalar>>
132class ConstrainedConjugateGradient;
133
134namespace internal {
135
136template<typename _MatrixType, int _UpLo, typename _FilterMatrixType, typename _Preconditioner>
137struct traits<
138 ConstrainedConjugateGradient<_MatrixType, _UpLo, _FilterMatrixType, _Preconditioner>> {
139 typedef _MatrixType MatrixType;
140 typedef _FilterMatrixType FilterMatrixType;
141 typedef _Preconditioner Preconditioner;
142};
143
144} // namespace internal
145
195template<typename _MatrixType, int _UpLo, typename _FilterMatrixType, typename _Preconditioner>
197 : public IterativeSolverBase<
198 ConstrainedConjugateGradient<_MatrixType, _UpLo, _FilterMatrixType, _Preconditioner>> {
199 typedef IterativeSolverBase<ConstrainedConjugateGradient> Base;
200 using Base::m_error;
201 using Base::m_info;
202 using Base::m_isInitialized;
203 using Base::m_iterations;
204 using Base::mp_matrix;
205
206 public:
207 typedef _MatrixType MatrixType;
208 typedef typename MatrixType::Scalar Scalar;
209 typedef typename MatrixType::Index Index;
210 typedef typename MatrixType::RealScalar RealScalar;
211 typedef _FilterMatrixType FilterMatrixType;
212 typedef _Preconditioner Preconditioner;
213
214 enum { UpLo = _UpLo };
215
216 public:
219
232
234
236 {
237 return m_filter;
238 }
240 {
241 return m_filter;
242 }
243
250 template<typename Rhs, typename Guess>
251 inline const internal::solve_retval_with_guess<ConstrainedConjugateGradient, Rhs, Guess>
252 solveWithGuess(const MatrixBase<Rhs> &b, const Guess &x0) const
253 {
254 eigen_assert(m_isInitialized && "ConjugateGradient is not initialized.");
255 eigen_assert(
256 Base::rows() == b.rows() &&
257 "ConjugateGradient::solve(): invalid number of rows of the right hand side matrix b");
258 return internal::solve_retval_with_guess<ConstrainedConjugateGradient, Rhs, Guess>(
259 *this, b.derived(), x0);
260 }
261
263 template<typename Rhs, typename Dest> void _solveWithGuess(const Rhs &b, Dest &x) const
264 {
265 m_iterations = Base::maxIterations();
266 m_error = Base::m_tolerance;
267
268 for (int j = 0; j < b.cols(); j++) {
269 m_iterations = Base::maxIterations();
270 m_error = Base::m_tolerance;
271
272 typename Dest::ColXpr xj(x, j);
273 internal::constrained_conjugate_gradient(mp_matrix->template selfadjointView<UpLo>(),
274 b.col(j),
275 xj,
276 m_filter,
277 Base::m_preconditioner,
278 m_iterations,
279 m_error);
280 }
281
282 m_isInitialized = true;
283 m_info = m_error <= Base::m_tolerance ? Success : NoConvergence;
284 }
285
287 template<typename Rhs, typename Dest> void _solve(const Rhs &b, Dest &x) const
288 {
289 x.setOnes();
290 _solveWithGuess(b, x);
291 }
292
293 protected:
295};
296
297namespace internal {
298
299template<typename _MatrixType, int _UpLo, typename _Filter, typename _Preconditioner, typename Rhs>
300struct solve_retval<ConstrainedConjugateGradient<_MatrixType, _UpLo, _Filter, _Preconditioner>,
301 Rhs>
302 : solve_retval_base<ConstrainedConjugateGradient<_MatrixType, _UpLo, _Filter, _Preconditioner>,
303 Rhs> {
305 EIGEN_MAKE_SOLVE_HELPERS(Dec, Rhs)
306
307 template<typename Dest> void evalTo(Dest &dst) const
308 {
309 dec()._solve(rhs(), dst);
310 }
311};
312
313} // end namespace internal
314
315} // end namespace Eigen
sqrt(x)+1/max(0
ATTR_WARN_UNUSED_RESULT const BMVert * v
SIMD_FORCE_INLINE const btScalar & z() const
Return the z value.
Definition btQuadWord.h:117
A conjugate gradient solver for sparse self-adjoint problems with additional constraints.
const 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
local_group_size(16, 16) .push_constant(Type b
local_group_size(16, 16) .push_constant(Type rhs
#define NULL
float Scalar
Definition eigen_utils.h:27
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)
ccl_device_inline float beta(float x, float y)
Definition util/math.h:833