15#include "./InternalHeaderCheck.h"
37template <
typename Scalar>
40 typedef typename NumTraits<Scalar>::Real RealScalar;
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; }
72 template <
typename Derived>
74 EIGEN_DEVICE_FUNC
bool makeJacobi(
const RealScalar& x,
const Scalar& y,
const RealScalar& z);
92template <
typename Scalar>
97 RealScalar deno = RealScalar(2) *
abs(y);
98 if (deno < (std::numeric_limits<RealScalar>::min)()) {
103 RealScalar tau = (x - z) / deno;
104 RealScalar w =
sqrt(numext::abs2(tau) + RealScalar(1));
106 if (tau > RealScalar(0)) {
107 t = RealScalar(1) / (tau + w);
109 t = RealScalar(1) / (tau - w);
111 RealScalar sign_t = t > RealScalar(0) ? RealScalar(1) : RealScalar(-1);
112 RealScalar n = RealScalar(1) /
sqrt(numext::abs2(t) + RealScalar(1));
113 m_s = -sign_t * (numext::conj(y) /
abs(y)) *
abs(t) * n;
129template <
typename Scalar>
130template <
typename Derived>
132 return makeJacobi(numext::real(m.coeff(p, p)), m.coeff(p, q), numext::real(m.coeff(q, q)));
151template <
typename Scalar>
153 makeGivens(p, q, r, std::conditional_t<NumTraits<Scalar>::IsComplex, internal::true_type, internal::false_type>());
157template <
typename Scalar>
159 internal::true_type) {
168 }
else if (p == Scalar(0)) {
173 RealScalar p1 = numext::norm1(p);
174 RealScalar q1 = numext::norm1(q);
177 RealScalar p2 = numext::abs2(ps);
179 RealScalar q2 = numext::abs2(qs);
181 RealScalar u =
sqrt(RealScalar(1) + q2 / p2);
182 if (numext::real(p) < RealScalar(0)) u = -u;
185 m_s = -qs *
conj(ps) * (m_c / p2);
189 RealScalar p2 = numext::abs2(ps);
191 RealScalar q2 = numext::abs2(qs);
193 RealScalar u = q1 *
sqrt(p2 + q2);
194 if (numext::real(p) < RealScalar(0)) u = -u;
199 m_s = -
conj(ps) * (q / u);
206template <
typename Scalar>
208 internal::false_type) {
211 if (numext::is_exactly_zero(q)) {
215 }
else if (numext::is_exactly_zero(p)) {
219 }
else if (
abs(p) >
abs(q)) {
222 if (p <
Scalar(0)) u = -u;
229 if (q <
Scalar(0)) u = -u;
248template <
typename VectorX,
typename VectorY,
typename OtherScalar>
249EIGEN_DEVICE_FUNC
void apply_rotation_in_the_plane(DenseBase<VectorX>& xpr_x, DenseBase<VectorY>& xpr_y,
250 const JacobiRotation<OtherScalar>& j);
259template <
typename Derived>
260template <
typename OtherScalar>
263 RowXpr x(this->
row(p));
264 RowXpr y(this->
row(q));
265 internal::apply_rotation_in_the_plane(x, y, j);
274template <
typename Derived>
275template <
typename OtherScalar>
278 ColXpr x(this->
col(p));
279 ColXpr y(this->
col(q));
280 internal::apply_rotation_in_the_plane(x, y, j.
transpose());
285template <
typename Scalar,
typename OtherScalar,
int SizeAtCompileTime,
int MinAlignment,
bool Vectorizable>
286struct apply_rotation_in_the_plane_selector {
289 for (
Index i = 0; i < size; ++i) {
292 *x = c * xi + numext::conj(s) * yi;
293 *y = -s * xi + numext::conj(c) * yi;
300template <
typename Scalar,
typename OtherScalar,
int SizeAtCompileTime,
int MinAlignment>
301struct apply_rotation_in_the_plane_selector<Scalar, OtherScalar, SizeAtCompileTime, MinAlignment,
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;
307 constexpr int RequiredAlignment =
308 (std::max)(unpacket_traits<Packet>::alignment, unpacket_traits<OtherPacket>::alignment);
309 constexpr Index PacketSize = packet_traits<Scalar>::size;
312 if (size >= 2 * PacketSize && SizeAtCompileTime ==
Dynamic && ((incrx == 1 && incry == 1) || PacketSize == 1)) {
314 constexpr Index Peeling = 2;
316 Index alignedStart = internal::first_default_aligned(y, size);
317 Index alignedEnd = alignedStart + ((size - alignedStart) / PacketSize) * PacketSize;
319 const OtherPacket pc = pset1<OtherPacket>(c);
320 const OtherPacket ps = pset1<OtherPacket>(s);
321 conj_helper<OtherPacket, Packet, NumTraits<OtherScalar>::IsComplex,
false> pcj;
322 conj_helper<OtherPacket, Packet, false, false> pm;
324 for (
Index i = 0; i < alignedStart; ++i) {
327 x[i] = c * xi + numext::conj(s) * yi;
328 y[i] = -s * xi + numext::conj(c) * yi;
331 Scalar* EIGEN_RESTRICT px = x + alignedStart;
332 Scalar* EIGEN_RESTRICT py = y + alignedStart;
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)));
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;
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)));
365 for (
Index i = alignedEnd; i < size; ++i) {
368 x[i] = c * xi + numext::conj(s) * yi;
369 y[i] = -s * xi + numext::conj(c) * yi;
374 else if (SizeAtCompileTime !=
Dynamic && MinAlignment >= RequiredAlignment) {
375 const OtherPacket pc = pset1<OtherPacket>(c);
376 const OtherPacket ps = pset1<OtherPacket>(s);
377 conj_helper<OtherPacket, Packet, NumTraits<OtherScalar>::IsComplex,
false> pcj;
378 conj_helper<OtherPacket, Packet, false, false> pm;
379 Scalar* EIGEN_RESTRICT px = x;
380 Scalar* EIGEN_RESTRICT py = y;
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)));
393 apply_rotation_in_the_plane_selector<Scalar, OtherScalar, SizeAtCompileTime, MinAlignment, false>::run(
394 x, incrx, y, incry, size, c, s);
399template <
typename VectorX,
typename VectorY,
typename OtherScalar>
400EIGEN_DEVICE_FUNC
void inline apply_rotation_in_the_plane(DenseBase<VectorX>& xpr_x, DenseBase<VectorY>& xpr_y,
401 const JacobiRotation<OtherScalar>& j) {
402 typedef typename VectorX::Scalar Scalar;
403 constexpr bool Vectorizable = (int(evaluator<VectorX>::Flags) & int(evaluator<VectorY>::Flags) &
PacketAccessBit) &&
404 (
int(packet_traits<Scalar>::size) ==
int(packet_traits<OtherScalar>::size));
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();
411 Scalar* EIGEN_RESTRICT x = &xpr_x.derived().coeffRef(0);
412 Scalar* EIGEN_RESTRICT y = &xpr_y.derived().coeffRef(0);
414 OtherScalar c = j.c();
415 OtherScalar s = j.s();
416 if (numext::is_exactly_one(c) && numext::is_exactly_zero(s))
return;
418 constexpr int Alignment = (std::min)(
int(evaluator<VectorX>::Alignment),
int(evaluator<VectorY>::Alignment));
419 apply_rotation_in_the_plane_selector<Scalar, OtherScalar, VectorX::SizeAtCompileTime, Alignment, Vectorizable>::run(
420 x, incrx, y, incry, size, c, s);
RowXpr row(Index i)
Definition DenseBase.h:1085
ColXpr col(Index i)
Definition DenseBase.h:1072
Rotation given by a cosine-sine pair.
Definition Jacobi.h:38
JacobiRotation()
Definition Jacobi.h:43
JacobiRotation(const Scalar &c, const Scalar &s)
Definition Jacobi.h:46
bool makeJacobi(const MatrixBase< Derived > &, Index p, Index q)
Definition Jacobi.h:131
JacobiRotation adjoint() const
Definition Jacobi.h:67
JacobiRotation transpose() const
Definition Jacobi.h:61
JacobiRotation operator*(const JacobiRotation &other)
Definition Jacobi.h:54
void makeGivens(const Scalar &p, const Scalar &q, Scalar *r=0)
Definition Jacobi.h:152
Base class for all dense matrices, vectors, and expressions.
Definition ForwardDeclarations.h:73
void applyOnTheLeft(const EigenBase< OtherDerived > &other)
Definition MatrixBase.h:532
void applyOnTheRight(const EigenBase< OtherDerived > &other)
Definition MatrixBase.h:521
const unsigned int PacketAccessBit
Definition Constants.h:97
Namespace containing all symbols from the Eigen library.
Definition B01_Experimental.dox:1
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_conjugate_op< typename Derived::Scalar >, const Derived > conj(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:82
const int Dynamic
Definition Constants.h:25