Prusa Slicer 2.6.0
Loading...
Searching...
No Matches
Eigen::JacobiRotation< Scalar > Class Template Reference

Rotation given by a cosine-sine pair. More...

#include <src/eigen/Eigen/src/Jacobi/Jacobi.h>

Public Types

typedef NumTraits< Scalar >::Real RealScalar
 

Public Member Functions

 JacobiRotation ()
 
 JacobiRotation (const Scalar &c, const Scalar &s)
 
Scalar & c ()
 
Scalar c () const
 
Scalar & s ()
 
Scalar s () const
 
JacobiRotation operator* (const JacobiRotation &other)
 
JacobiRotation transpose () const
 
JacobiRotation adjoint () const
 
template<typename Derived >
bool makeJacobi (const MatrixBase< Derived > &, Index p, Index q)
 
bool makeJacobi (const RealScalar &x, const Scalar &y, const RealScalar &z)
 
void makeGivens (const Scalar &p, const Scalar &q, Scalar *r=0)
 

Protected Member Functions

void makeGivens (const Scalar &p, const Scalar &q, Scalar *r, internal::true_type)
 
void makeGivens (const Scalar &p, const Scalar &q, Scalar *r, internal::false_type)
 

Protected Attributes

Scalar m_c
 
Scalar m_s
 

Detailed Description

template<typename Scalar>
class Eigen::JacobiRotation< Scalar >

Rotation given by a cosine-sine pair.

\jacobi_module

This class represents a Jacobi or Givens rotation. This is a 2D rotation in the plane J of angle $ \theta $ defined by its cosine c and sine s as follow: $ J = \left ( \begin{array}{cc} c & \overline s \\ -s  & \overline c \end{array} \right ) $

You can apply the respective counter-clockwise rotation to a column vector v by applying its adjoint on the left: $ v = J^* v $ that translates to the following Eigen code:

v.applyOnTheLeft(J.adjoint());
See also
MatrixBase::applyOnTheLeft(), MatrixBase::applyOnTheRight()

Member Typedef Documentation

◆ RealScalar

template<typename Scalar >
typedef NumTraits<Scalar>::Real Eigen::JacobiRotation< Scalar >::RealScalar

Constructor & Destructor Documentation

◆ JacobiRotation() [1/2]

template<typename Scalar >
Eigen::JacobiRotation< Scalar >::JacobiRotation ( )
inline

Default constructor without any initialization.

40{}

Referenced by Eigen::JacobiRotation< Scalar >::adjoint(), Eigen::JacobiRotation< Scalar >::operator*(), and Eigen::JacobiRotation< Scalar >::transpose().

+ Here is the caller graph for this function:

◆ JacobiRotation() [2/2]

template<typename Scalar >
Eigen::JacobiRotation< Scalar >::JacobiRotation ( const Scalar &  c,
const Scalar &  s 
)
inline

Construct a planar rotation from a cosine-sine pair (c, s).

43: m_c(c), m_s(s) {}
Scalar m_c
Definition Jacobi.h:74
Scalar & s()
Definition Jacobi.h:47
Scalar & c()
Definition Jacobi.h:45
Scalar m_s
Definition Jacobi.h:74

Member Function Documentation

◆ adjoint()

template<typename Scalar >
JacobiRotation Eigen::JacobiRotation< Scalar >::adjoint ( ) const
inline

Returns the adjoint transformation

62{ using numext::conj; return JacobiRotation(conj(m_c), -m_s); }
JacobiRotation()
Definition Jacobi.h:40

References Eigen::JacobiRotation< Scalar >::JacobiRotation(), Eigen::JacobiRotation< Scalar >::m_c, and Eigen::JacobiRotation< Scalar >::m_s.

Referenced by Eigen::RealQZ< _MatrixType >::hessenbergTriangular(), Eigen::RealQZ< _MatrixType >::pushDownZero(), Eigen::ComplexSchur< _MatrixType >::reduceToTriangularForm(), Eigen::internal::svd_precondition_2x2_block_to_be_real< MatrixType, QRPreconditioner, true >::run(), Eigen::RealQZ< _MatrixType >::splitOffTwoRows(), Eigen::RealSchur< _MatrixType >::splitOffTwoRows(), and Eigen::RealQZ< _MatrixType >::step().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ c() [1/2]

template<typename Scalar >
Scalar & Eigen::JacobiRotation< Scalar >::c ( )
inline
45{ return m_c; }

References Eigen::JacobiRotation< Scalar >::m_c.

Referenced by Eigen::MatrixBase< Homogeneous< MatrixType, _Direction > >::applyOnTheRight(), Eigen::internal::real_2x2_jacobi_svd(), Eigen::internal::svd_precondition_2x2_block_to_be_real< MatrixType, QRPreconditioner, true >::run(), and Eigen::internal::tridiagonal_qr_step().

+ Here is the caller graph for this function:

◆ c() [2/2]

template<typename Scalar >
Scalar Eigen::JacobiRotation< Scalar >::c ( ) const
inline
46{ return m_c; }

References Eigen::JacobiRotation< Scalar >::m_c.

◆ makeGivens() [1/3]

template<typename Scalar >
void Eigen::JacobiRotation< Scalar >::makeGivens ( const Scalar &  p,
const Scalar &  q,
Scalar *  r,
internal::false_type   
)
protected
216{
217 using std::sqrt;
218 using std::abs;
219 if(q==Scalar(0))
220 {
221 m_c = p<Scalar(0) ? Scalar(-1) : Scalar(1);
222 m_s = Scalar(0);
223 if(r) *r = abs(p);
224 }
225 else if(p==Scalar(0))
226 {
227 m_c = Scalar(0);
228 m_s = q<Scalar(0) ? Scalar(1) : Scalar(-1);
229 if(r) *r = abs(q);
230 }
231 else if(abs(p) > abs(q))
232 {
233 Scalar t = q/p;
234 Scalar u = sqrt(Scalar(1) + numext::abs2(t));
235 if(p<Scalar(0))
236 u = -u;
237 m_c = Scalar(1)/u;
238 m_s = -t * m_c;
239 if(r) *r = p * u;
240 }
241 else
242 {
243 Scalar t = p/q;
244 Scalar u = sqrt(Scalar(1) + numext::abs2(t));
245 if(q<Scalar(0))
246 u = -u;
247 m_s = -Scalar(1)/u;
248 m_c = -t * m_s;
249 if(r) *r = q * u;
250 }
251
252}
EIGEN_DEVICE_FUNC const SqrtReturnType sqrt() const
Definition ArrayCwiseUnaryOps.h:152
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC half abs(const half &a)
Definition Half.h:445
typename Traits< remove_cvref_t< L > >::Scalar Scalar
Definition Line.hpp:36

References sqrt().

+ Here is the call graph for this function:

◆ makeGivens() [2/3]

template<typename Scalar >
void Eigen::JacobiRotation< Scalar >::makeGivens ( const Scalar &  p,
const Scalar &  q,
Scalar *  r,
internal::true_type   
)
protected
157{
158 using std::sqrt;
159 using std::abs;
160 using numext::conj;
161
162 if(q==Scalar(0))
163 {
164 m_c = numext::real(p)<0 ? Scalar(-1) : Scalar(1);
165 m_s = 0;
166 if(r) *r = m_c * p;
167 }
168 else if(p==Scalar(0))
169 {
170 m_c = 0;
171 m_s = -q/abs(q);
172 if(r) *r = abs(q);
173 }
174 else
175 {
176 RealScalar p1 = numext::norm1(p);
177 RealScalar q1 = numext::norm1(q);
178 if(p1>=q1)
179 {
180 Scalar ps = p / p1;
181 RealScalar p2 = numext::abs2(ps);
182 Scalar qs = q / p1;
183 RealScalar q2 = numext::abs2(qs);
184
185 RealScalar u = sqrt(RealScalar(1) + q2/p2);
186 if(numext::real(p)<RealScalar(0))
187 u = -u;
188
189 m_c = Scalar(1)/u;
190 m_s = -qs*conj(ps)*(m_c/p2);
191 if(r) *r = p * u;
192 }
193 else
194 {
195 Scalar ps = p / q1;
196 RealScalar p2 = numext::abs2(ps);
197 Scalar qs = q / q1;
198 RealScalar q2 = numext::abs2(qs);
199
200 RealScalar u = q1 * sqrt(p2 + q2);
201 if(numext::real(p)<RealScalar(0))
202 u = -u;
203
204 p1 = abs(p);
205 ps = p/p1;
206 m_c = p1/u;
207 m_s = -conj(ps) * (q/u);
208 if(r) *r = ps * u;
209 }
210 }
211}
NumTraits< Scalar >::Real RealScalar
Definition Jacobi.h:37

References sqrt().

+ Here is the call graph for this function:

◆ makeGivens() [3/3]

template<typename Scalar >
void Eigen::JacobiRotation< Scalar >::makeGivens ( const Scalar &  p,
const Scalar &  q,
Scalar *  r = 0 
)

Makes *this as a Givens rotation G such that applying $ G^* $ to the left of the vector $ V = \left ( \begin{array}{c} p \\ q \end{array} \right )$ yields: $ G^* V = \left ( \begin{array}{c} r \\ 0 \end{array} \right )$.

The value of r is returned if r is not null (the default is null). Also note that G is built such that the cosine is always real.

Example:

Output:

This function implements the continuous Givens rotation generation algorithm found in Anderson (2000), Discontinuous Plane Rotations and the Symmetric Eigenvalue Problem. LAPACK Working Note 150, University of Tennessee, UT-CS-00-454, December 4, 2000.

See also
MatrixBase::applyOnTheLeft(), MatrixBase::applyOnTheRight()
149{
150 makeGivens(p, q, r, typename internal::conditional<NumTraits<Scalar>::IsComplex, internal::true_type, internal::false_type>::type());
151}
void makeGivens(const Scalar &p, const Scalar &q, Scalar *r=0)
Definition Jacobi.h:148

Referenced by Eigen::RealQZ< _MatrixType >::hessenbergTriangular(), Eigen::internal::llt_rank_update_lower(), Eigen::RealQZ< _MatrixType >::pushDownZero(), Eigen::ComplexSchur< _MatrixType >::reduceToTriangularForm(), Eigen::RealQZ< _MatrixType >::splitOffTwoRows(), Eigen::RealSchur< _MatrixType >::splitOffTwoRows(), Eigen::RealQZ< _MatrixType >::step(), and Eigen::internal::tridiagonal_qr_step().

+ Here is the caller graph for this function:

◆ makeJacobi() [1/2]

template<typename Scalar >
template<typename Derived >
bool Eigen::JacobiRotation< Scalar >::makeJacobi ( const MatrixBase< Derived > &  m,
Index  p,
Index  q 
)
inline

Makes *this as a Jacobi rotation J such that applying J on both the right and left sides of the 2x2 selfadjoint matrix $ B = \left ( \begin{array}{cc} \text{this}_{pp} & \text{this}_{pq} \\ (\text{this}_{pq})^* & \text{this}_{qq} \end{array} \right )$ yields a diagonal matrix $ A = J^* B J $

Example:

Output:

See also
JacobiRotation::makeJacobi(RealScalar, Scalar, RealScalar), MatrixBase::applyOnTheLeft(), MatrixBase::applyOnTheRight()
127{
128 return makeJacobi(numext::real(m.coeff(p,p)), m.coeff(p,q), numext::real(m.coeff(q,q)));
129}
bool makeJacobi(const MatrixBase< Derived > &, Index p, Index q)
Definition Jacobi.h:126

Referenced by Eigen::internal::real_2x2_jacobi_svd().

+ Here is the caller graph for this function:

◆ makeJacobi() [2/2]

template<typename Scalar >
bool Eigen::JacobiRotation< Scalar >::makeJacobi ( const RealScalar x,
const Scalar &  y,
const RealScalar z 
)

Makes *this as a Jacobi rotation J such that applying J on both the right and left sides of the selfadjoint 2x2 matrix $ B = \left ( \begin{array}{cc} x & y \\ \overline y & z \end{array} \right )$ yields a diagonal matrix $ A = J^* B J $

See also
MatrixBase::makeJacobi(const MatrixBase<Derived>&, Index, Index), MatrixBase::applyOnTheLeft(), MatrixBase::applyOnTheRight()
84{
85 using std::sqrt;
86 using std::abs;
87 RealScalar deno = RealScalar(2)*abs(y);
88 if(deno < (std::numeric_limits<RealScalar>::min)())
89 {
90 m_c = Scalar(1);
91 m_s = Scalar(0);
92 return false;
93 }
94 else
95 {
96 RealScalar tau = (x-z)/deno;
97 RealScalar w = sqrt(numext::abs2(tau) + RealScalar(1));
98 RealScalar t;
99 if(tau>RealScalar(0))
100 {
101 t = RealScalar(1) / (tau + w);
102 }
103 else
104 {
105 t = RealScalar(1) / (tau - w);
106 }
107 RealScalar sign_t = t > RealScalar(0) ? RealScalar(1) : RealScalar(-1);
108 RealScalar n = RealScalar(1) / sqrt(numext::abs2(t)+RealScalar(1));
109 m_s = - sign_t * (numext::conj(y) / abs(y)) * abs(t) * n;
110 m_c = n;
111 return true;
112 }
113}
TCoord< P > x(const P &p)
Definition geometry_traits.hpp:297

References sqrt().

+ Here is the call graph for this function:

◆ operator*()

template<typename Scalar >
JacobiRotation Eigen::JacobiRotation< Scalar >::operator* ( const JacobiRotation< Scalar > &  other)
inline

Concatenates two planar rotation

52 {
53 using numext::conj;
54 return JacobiRotation(m_c * other.m_c - conj(m_s) * other.m_s,
55 conj(m_c * conj(other.m_s) + conj(m_s) * conj(other.m_c)));
56 }

References Eigen::JacobiRotation< Scalar >::JacobiRotation(), Eigen::JacobiRotation< Scalar >::m_c, and Eigen::JacobiRotation< Scalar >::m_s.

+ Here is the call graph for this function:

◆ s() [1/2]

template<typename Scalar >
Scalar & Eigen::JacobiRotation< Scalar >::s ( )
inline
47{ return m_s; }

References Eigen::JacobiRotation< Scalar >::m_s.

Referenced by Eigen::MatrixBase< Homogeneous< MatrixType, _Direction > >::applyOnTheRight(), Eigen::internal::real_2x2_jacobi_svd(), Eigen::internal::svd_precondition_2x2_block_to_be_real< MatrixType, QRPreconditioner, true >::run(), and Eigen::internal::tridiagonal_qr_step().

+ Here is the caller graph for this function:

◆ s() [2/2]

template<typename Scalar >
Scalar Eigen::JacobiRotation< Scalar >::s ( ) const
inline
48{ return m_s; }

References Eigen::JacobiRotation< Scalar >::m_s.

◆ transpose()

template<typename Scalar >
JacobiRotation Eigen::JacobiRotation< Scalar >::transpose ( ) const
inline

Returns the transposed transformation

59{ using numext::conj; return JacobiRotation(m_c, -conj(m_s)); }

References Eigen::JacobiRotation< Scalar >::JacobiRotation(), Eigen::JacobiRotation< Scalar >::m_c, and Eigen::JacobiRotation< Scalar >::m_s.

Referenced by Eigen::MatrixBase< Derived >::applyOnTheRight(), Eigen::RealQZ< _MatrixType >::compute(), Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::compute(), and Eigen::internal::real_2x2_jacobi_svd().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

Member Data Documentation

◆ m_c

◆ m_s


The documentation for this class was generated from the following files: