10#ifndef EIGEN_SELFADJOINT_MATRIX_VECTOR_H
11#define EIGEN_SELFADJOINT_MATRIX_VECTOR_H
23template<
typename Scalar,
typename Index,
int StorageOrder,
int UpLo,
bool ConjugateLhs,
bool ConjugateRhs,
int Version=Specialized>
24struct selfadjoint_matrix_vector_product;
26template<
typename Scalar,
typename Index,
int StorageOrder,
int UpLo,
bool ConjugateLhs,
bool ConjugateRhs,
int Version>
27struct selfadjoint_matrix_vector_product
30static EIGEN_DONT_INLINE EIGEN_DEVICE_FUNC
33 const Scalar* lhs,
Index lhsStride,
39template<
typename Scalar,
typename Index,
int StorageOrder,
int UpLo,
bool ConjugateLhs,
bool ConjugateRhs,
int Version>
40EIGEN_DONT_INLINE EIGEN_DEVICE_FUNC
41void selfadjoint_matrix_vector_product<Scalar,Index,StorageOrder,UpLo,ConjugateLhs,ConjugateRhs,Version>::run(
43 const Scalar* lhs,
Index lhsStride,
48 typedef typename packet_traits<Scalar>::type Packet;
49 typedef typename NumTraits<Scalar>::Real RealScalar;
50 const Index PacketSize =
sizeof(Packet)/
sizeof(Scalar);
53 IsRowMajor = StorageOrder==
RowMajor ? 1 : 0,
54 IsLower = UpLo ==
Lower ? 1 : 0,
55 FirstTriangular = IsRowMajor == IsLower
58 conj_helper<Scalar,Scalar,NumTraits<Scalar>::IsComplex && EIGEN_LOGICAL_XOR(ConjugateLhs, IsRowMajor), ConjugateRhs> cj0;
59 conj_helper<Scalar,Scalar,NumTraits<Scalar>::IsComplex && EIGEN_LOGICAL_XOR(ConjugateLhs, !IsRowMajor), ConjugateRhs> cj1;
60 conj_helper<RealScalar,Scalar,false, ConjugateRhs> cjd;
62 conj_helper<Packet,Packet,NumTraits<Scalar>::IsComplex && EIGEN_LOGICAL_XOR(ConjugateLhs, IsRowMajor), ConjugateRhs> pcj0;
63 conj_helper<Packet,Packet,NumTraits<Scalar>::IsComplex && EIGEN_LOGICAL_XOR(ConjugateLhs, !IsRowMajor), ConjugateRhs> pcj1;
65 Scalar cjAlpha = ConjugateRhs ? numext::conj(alpha) : alpha;
67 Index bound = numext::maxi(
Index(0), size-8) & 0xfffffffe;
71 for (
Index j=FirstTriangular ? bound : 0;
72 j<(FirstTriangular ? size : bound);j+=2)
74 const Scalar* EIGEN_RESTRICT A0 = lhs + j*lhsStride;
75 const Scalar* EIGEN_RESTRICT A1 = lhs + (j+1)*lhsStride;
77 Scalar t0 = cjAlpha * rhs[j];
78 Packet ptmp0 = pset1<Packet>(t0);
79 Scalar t1 = cjAlpha * rhs[j+1];
80 Packet ptmp1 = pset1<Packet>(t1);
83 Packet ptmp2 = pset1<Packet>(t2);
85 Packet ptmp3 = pset1<Packet>(t3);
87 Index starti = FirstTriangular ? 0 : j+2;
88 Index endi = FirstTriangular ? j : size;
89 Index alignedStart = (starti) + internal::first_default_aligned(&res[starti], endi-starti);
90 Index alignedEnd = alignedStart + ((endi-alignedStart)/(PacketSize))*(PacketSize);
92 res[j] += cjd.pmul(numext::real(A0[j]), t0);
93 res[j+1] += cjd.pmul(numext::real(A1[j+1]), t1);
96 res[j] += cj0.pmul(A1[j], t1);
97 t3 += cj1.pmul(A1[j], rhs[j]);
101 res[j+1] += cj0.pmul(A0[j+1],t0);
102 t2 += cj1.pmul(A0[j+1], rhs[j+1]);
105 for (
Index i=starti; i<alignedStart; ++i)
107 res[i] += cj0.pmul(A0[i], t0) + cj0.pmul(A1[i],t1);
108 t2 += cj1.pmul(A0[i], rhs[i]);
109 t3 += cj1.pmul(A1[i], rhs[i]);
113 const Scalar* EIGEN_RESTRICT a0It = A0 + alignedStart;
114 const Scalar* EIGEN_RESTRICT a1It = A1 + alignedStart;
115 const Scalar* EIGEN_RESTRICT rhsIt = rhs + alignedStart;
116 Scalar* EIGEN_RESTRICT resIt = res + alignedStart;
117 for (
Index i=alignedStart; i<alignedEnd; i+=PacketSize)
119 Packet A0i = ploadu<Packet>(a0It); a0It += PacketSize;
120 Packet A1i = ploadu<Packet>(a1It); a1It += PacketSize;
121 Packet Bi = ploadu<Packet>(rhsIt); rhsIt += PacketSize;
122 Packet Xi = pload <Packet>(resIt);
124 Xi = pcj0.pmadd(A0i,ptmp0, pcj0.pmadd(A1i,ptmp1,Xi));
125 ptmp2 = pcj1.pmadd(A0i, Bi, ptmp2);
126 ptmp3 = pcj1.pmadd(A1i, Bi, ptmp3);
127 pstore(resIt,Xi); resIt += PacketSize;
129 for (
Index i=alignedEnd; i<endi; i++)
131 res[i] += cj0.pmul(A0[i], t0) + cj0.pmul(A1[i],t1);
132 t2 += cj1.pmul(A0[i], rhs[i]);
133 t3 += cj1.pmul(A1[i], rhs[i]);
136 res[j] += alpha * (t2 + predux(ptmp2));
137 res[j+1] += alpha * (t3 + predux(ptmp3));
139 for (
Index j=FirstTriangular ? 0 : bound;j<(FirstTriangular ? bound : size);j++)
141 const Scalar* EIGEN_RESTRICT A0 = lhs + j*lhsStride;
143 Scalar t1 = cjAlpha * rhs[j];
145 res[j] += cjd.pmul(numext::real(A0[j]), t1);
146 for (
Index i=FirstTriangular ? 0 : j+1; i<(FirstTriangular ? j : size); i++)
148 res[i] += cj0.pmul(A0[i], t1);
149 t2 += cj1.pmul(A0[i], rhs[i]);
151 res[j] += alpha * t2;
163template<
typename Lhs,
int LhsMode,
typename Rhs>
164struct selfadjoint_product_impl<Lhs,LhsMode,false,Rhs,0,true>
166 typedef typename Product<Lhs,Rhs>::Scalar Scalar;
168 typedef internal::blas_traits<Lhs> LhsBlasTraits;
169 typedef typename LhsBlasTraits::DirectLinearAccessType ActualLhsType;
170 typedef typename internal::remove_all<ActualLhsType>::type ActualLhsTypeCleaned;
172 typedef internal::blas_traits<Rhs> RhsBlasTraits;
173 typedef typename RhsBlasTraits::DirectLinearAccessType ActualRhsType;
174 typedef typename internal::remove_all<ActualRhsType>::type ActualRhsTypeCleaned;
178 template<
typename Dest>
179 static EIGEN_DEVICE_FUNC
180 void run(Dest& dest,
const Lhs &a_lhs,
const Rhs &a_rhs,
const Scalar& alpha)
182 typedef typename Dest::Scalar ResScalar;
183 typedef typename Rhs::Scalar RhsScalar;
184 typedef Map<Matrix<ResScalar,Dynamic,1>, EIGEN_PLAIN_ENUM_MIN(AlignedMax,internal::packet_traits<ResScalar>::size)> MappedDest;
186 eigen_assert(dest.rows()==a_lhs.rows() && dest.cols()==a_rhs.cols());
188 typename internal::add_const_on_value_type<ActualLhsType>::type lhs = LhsBlasTraits::extract(a_lhs);
189 typename internal::add_const_on_value_type<ActualRhsType>::type rhs = RhsBlasTraits::extract(a_rhs);
191 Scalar actualAlpha = alpha * LhsBlasTraits::extractScalarFactor(a_lhs)
192 * RhsBlasTraits::extractScalarFactor(a_rhs);
195 EvalToDest = (Dest::InnerStrideAtCompileTime==1),
196 UseRhs = (ActualRhsTypeCleaned::InnerStrideAtCompileTime==1)
199 internal::gemv_static_vector_if<ResScalar,Dest::SizeAtCompileTime,Dest::MaxSizeAtCompileTime,!EvalToDest> static_dest;
200 internal::gemv_static_vector_if<RhsScalar,ActualRhsTypeCleaned::SizeAtCompileTime,ActualRhsTypeCleaned::MaxSizeAtCompileTime,!UseRhs> static_rhs;
202 ei_declare_aligned_stack_constructed_variable(ResScalar,actualDestPtr,dest.size(),
203 EvalToDest ? dest.data() : static_dest.data());
205 ei_declare_aligned_stack_constructed_variable(RhsScalar,actualRhsPtr,rhs.size(),
206 UseRhs ?
const_cast<RhsScalar*
>(rhs.data()) : static_rhs.data());
210 #ifdef EIGEN_DENSE_STORAGE_CTOR_PLUGIN
211 Index size = dest.size();
212 EIGEN_DENSE_STORAGE_CTOR_PLUGIN
214 MappedDest(actualDestPtr, dest.size()) = dest;
219 #ifdef EIGEN_DENSE_STORAGE_CTOR_PLUGIN
220 Index size = rhs.size();
221 EIGEN_DENSE_STORAGE_CTOR_PLUGIN
223 Map<typename ActualRhsTypeCleaned::PlainObject>(actualRhsPtr, rhs.size()) = rhs;
228 int(LhsUpLo), bool(LhsBlasTraits::NeedToConjugate), bool(RhsBlasTraits::NeedToConjugate)>::run
231 &lhs.coeffRef(0,0), lhs.outerStride(),
238 dest = MappedDest(actualDestPtr, dest.size());
242template<
typename Lhs,
typename Rhs,
int RhsMode>
243struct selfadjoint_product_impl<Lhs,0,true,Rhs,RhsMode,false>
245 typedef typename Product<Lhs,Rhs>::Scalar Scalar;
248 template<
typename Dest>
249 static void run(Dest& dest,
const Lhs &a_lhs,
const Rhs &a_rhs,
const Scalar& alpha)
252 Transpose<Dest> destT(dest);
253 selfadjoint_product_impl<Transpose<const Rhs>, int(RhsUpLo)==
Upper ?
Lower :
Upper,
false,
254 Transpose<const Lhs>, 0,
true>::run(destT, a_rhs.transpose(), a_lhs.transpose(), alpha);
@ Lower
Definition Constants.h:209
@ Upper
Definition Constants.h:211
@ ColMajor
Definition Constants.h:319
@ RowMajor
Definition Constants.h:321
const unsigned int RowMajorBit
Definition Constants.h:66
Namespace containing all symbols from the Eigen library.
Definition B01_Experimental.dox:1
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:74