11 #ifndef EIGEN_JACOBI_H
12 #define EIGEN_JACOBI_H
14 #include "./InternalHeaderCheck.h"
49 EIGEN_DEVICE_FUNC Scalar& c() {
return m_c; }
50 EIGEN_DEVICE_FUNC Scalar c()
const {
return m_c; }
51 EIGEN_DEVICE_FUNC Scalar& s() {
return m_s; }
52 EIGEN_DEVICE_FUNC Scalar s()
const {
return m_s; }
71 template<
typename Derived>
75 bool makeJacobi(
const RealScalar& x,
const Scalar& y,
const RealScalar& z);
78 void makeGivens(
const Scalar& p,
const Scalar& q, Scalar* r=0);
82 void makeGivens(
const Scalar& p,
const Scalar& q, Scalar* r, internal::true_type);
84 void makeGivens(
const Scalar& p,
const Scalar& q, Scalar* r, internal::false_type);
94 template<
typename Scalar>
101 RealScalar deno = RealScalar(2)*
abs(y);
102 if(deno < (std::numeric_limits<RealScalar>::min)())
110 RealScalar tau = (x-z)/deno;
111 RealScalar w =
sqrt(numext::abs2(tau) + RealScalar(1));
113 if(tau>RealScalar(0))
115 t = RealScalar(1) / (tau + w);
119 t = RealScalar(1) / (tau - w);
121 RealScalar sign_t = t > RealScalar(0) ? RealScalar(1) : RealScalar(-1);
122 RealScalar n = RealScalar(1) /
sqrt(numext::abs2(t)+RealScalar(1));
123 m_s = - sign_t * (numext::conj(y) /
abs(y)) *
abs(t) * n;
138 template<
typename Scalar>
139 template<
typename Derived>
143 return makeJacobi(numext::real(m.
coeff(p,p)), m.
coeff(p,q), numext::real(m.
coeff(q,q)));
162 template<
typename Scalar>
171 template<
typename Scalar>
181 m_c = numext::real(p)<0 ? Scalar(-1) : Scalar(1);
185 else if(p==Scalar(0))
193 RealScalar p1 = numext::norm1(p);
194 RealScalar q1 = numext::norm1(q);
198 RealScalar p2 = numext::abs2(ps);
200 RealScalar q2 = numext::abs2(qs);
202 RealScalar u =
sqrt(RealScalar(1) + q2/p2);
203 if(numext::real(p)<RealScalar(0))
207 m_s = -qs*
conj(ps)*(m_c/p2);
213 RealScalar p2 = numext::abs2(ps);
215 RealScalar q2 = numext::abs2(qs);
217 RealScalar u = q1 *
sqrt(p2 + q2);
218 if(numext::real(p)<RealScalar(0))
224 m_s = -
conj(ps) * (q/u);
231 template<
typename Scalar>
237 if(numext::is_exactly_zero(q))
239 m_c = p<Scalar(0) ? Scalar(-1) : Scalar(1);
243 else if(numext::is_exactly_zero(p))
246 m_s = q<Scalar(0) ? Scalar(1) : Scalar(-1);
252 Scalar u =
sqrt(Scalar(1) + numext::abs2(t));
262 Scalar u =
sqrt(Scalar(1) + numext::abs2(t));
283 template<
typename VectorX,
typename VectorY,
typename OtherScalar>
285 void apply_rotation_in_the_plane(DenseBase<VectorX>& xpr_x, DenseBase<VectorY>& xpr_y,
const JacobiRotation<OtherScalar>& j);
294 template<
typename Derived>
295 template<
typename OtherScalar>
299 RowXpr x(this->row(p));
300 RowXpr y(this->row(q));
301 internal::apply_rotation_in_the_plane(x, y, j);
310 template<
typename Derived>
311 template<
typename OtherScalar>
315 ColXpr x(this->col(p));
316 ColXpr y(this->col(q));
317 internal::apply_rotation_in_the_plane(x, y, j.
transpose());
322 template<
typename Scalar,
typename OtherScalar,
323 int SizeAtCompileTime,
int MinAlignment,
bool Vectorizable>
324 struct apply_rotation_in_the_plane_selector
326 static EIGEN_DEVICE_FUNC
327 inline void run(Scalar *x,
Index incrx, Scalar *y,
Index incry,
Index size, OtherScalar c, OtherScalar s)
329 for(
Index i=0; i<size; ++i)
333 *x = c * xi + numext::conj(s) * yi;
334 *y = -s * xi + numext::conj(c) * yi;
341 template<
typename Scalar,
typename OtherScalar,
342 int SizeAtCompileTime,
int MinAlignment>
343 struct apply_rotation_in_the_plane_selector<Scalar,OtherScalar,SizeAtCompileTime,MinAlignment,true >
345 static inline void run(Scalar *x,
Index incrx, Scalar *y,
Index incry,
Index size, OtherScalar c, OtherScalar s)
348 PacketSize = packet_traits<Scalar>::size,
349 OtherPacketSize = packet_traits<OtherScalar>::size
351 typedef typename packet_traits<Scalar>::type Packet;
352 typedef typename packet_traits<OtherScalar>::type OtherPacket;
355 if(SizeAtCompileTime ==
Dynamic && ((incrx==1 && incry==1) || PacketSize == 1))
358 enum { Peeling = 2 };
360 Index alignedStart = internal::first_default_aligned(y, size);
361 Index alignedEnd = alignedStart + ((size-alignedStart)/PacketSize)*PacketSize;
363 const OtherPacket pc = pset1<OtherPacket>(c);
364 const OtherPacket ps = pset1<OtherPacket>(s);
365 conj_helper<OtherPacket,Packet,NumTraits<OtherScalar>::IsComplex,
false> pcj;
366 conj_helper<OtherPacket,Packet,false,false> pm;
368 for(
Index i=0; i<alignedStart; ++i)
372 x[i] = c * xi + numext::conj(s) * yi;
373 y[i] = -s * xi + numext::conj(c) * yi;
376 Scalar* EIGEN_RESTRICT px = x + alignedStart;
377 Scalar* EIGEN_RESTRICT py = y + alignedStart;
379 if(internal::first_default_aligned(x, size)==alignedStart)
381 for(
Index i=alignedStart; i<alignedEnd; i+=PacketSize)
383 Packet xi = pload<Packet>(px);
384 Packet yi = pload<Packet>(py);
385 pstore(px, padd(pm.pmul(pc,xi),pcj.pmul(ps,yi)));
386 pstore(py, psub(pcj.pmul(pc,yi),pm.pmul(ps,xi)));
393 Index peelingEnd = alignedStart + ((size-alignedStart)/(Peeling*PacketSize))*(Peeling*PacketSize);
394 for(
Index i=alignedStart; i<peelingEnd; i+=Peeling*PacketSize)
396 Packet xi = ploadu<Packet>(px);
397 Packet xi1 = ploadu<Packet>(px+PacketSize);
398 Packet yi = pload <Packet>(py);
399 Packet yi1 = pload <Packet>(py+PacketSize);
400 pstoreu(px, padd(pm.pmul(pc,xi),pcj.pmul(ps,yi)));
401 pstoreu(px+PacketSize, padd(pm.pmul(pc,xi1),pcj.pmul(ps,yi1)));
402 pstore (py, psub(pcj.pmul(pc,yi),pm.pmul(ps,xi)));
403 pstore (py+PacketSize, psub(pcj.pmul(pc,yi1),pm.pmul(ps,xi1)));
404 px += Peeling*PacketSize;
405 py += Peeling*PacketSize;
407 if(alignedEnd!=peelingEnd)
409 Packet xi = ploadu<Packet>(x+peelingEnd);
410 Packet yi = pload <Packet>(y+peelingEnd);
411 pstoreu(x+peelingEnd, padd(pm.pmul(pc,xi),pcj.pmul(ps,yi)));
412 pstore (y+peelingEnd, psub(pcj.pmul(pc,yi),pm.pmul(ps,xi)));
416 for(
Index i=alignedEnd; i<size; ++i)
420 x[i] = c * xi + numext::conj(s) * yi;
421 y[i] = -s * xi + numext::conj(c) * yi;
426 else if(SizeAtCompileTime !=
Dynamic && MinAlignment>0)
428 const OtherPacket pc = pset1<OtherPacket>(c);
429 const OtherPacket ps = pset1<OtherPacket>(s);
430 conj_helper<OtherPacket,Packet,NumTraits<OtherPacket>::IsComplex,
false> pcj;
431 conj_helper<OtherPacket,Packet,false,false> pm;
432 Scalar* EIGEN_RESTRICT px = x;
433 Scalar* EIGEN_RESTRICT py = y;
434 for(
Index i=0; i<size; i+=PacketSize)
436 Packet xi = pload<Packet>(px);
437 Packet yi = pload<Packet>(py);
438 pstore(px, padd(pm.pmul(pc,xi),pcj.pmul(ps,yi)));
439 pstore(py, psub(pcj.pmul(pc,yi),pm.pmul(ps,xi)));
448 apply_rotation_in_the_plane_selector<Scalar,OtherScalar,SizeAtCompileTime,MinAlignment,false>::run(x,incrx,y,incry,size,c,s);
453 template<
typename VectorX,
typename VectorY,
typename OtherScalar>
455 void apply_rotation_in_the_plane(DenseBase<VectorX>& xpr_x, DenseBase<VectorY>& xpr_y,
const JacobiRotation<OtherScalar>& j)
457 typedef typename VectorX::Scalar Scalar;
459 && (
int(packet_traits<Scalar>::size) == int(packet_traits<OtherScalar>::size));
461 eigen_assert(xpr_x.size() == xpr_y.size());
462 Index size = xpr_x.size();
463 Index incrx = xpr_x.derived().innerStride();
464 Index incry = xpr_y.derived().innerStride();
466 Scalar* EIGEN_RESTRICT x = &xpr_x.derived().coeffRef(0);
467 Scalar* EIGEN_RESTRICT y = &xpr_y.derived().coeffRef(0);
469 OtherScalar c = j.c();
470 OtherScalar s = j.s();
471 if (numext::is_exactly_one(c) && numext::is_exactly_zero(s))
474 apply_rotation_in_the_plane_selector<
477 plain_enum_min(evaluator<VectorX>::Alignment, evaluator<VectorY>::Alignment),
478 Vectorizable>::run(x,incrx,y,incry,size,c,s);
@ SizeAtCompileTime
Definition: DenseBase.h:108
@ Flags
Definition: DenseBase.h:159
CoeffReturnType coeff(Index row, Index col) const
Definition: DenseCoeffsBase.h:99
Rotation given by a cosine-sine pair.
Definition: Jacobi.h:37
JacobiRotation()
Definition: Jacobi.h:43
JacobiRotation(const Scalar &c, const Scalar &s)
Definition: Jacobi.h:47
bool makeJacobi(const MatrixBase< Derived > &, Index p, Index q)
Definition: Jacobi.h:141
JacobiRotation adjoint() const
Definition: Jacobi.h:69
JacobiRotation transpose() const
Definition: Jacobi.h:65
JacobiRotation operator*(const JacobiRotation &other)
Definition: Jacobi.h:56
void makeGivens(const Scalar &p, const Scalar &q, Scalar *r=0)
Definition: Jacobi.h:164
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
void applyOnTheLeft(const EigenBase< OtherDerived > &other)
Definition: MatrixBase.h:545
void applyOnTheRight(const EigenBase< OtherDerived > &other)
Definition: MatrixBase.h:533
const unsigned int PacketAccessBit
Definition: Constants.h:96
Namespace containing all symbols from the Eigen library.
Definition: Core:139
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:59
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_conjugate_op< typename Derived::Scalar >, const Derived > conj(const Eigen::ArrayBase< Derived > &x)
const int Dynamic
Definition: Constants.h:24
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:41
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:231