Jacobi.h
Go to the documentation of this file.
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2009 Benoit Jacob <jacob.benoit.1@gmail.com>
5 // Copyright (C) 2009 Gael Guennebaud <gael.guennebaud@inria.fr>
6 //
7 // This Source Code Form is subject to the terms of the Mozilla
8 // Public License v. 2.0. If a copy of the MPL was not distributed
9 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
10 
11 #ifndef EIGEN_JACOBI_H
12 #define EIGEN_JACOBI_H
13 
14 // IWYU pragma: private
15 #include "./InternalHeaderCheck.h"
16 
17 namespace Eigen {
18 
37 template <typename Scalar>
39  public:
41 
44 
47 
48  EIGEN_DEVICE_FUNC Scalar& c() { return m_c; }
49  EIGEN_DEVICE_FUNC Scalar c() const { return m_c; }
50  EIGEN_DEVICE_FUNC Scalar& s() { return m_s; }
51  EIGEN_DEVICE_FUNC Scalar s() const { return m_s; }
52 
55  using numext::conj;
56  return JacobiRotation(m_c * other.m_c - conj(m_s) * other.m_s,
57  conj(m_c * conj(other.m_s) + conj(m_s) * conj(other.m_c)));
58  }
59 
62  using numext::conj;
63  return JacobiRotation(m_c, -conj(m_s));
64  }
65 
68  using numext::conj;
69  return JacobiRotation(conj(m_c), -m_s);
70  }
71 
72  template <typename Derived>
74  EIGEN_DEVICE_FUNC bool makeJacobi(const RealScalar& x, const Scalar& y, const RealScalar& z);
75 
76  EIGEN_DEVICE_FUNC void makeGivens(const Scalar& p, const Scalar& q, Scalar* r = 0);
77 
78  protected:
81 
83 };
84 
92 template <typename Scalar>
94  using std::abs;
95  using std::sqrt;
96 
97  RealScalar deno = RealScalar(2) * abs(y);
98  if (deno < (std::numeric_limits<RealScalar>::min)()) {
99  m_c = Scalar(1);
100  m_s = Scalar(0);
101  return false;
102  } else {
103  RealScalar tau = (x - z) / deno;
104  RealScalar w = sqrt(numext::abs2(tau) + RealScalar(1));
105  RealScalar t;
106  if (tau > RealScalar(0)) {
107  t = RealScalar(1) / (tau + w);
108  } else {
109  t = RealScalar(1) / (tau - w);
110  }
111  RealScalar sign_t = t > RealScalar(0) ? RealScalar(1) : RealScalar(-1);
113  m_s = -sign_t * (numext::conj(y) / abs(y)) * abs(t) * n;
114  m_c = n;
115  return true;
116  }
117 }
118 
129 template <typename Scalar>
130 template <typename Derived>
132  return makeJacobi(numext::real(m.coeff(p, p)), m.coeff(p, q), numext::real(m.coeff(q, q)));
133 }
134 
151 template <typename Scalar>
154 }
155 
156 // specialization for complexes
157 template <typename Scalar>
160  using numext::conj;
161  using std::abs;
162  using std::sqrt;
163 
164  if (q == Scalar(0)) {
165  m_c = numext::real(p) < 0 ? Scalar(-1) : Scalar(1);
166  m_s = 0;
167  if (r) *r = m_c * p;
168  } else if (p == Scalar(0)) {
169  m_c = 0;
170  m_s = -q / abs(q);
171  if (r) *r = abs(q);
172  } else {
173  RealScalar p1 = numext::norm1(p);
174  RealScalar q1 = numext::norm1(q);
175  if (p1 >= q1) {
176  Scalar ps = p / p1;
177  RealScalar p2 = numext::abs2(ps);
178  Scalar qs = q / p1;
179  RealScalar q2 = numext::abs2(qs);
180 
181  RealScalar u = sqrt(RealScalar(1) + q2 / p2);
182  if (numext::real(p) < RealScalar(0)) u = -u;
183 
184  m_c = Scalar(1) / u;
185  m_s = -qs * conj(ps) * (m_c / p2);
186  if (r) *r = p * u;
187  } else {
188  Scalar ps = p / q1;
189  RealScalar p2 = numext::abs2(ps);
190  Scalar qs = q / q1;
191  RealScalar q2 = numext::abs2(qs);
192 
193  RealScalar u = q1 * sqrt(p2 + q2);
194  if (numext::real(p) < RealScalar(0)) u = -u;
195 
196  p1 = abs(p);
197  ps = p / p1;
198  m_c = p1 / u;
199  m_s = -conj(ps) * (q / u);
200  if (r) *r = ps * u;
201  }
202  }
203 }
204 
205 // specialization for reals
206 template <typename Scalar>
209  using std::abs;
210  using std::sqrt;
211  if (numext::is_exactly_zero(q)) {
212  m_c = p < Scalar(0) ? Scalar(-1) : Scalar(1);
213  m_s = Scalar(0);
214  if (r) *r = abs(p);
215  } else if (numext::is_exactly_zero(p)) {
216  m_c = Scalar(0);
217  m_s = q < Scalar(0) ? Scalar(1) : Scalar(-1);
218  if (r) *r = abs(q);
219  } else if (abs(p) > abs(q)) {
220  Scalar t = q / p;
221  Scalar u = sqrt(Scalar(1) + numext::abs2(t));
222  if (p < Scalar(0)) u = -u;
223  m_c = Scalar(1) / u;
224  m_s = -t * m_c;
225  if (r) *r = p * u;
226  } else {
227  Scalar t = p / q;
228  Scalar u = sqrt(Scalar(1) + numext::abs2(t));
229  if (q < Scalar(0)) u = -u;
230  m_s = -Scalar(1) / u;
231  m_c = -t * m_s;
232  if (r) *r = q * u;
233  }
234 }
235 
236 /****************************************************************************************
237  * Implementation of MatrixBase methods
238  ****************************************************************************************/
239 
240 namespace internal {
248 template <typename VectorX, typename VectorY, typename OtherScalar>
251 } // namespace internal
252 
259 template <typename Derived>
260 template <typename OtherScalar>
263  RowXpr x(this->row(p));
264  RowXpr y(this->row(q));
266 }
267 
274 template <typename Derived>
275 template <typename OtherScalar>
278  ColXpr x(this->col(p));
279  ColXpr y(this->col(q));
281 }
282 
283 namespace internal {
284 
285 template <typename Scalar, typename OtherScalar, int SizeAtCompileTime, int MinAlignment, bool Vectorizable>
287  static EIGEN_DEVICE_FUNC inline void run(Scalar* x, Index incrx, Scalar* y, Index incry, Index size, OtherScalar c,
288  OtherScalar s) {
289  for (Index i = 0; i < size; ++i) {
290  Scalar xi = *x;
291  Scalar yi = *y;
292  *x = c * xi + numext::conj(s) * yi;
293  *y = -s * xi + numext::conj(c) * yi;
294  x += incrx;
295  y += incry;
296  }
297  }
298 };
299 
300 template <typename Scalar, typename OtherScalar, int SizeAtCompileTime, int MinAlignment>
301 struct apply_rotation_in_the_plane_selector<Scalar, OtherScalar, SizeAtCompileTime, MinAlignment,
302  true /* vectorizable */> {
303  static inline void run(Scalar* x, Index incrx, Scalar* y, Index incry, Index size, OtherScalar c, OtherScalar s) {
304  typedef typename packet_traits<Scalar>::type Packet;
305  typedef typename packet_traits<OtherScalar>::type OtherPacket;
306 
307  constexpr int RequiredAlignment =
309  constexpr Index PacketSize = packet_traits<Scalar>::size;
310 
311  /*** dynamic-size vectorized paths ***/
312  if (size >= 2 * PacketSize && SizeAtCompileTime == Dynamic && ((incrx == 1 && incry == 1) || PacketSize == 1)) {
313  // both vectors are sequentially stored in memory => vectorization
314  constexpr Index Peeling = 2;
315 
316  Index alignedStart = internal::first_default_aligned(y, size);
317  Index alignedEnd = alignedStart + ((size - alignedStart) / PacketSize) * PacketSize;
318 
319  const OtherPacket pc = pset1<OtherPacket>(c);
320  const OtherPacket ps = pset1<OtherPacket>(s);
323 
324  for (Index i = 0; i < alignedStart; ++i) {
325  Scalar xi = x[i];
326  Scalar yi = y[i];
327  x[i] = c * xi + numext::conj(s) * yi;
328  y[i] = -s * xi + numext::conj(c) * yi;
329  }
330 
331  Scalar* EIGEN_RESTRICT px = x + alignedStart;
332  Scalar* EIGEN_RESTRICT py = y + alignedStart;
333 
334  if (internal::first_default_aligned(x, size) == alignedStart) {
335  for (Index i = alignedStart; i < alignedEnd; i += PacketSize) {
336  Packet xi = pload<Packet>(px);
337  Packet yi = pload<Packet>(py);
338  pstore(px, padd(pm.pmul(pc, xi), pcj.pmul(ps, yi)));
339  pstore(py, psub(pcj.pmul(pc, yi), pm.pmul(ps, xi)));
340  px += PacketSize;
341  py += PacketSize;
342  }
343  } else {
344  Index peelingEnd = alignedStart + ((size - alignedStart) / (Peeling * PacketSize)) * (Peeling * PacketSize);
345  for (Index i = alignedStart; i < peelingEnd; i += Peeling * PacketSize) {
346  Packet xi = ploadu<Packet>(px);
347  Packet xi1 = ploadu<Packet>(px + PacketSize);
348  Packet yi = pload<Packet>(py);
349  Packet yi1 = pload<Packet>(py + PacketSize);
350  pstoreu(px, padd(pm.pmul(pc, xi), pcj.pmul(ps, yi)));
351  pstoreu(px + PacketSize, padd(pm.pmul(pc, xi1), pcj.pmul(ps, yi1)));
352  pstore(py, psub(pcj.pmul(pc, yi), pm.pmul(ps, xi)));
353  pstore(py + PacketSize, psub(pcj.pmul(pc, yi1), pm.pmul(ps, xi1)));
354  px += Peeling * PacketSize;
355  py += Peeling * PacketSize;
356  }
357  if (alignedEnd != peelingEnd) {
358  Packet xi = ploadu<Packet>(x + peelingEnd);
359  Packet yi = pload<Packet>(y + peelingEnd);
360  pstoreu(x + peelingEnd, padd(pm.pmul(pc, xi), pcj.pmul(ps, yi)));
361  pstore(y + peelingEnd, psub(pcj.pmul(pc, yi), pm.pmul(ps, xi)));
362  }
363  }
364 
365  for (Index i = alignedEnd; i < size; ++i) {
366  Scalar xi = x[i];
367  Scalar yi = y[i];
368  x[i] = c * xi + numext::conj(s) * yi;
369  y[i] = -s * xi + numext::conj(c) * yi;
370  }
371  }
372 
373  /*** fixed-size vectorized path ***/
374  else if (SizeAtCompileTime != Dynamic && MinAlignment >= RequiredAlignment) {
375  const OtherPacket pc = pset1<OtherPacket>(c);
376  const OtherPacket ps = pset1<OtherPacket>(s);
381  for (Index i = 0; i < size; i += PacketSize) {
382  Packet xi = pload<Packet>(px);
383  Packet yi = pload<Packet>(py);
384  pstore(px, padd(pm.pmul(pc, xi), pcj.pmul(ps, yi)));
385  pstore(py, psub(pcj.pmul(pc, yi), pm.pmul(ps, xi)));
386  px += PacketSize;
387  py += PacketSize;
388  }
389  }
390 
391  /*** non-vectorized path ***/
392  else {
394  x, incrx, y, incry, size, c, s);
395  }
396  }
397 };
398 
399 template <typename VectorX, typename VectorY, typename OtherScalar>
402  typedef typename VectorX::Scalar Scalar;
403  constexpr bool Vectorizable = (int(evaluator<VectorX>::Flags) & int(evaluator<VectorY>::Flags) & PacketAccessBit) &&
405 
406  eigen_assert(xpr_x.size() == xpr_y.size());
407  Index size = xpr_x.size();
408  Index incrx = xpr_x.derived().innerStride();
409  Index incry = xpr_y.derived().innerStride();
410 
411  Scalar* EIGEN_RESTRICT x = &xpr_x.derived().coeffRef(0);
412  Scalar* EIGEN_RESTRICT y = &xpr_y.derived().coeffRef(0);
413 
414  OtherScalar c = j.c();
415  OtherScalar s = j.s();
417 
418  constexpr int Alignment = (std::min)(int(evaluator<VectorX>::Alignment), int(evaluator<VectorY>::Alignment));
420  x, incrx, y, incry, size, c, s);
421 }
422 
423 } // end namespace internal
424 
425 } // end namespace Eigen
426 
427 #endif // EIGEN_JACOBI_H
AnnoyingScalar abs(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:135
AnnoyingScalar conj(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:133
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
Matrix4d pm
Definition: HessenbergDecomposition_packedMatrix.cpp:4
G makeGivens(v.x(), v.y())
J makeJacobi(m, 0, 1)
#define EIGEN_RESTRICT
Definition: Macros.h:1067
#define EIGEN_DEVICE_FUNC
Definition: Macros.h:892
#define eigen_assert(x)
Definition: Macros.h:910
Vector3f p1
Definition: MatrixBase_all.cpp:2
m col(1)
m row(1)
RowVector3d w
Definition: Matrix_resize_int.cpp:3
float * p
Definition: Tutorial_Map_using.cpp:9
Scalar Scalar int size
Definition: benchVecAdd.cpp:17
SCALAR Scalar
Definition: bench_gemm.cpp:45
NumTraits< Scalar >::Real RealScalar
Definition: bench_gemm.cpp:46
internal::packet_traits< Scalar >::type Packet
Definition: benchmark-blocking-sizes.cpp:54
Base class for all dense matrices, vectors, and arrays.
Definition: DenseBase.h:44
Rotation given by a cosine-sine pair.
Definition: Jacobi.h:38
EIGEN_DEVICE_FUNC Scalar & s()
Definition: Jacobi.h:50
EIGEN_DEVICE_FUNC void makeGivens(const Scalar &p, const Scalar &q, Scalar *r=0)
Definition: Jacobi.h:152
EIGEN_DEVICE_FUNC JacobiRotation transpose() const
Definition: Jacobi.h:61
Scalar m_c
Definition: Jacobi.h:82
EIGEN_DEVICE_FUNC Scalar c() const
Definition: Jacobi.h:49
EIGEN_DEVICE_FUNC JacobiRotation operator*(const JacobiRotation &other)
Definition: Jacobi.h:54
NumTraits< Scalar >::Real RealScalar
Definition: Jacobi.h:40
EIGEN_DEVICE_FUNC JacobiRotation()
Definition: Jacobi.h:43
EIGEN_DEVICE_FUNC Scalar & c()
Definition: Jacobi.h:48
EIGEN_DEVICE_FUNC JacobiRotation adjoint() const
Definition: Jacobi.h:67
EIGEN_DEVICE_FUNC Scalar s() const
Definition: Jacobi.h:51
EIGEN_DEVICE_FUNC JacobiRotation(const Scalar &c, const Scalar &s)
Definition: Jacobi.h:46
Scalar m_s
Definition: Jacobi.h:82
EIGEN_DEVICE_FUNC bool makeJacobi(const MatrixBase< Derived > &, Index p, Index q)
Definition: Jacobi.h:131
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
void applyOnTheLeft(const EigenBase< OtherDerived > &other)
Definition: MatrixBase.h:532
void applyOnTheRight(const EigenBase< OtherDerived > &other)
Definition: MatrixBase.h:521
Base::RowXpr RowXpr
Definition: MatrixBase.h:89
Base::ColXpr ColXpr
Definition: MatrixBase.h:90
internal::traits< Derived >::Scalar Scalar
Definition: PlainObjectBase.h:127
float real
Definition: datatypes.h:10
#define min(a, b)
Definition: datatypes.h:22
#define max(a, b)
Definition: datatypes.h:23
const unsigned int PacketAccessBit
Definition: Constants.h:97
RealScalar RealScalar * px
Definition: level1_cplx_impl.h:27
Scalar * y
Definition: level1_cplx_impl.h:128
RealScalar s
Definition: level1_cplx_impl.h:130
return int(ret)+1
int RealScalar int RealScalar int RealScalar RealScalar * ps
Definition: level1_cplx_impl.h:124
int RealScalar int RealScalar * py
Definition: level1_cplx_impl.h:124
int RealScalar int RealScalar int RealScalar * pc
Definition: level1_cplx_impl.h:124
int * m
Definition: level2_cplx_impl.h:294
EIGEN_DEVICE_FUNC Packet padd(const Packet &a, const Packet &b)
Definition: GenericPacketMath.h:318
const Scalar & y
Definition: RandomImpl.h:36
static Index first_default_aligned(const DenseBase< Derived > &m)
Definition: DenseCoeffsBase.h:539
EIGEN_DEVICE_FUNC void pstore(Scalar *to, const Packet &from)
Definition: GenericPacketMath.h:891
EIGEN_DEVICE_FUNC void apply_rotation_in_the_plane(DenseBase< VectorX > &xpr_x, DenseBase< VectorY > &xpr_y, const JacobiRotation< OtherScalar > &j)
Definition: Jacobi.h:400
EIGEN_DEVICE_FUNC void pstoreu(Scalar *to, const Packet &from)
Definition: GenericPacketMath.h:911
EIGEN_DEVICE_FUNC Packet psub(const Packet &a, const Packet &b)
Definition: GenericPacketMath.h:337
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC bool is_exactly_zero(const X &x)
Definition: Meta.h:592
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC bool is_exactly_one(const X &x)
Definition: Meta.h:601
EIGEN_DEVICE_FUNC const Scalar & q
Definition: SpecialFunctionsImpl.h:2019
EIGEN_DEVICE_FUNC bool abs2(bool x)
Definition: MathFunctions.h:1102
Namespace containing all symbols from the Eigen library.
Definition: bench_norm.cpp:70
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:83
const AutoDiffScalar< DerType > & conj(const AutoDiffScalar< DerType > &x)
Definition: AutoDiffScalar.h:482
const int Dynamic
Definition: Constants.h:25
r
Definition: UniformPSDSelfTest.py:20
int c
Definition: calibrate.py:100
Definition: Eigen_Colamd.h:49
list x
Definition: plotDoE.py:28
t
Definition: plotPSD.py:36
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:217
static void run(Scalar *x, Index incrx, Scalar *y, Index incry, Index size, OtherScalar c, OtherScalar s)
Definition: Jacobi.h:303
static EIGEN_DEVICE_FUNC void run(Scalar *x, Index incrx, Scalar *y, Index incry, Index size, OtherScalar c, OtherScalar s)
Definition: Jacobi.h:287
Definition: ConjHelper.h:71
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE ResultType pmul(const LhsType &x, const RhsType &y) const
Definition: ConjHelper.h:79
Definition: CoreEvaluators.h:104
Definition: Meta.h:97
Definition: GenericPacketMath.h:108
Definition: Meta.h:94
Definition: GenericPacketMath.h:134
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2
Definition: ZVector/PacketMath.h:50