Blender  V2.93
ConstrainedConjugateGradient.h
Go to the documentation of this file.
1 /*
2  * This program is free software; you can redistribute it and/or
3  * modify it under the terms of the GNU General Public License
4  * as published by the Free Software Foundation; either version 2
5  * of the License, or (at your option) any later version.
6  *
7  * This program is distributed in the hope that it will be useful,
8  * but WITHOUT ANY WARRANTY; without even the implied warranty of
9  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
10  * GNU General Public License for more details.
11  *
12  * You should have received a copy of the GNU General Public License
13  * along with this program; if not, write to the Free Software Foundation,
14  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
15  *
16  * The Original Code is Copyright (C) Blender Foundation
17  * All rights reserved.
18  */
19 
20 #pragma once
21 
26 #include <Eigen/Core>
27 
28 namespace Eigen {
29 
30 namespace internal {
31 
43 template<typename MatrixType,
44  typename Rhs,
45  typename Dest,
46  typename FilterMatrixType,
47  typename Preconditioner>
48 EIGEN_DONT_INLINE void constrained_conjugate_gradient(const MatrixType &mat,
49  const Rhs &rhs,
50  Dest &x,
51  const FilterMatrixType &filter,
52  const Preconditioner &precond,
53  int &iters,
54  typename Dest::RealScalar &tol_error)
55 {
56  using std::abs;
57  using std::sqrt;
58  typedef typename Dest::RealScalar RealScalar;
59  typedef typename Dest::Scalar Scalar;
60  typedef Matrix<Scalar, Dynamic, 1> VectorType;
61 
62  RealScalar tol = tol_error;
63  int maxIters = iters;
64 
65  int n = mat.cols();
66 
67  VectorType residual = filter * (rhs - mat * x); /* initial residual */
68 
69  RealScalar rhsNorm2 = (filter * rhs).squaredNorm();
70  if (rhsNorm2 == 0) {
71  /* XXX TODO set constrained result here */
72  x.setZero();
73  iters = 0;
74  tol_error = 0;
75  return;
76  }
77  RealScalar threshold = tol * tol * rhsNorm2;
78  RealScalar residualNorm2 = residual.squaredNorm();
79  if (residualNorm2 < threshold) {
80  iters = 0;
81  tol_error = sqrt(residualNorm2 / rhsNorm2);
82  return;
83  }
84 
85  VectorType p(n);
86  p = filter * precond.solve(residual); /* initial search direction */
87 
88  VectorType z(n), tmp(n);
89  RealScalar absNew = numext::real(
90  residual.dot(p)); /* the square of the absolute value of r scaled by invM */
91  int i = 0;
92  while (i < maxIters) {
93  tmp.noalias() = filter * (mat * p); /* the bottleneck of the algorithm */
94 
95  Scalar alpha = absNew / p.dot(tmp); /* the amount we travel on dir */
96  x += alpha * p; /* update solution */
97  residual -= alpha * tmp; /* update residue */
98 
99  residualNorm2 = residual.squaredNorm();
100  if (residualNorm2 < threshold) {
101  break;
102  }
103 
104  z = precond.solve(residual); /* approximately solve for "A z = residual" */
105 
106  RealScalar absOld = absNew;
107  absNew = numext::real(residual.dot(z)); /* update the absolute value of r */
108  RealScalar beta =
109  absNew /
110  absOld; /* calculate the Gram-Schmidt value used to create the new search direction */
111  p = filter * (z + beta * p); /* update search direction */
112  i++;
113  }
114  tol_error = sqrt(residualNorm2 / rhsNorm2);
115  iters = i;
116 }
117 
118 } // namespace internal
119 
120 #if 0 /* unused */
121 template<typename MatrixType> struct MatrixFilter {
122  MatrixFilter() : m_cmat(NULL)
123  {
124  }
125 
126  MatrixFilter(const MatrixType &cmat) : m_cmat(&cmat)
127  {
128  }
129 
130  void setMatrix(const MatrixType &cmat)
131  {
132  m_cmat = &cmat;
133  }
134 
135  template<typename VectorType> void apply(VectorType v) const
136  {
137  v = (*m_cmat) * v;
138  }
139 
140  protected:
141  const MatrixType *m_cmat;
142 };
143 #endif
144 
145 template<typename _MatrixType,
146  int _UpLo = Lower,
147  typename _FilterMatrixType = _MatrixType,
148  typename _Preconditioner = DiagonalPreconditioner<typename _MatrixType::Scalar>>
149 class ConstrainedConjugateGradient;
150 
151 namespace internal {
152 
153 template<typename _MatrixType, int _UpLo, typename _FilterMatrixType, typename _Preconditioner>
154 struct traits<
155  ConstrainedConjugateGradient<_MatrixType, _UpLo, _FilterMatrixType, _Preconditioner>> {
156  typedef _MatrixType MatrixType;
157  typedef _FilterMatrixType FilterMatrixType;
158  typedef _Preconditioner Preconditioner;
159 };
160 
161 } // namespace internal
162 
212 template<typename _MatrixType, int _UpLo, typename _FilterMatrixType, typename _Preconditioner>
214  : public IterativeSolverBase<
215  ConstrainedConjugateGradient<_MatrixType, _UpLo, _FilterMatrixType, _Preconditioner>> {
216  typedef IterativeSolverBase<ConstrainedConjugateGradient> Base;
217  using Base::m_error;
218  using Base::m_info;
219  using Base::m_isInitialized;
220  using Base::m_iterations;
221  using Base::mp_matrix;
222 
223  public:
224  typedef _MatrixType MatrixType;
225  typedef typename MatrixType::Scalar Scalar;
226  typedef typename MatrixType::Index Index;
227  typedef typename MatrixType::RealScalar RealScalar;
228  typedef _FilterMatrixType FilterMatrixType;
229  typedef _Preconditioner Preconditioner;
230 
231  enum { UpLo = _UpLo };
232 
233  public:
236  {
237  }
238 
250  {
251  }
252 
254  {
255  }
256 
258  {
259  return m_filter;
260  }
261  const FilterMatrixType &filter() const
262  {
263  return m_filter;
264  }
265 
271  template<typename Rhs, typename Guess>
272  inline const internal::solve_retval_with_guess<ConstrainedConjugateGradient, Rhs, Guess>
273  solveWithGuess(const MatrixBase<Rhs> &b, const Guess &x0) const
274  {
275  eigen_assert(m_isInitialized && "ConjugateGradient is not initialized.");
276  eigen_assert(
277  Base::rows() == b.rows() &&
278  "ConjugateGradient::solve(): invalid number of rows of the right hand side matrix b");
279  return internal::solve_retval_with_guess<ConstrainedConjugateGradient, Rhs, Guess>(
280  *this, b.derived(), x0);
281  }
282 
284  template<typename Rhs, typename Dest> void _solveWithGuess(const Rhs &b, Dest &x) const
285  {
286  m_iterations = Base::maxIterations();
287  m_error = Base::m_tolerance;
288 
289  for (int j = 0; j < b.cols(); j++) {
290  m_iterations = Base::maxIterations();
291  m_error = Base::m_tolerance;
292 
293  typename Dest::ColXpr xj(x, j);
294  internal::constrained_conjugate_gradient(mp_matrix->template selfadjointView<UpLo>(),
295  b.col(j),
296  xj,
297  m_filter,
298  Base::m_preconditioner,
299  m_iterations,
300  m_error);
301  }
302 
303  m_isInitialized = true;
304  m_info = m_error <= Base::m_tolerance ? Success : NoConvergence;
305  }
306 
308  template<typename Rhs, typename Dest> void _solve(const Rhs &b, Dest &x) const
309  {
310  x.setOnes();
311  _solveWithGuess(b, x);
312  }
313 
314  protected:
316 };
317 
318 namespace internal {
319 
320 template<typename _MatrixType, int _UpLo, typename _Filter, typename _Preconditioner, typename Rhs>
321 struct solve_retval<ConstrainedConjugateGradient<_MatrixType, _UpLo, _Filter, _Preconditioner>,
322  Rhs>
323  : solve_retval_base<ConstrainedConjugateGradient<_MatrixType, _UpLo, _Filter, _Preconditioner>,
324  Rhs> {
326  EIGEN_MAKE_SOLVE_HELPERS(Dec, Rhs)
327 
328  template<typename Dest> void evalTo(Dest &dst) const
329  {
330  dec()._solve(rhs(), dst);
331  }
332 };
333 
334 } // end namespace internal
335 
336 } // end namespace Eigen
sqrt(x)+1/max(0
_GL_VOID GLfloat value _GL_VOID_RET _GL_VOID const GLuint GLboolean *residences _GL_BOOL_RET _GL_VOID GLsizei GLfloat GLfloat GLfloat GLfloat const GLubyte *bitmap _GL_VOID_RET _GL_VOID GLenum const void *lists _GL_VOID_RET _GL_VOID const GLdouble *equation _GL_VOID_RET _GL_VOID GLdouble GLdouble blue _GL_VOID_RET _GL_VOID GLfloat GLfloat blue _GL_VOID_RET _GL_VOID GLint GLint blue _GL_VOID_RET _GL_VOID GLshort GLshort blue _GL_VOID_RET _GL_VOID GLubyte GLubyte blue _GL_VOID_RET _GL_VOID GLuint GLuint blue _GL_VOID_RET _GL_VOID GLushort GLushort blue _GL_VOID_RET _GL_VOID GLbyte GLbyte GLbyte alpha _GL_VOID_RET _GL_VOID GLdouble GLdouble GLdouble alpha _GL_VOID_RET _GL_VOID GLfloat GLfloat GLfloat alpha _GL_VOID_RET _GL_VOID GLint GLint GLint alpha _GL_VOID_RET _GL_VOID GLshort GLshort GLshort alpha _GL_VOID_RET _GL_VOID GLubyte GLubyte GLubyte alpha _GL_VOID_RET _GL_VOID GLuint GLuint GLuint alpha _GL_VOID_RET _GL_VOID GLushort GLushort GLushort alpha _GL_VOID_RET _GL_VOID GLenum mode _GL_VOID_RET _GL_VOID GLint GLsizei GLsizei GLenum type _GL_VOID_RET _GL_VOID GLsizei GLenum GLenum const void *pixels _GL_VOID_RET _GL_VOID const void *pointer _GL_VOID_RET _GL_VOID GLdouble v _GL_VOID_RET _GL_VOID GLfloat v _GL_VOID_RET _GL_VOID GLint GLint i2 _GL_VOID_RET _GL_VOID GLint j _GL_VOID_RET _GL_VOID GLfloat param _GL_VOID_RET _GL_VOID GLint param _GL_VOID_RET _GL_VOID GLdouble GLdouble GLdouble GLdouble GLdouble zFar _GL_VOID_RET _GL_UINT GLdouble *equation _GL_VOID_RET _GL_VOID GLenum GLint *params _GL_VOID_RET _GL_VOID GLenum GLfloat *v _GL_VOID_RET _GL_VOID GLenum GLfloat *params _GL_VOID_RET _GL_VOID GLfloat *values _GL_VOID_RET _GL_VOID GLushort *values _GL_VOID_RET _GL_VOID GLenum GLfloat *params _GL_VOID_RET _GL_VOID GLenum GLdouble *params _GL_VOID_RET _GL_VOID GLenum GLint *params _GL_VOID_RET _GL_VOID GLsizei const void *pointer _GL_VOID_RET _GL_VOID GLsizei const void *pointer _GL_VOID_RET _GL_BOOL GLfloat param _GL_VOID_RET _GL_VOID GLint param _GL_VOID_RET _GL_VOID GLenum GLfloat param _GL_VOID_RET _GL_VOID GLenum GLint param _GL_VOID_RET _GL_VOID GLushort pattern _GL_VOID_RET _GL_VOID GLdouble GLdouble GLint GLint const GLdouble *points _GL_VOID_RET _GL_VOID GLdouble GLdouble GLint GLint GLdouble GLdouble GLint GLint const GLdouble *points _GL_VOID_RET _GL_VOID GLdouble GLdouble u2 _GL_VOID_RET _GL_VOID GLdouble GLdouble GLint GLdouble GLdouble v2 _GL_VOID_RET _GL_VOID GLenum GLfloat param _GL_VOID_RET _GL_VOID GLenum GLint param _GL_VOID_RET _GL_VOID GLenum mode _GL_VOID_RET _GL_VOID GLdouble GLdouble nz _GL_VOID_RET _GL_VOID GLfloat GLfloat nz _GL_VOID_RET _GL_VOID GLint GLint nz _GL_VOID_RET _GL_VOID GLshort GLshort nz _GL_VOID_RET _GL_VOID GLsizei const void *pointer _GL_VOID_RET _GL_VOID GLsizei const GLfloat *values _GL_VOID_RET _GL_VOID GLsizei const GLushort *values _GL_VOID_RET _GL_VOID GLint param _GL_VOID_RET _GL_VOID const GLuint const GLclampf *priorities _GL_VOID_RET _GL_VOID GLdouble y _GL_VOID_RET _GL_VOID GLfloat y _GL_VOID_RET _GL_VOID GLint y _GL_VOID_RET _GL_VOID GLshort y _GL_VOID_RET _GL_VOID GLdouble GLdouble z _GL_VOID_RET _GL_VOID GLfloat GLfloat z _GL_VOID_RET _GL_VOID GLint GLint z _GL_VOID_RET _GL_VOID GLshort GLshort z _GL_VOID_RET _GL_VOID GLdouble GLdouble z
ATTR_WARN_UNUSED_RESULT const BMVert * v
#define A
A conjugate gradient solver for sparse self-adjoint problems with additional constraints.
void _solve(const Rhs &b, Dest &x) const
void _solveWithGuess(const Rhs &b, Dest &x) const
const internal::solve_retval_with_guess< ConstrainedConjugateGradient, Rhs, Guess > solveWithGuess(const MatrixBase< Rhs > &b, const Guess &x0) const
static CCL_NAMESPACE_BEGIN const double alpha
float Scalar
Definition: eigen_utils.h:42
DO_INLINE void filter(lfVector *V, fmatrix3x3 *S)
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)
double real
Definition: Precision.h:26
__forceinline const avxi abs(const avxi &a)
Definition: util_avxi.h:186
ccl_device_inline float beta(float x, float y)
Definition: util_math.h:666