TR-mbed 1.0
Loading...
Searching...
No Matches
TriangularMatrixVector.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 Gael Guennebaud <gael.guennebaud@inria.fr>
5//
6// This Source Code Form is subject to the terms of the Mozilla
7// Public License v. 2.0. If a copy of the MPL was not distributed
8// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9
10#ifndef EIGEN_TRIANGULARMATRIXVECTOR_H
11#define EIGEN_TRIANGULARMATRIXVECTOR_H
12
13namespace Eigen {
14
15namespace internal {
16
17template<typename Index, int Mode, typename LhsScalar, bool ConjLhs, typename RhsScalar, bool ConjRhs, int StorageOrder, int Version=Specialized>
19
20template<typename Index, int Mode, typename LhsScalar, bool ConjLhs, typename RhsScalar, bool ConjRhs, int Version>
21struct triangular_matrix_vector_product<Index,Mode,LhsScalar,ConjLhs,RhsScalar,ConjRhs,ColMajor,Version>
22{
24 enum {
25 IsLower = ((Mode&Lower)==Lower),
26 HasUnitDiag = (Mode & UnitDiag)==UnitDiag,
27 HasZeroDiag = (Mode & ZeroDiag)==ZeroDiag
28 };
29 static EIGEN_DONT_INLINE void run(Index _rows, Index _cols, const LhsScalar* _lhs, Index lhsStride,
30 const RhsScalar* _rhs, Index rhsIncr, ResScalar* _res, Index resIncr, const RhsScalar& alpha);
31};
32
33template<typename Index, int Mode, typename LhsScalar, bool ConjLhs, typename RhsScalar, bool ConjRhs, int Version>
34EIGEN_DONT_INLINE void triangular_matrix_vector_product<Index,Mode,LhsScalar,ConjLhs,RhsScalar,ConjRhs,ColMajor,Version>
35 ::run(Index _rows, Index _cols, const LhsScalar* _lhs, Index lhsStride,
36 const RhsScalar* _rhs, Index rhsIncr, ResScalar* _res, Index resIncr, const RhsScalar& alpha)
37 {
39 Index size = (std::min)(_rows,_cols);
40 Index rows = IsLower ? _rows : (std::min)(_rows,_cols);
41 Index cols = IsLower ? (std::min)(_rows,_cols) : _cols;
42
44 const LhsMap lhs(_lhs,rows,cols,OuterStride<>(lhsStride));
46
48 const RhsMap rhs(_rhs,cols,InnerStride<>(rhsIncr));
50
52 ResMap res(_res,rows);
53
56
57 for (Index pi=0; pi<size; pi+=PanelWidth)
58 {
60 for (Index k=0; k<actualPanelWidth; ++k)
61 {
62 Index i = pi + k;
63 Index s = IsLower ? ((HasUnitDiag||HasZeroDiag) ? i+1 : i ) : pi;
64 Index r = IsLower ? actualPanelWidth-k : k+1;
65 if ((!(HasUnitDiag||HasZeroDiag)) || (--r)>0)
66 res.segment(s,r) += (alpha * cjRhs.coeff(i)) * cjLhs.col(i).segment(s,r);
67 if (HasUnitDiag)
68 res.coeffRef(i) += alpha * cjRhs.coeff(i);
69 }
70 Index r = IsLower ? rows - pi - actualPanelWidth : pi;
71 if (r>0)
72 {
73 Index s = IsLower ? pi+actualPanelWidth : 0;
76 LhsMapper(&lhs.coeffRef(s,pi), lhsStride),
77 RhsMapper(&rhs.coeffRef(pi), rhsIncr),
78 &res.coeffRef(s), resIncr, alpha);
79 }
80 }
81 if((!IsLower) && cols>size)
82 {
84 rows, cols-size,
85 LhsMapper(&lhs.coeffRef(0,size), lhsStride),
86 RhsMapper(&rhs.coeffRef(size), rhsIncr),
87 _res, resIncr, alpha);
88 }
89 }
90
91template<typename Index, int Mode, typename LhsScalar, bool ConjLhs, typename RhsScalar, bool ConjRhs,int Version>
92struct triangular_matrix_vector_product<Index,Mode,LhsScalar,ConjLhs,RhsScalar,ConjRhs,RowMajor,Version>
93{
95 enum {
96 IsLower = ((Mode&Lower)==Lower),
97 HasUnitDiag = (Mode & UnitDiag)==UnitDiag,
98 HasZeroDiag = (Mode & ZeroDiag)==ZeroDiag
99 };
100 static EIGEN_DONT_INLINE void run(Index _rows, Index _cols, const LhsScalar* _lhs, Index lhsStride,
101 const RhsScalar* _rhs, Index rhsIncr, ResScalar* _res, Index resIncr, const ResScalar& alpha);
102};
103
104template<typename Index, int Mode, typename LhsScalar, bool ConjLhs, typename RhsScalar, bool ConjRhs,int Version>
105EIGEN_DONT_INLINE void triangular_matrix_vector_product<Index,Mode,LhsScalar,ConjLhs,RhsScalar,ConjRhs,RowMajor,Version>
106 ::run(Index _rows, Index _cols, const LhsScalar* _lhs, Index lhsStride,
107 const RhsScalar* _rhs, Index rhsIncr, ResScalar* _res, Index resIncr, const ResScalar& alpha)
108 {
110 Index diagSize = (std::min)(_rows,_cols);
111 Index rows = IsLower ? _rows : diagSize;
112 Index cols = IsLower ? diagSize : _cols;
113
115 const LhsMap lhs(_lhs,rows,cols,OuterStride<>(lhsStride));
117
119 const RhsMap rhs(_rhs,cols);
121
124
127
128 for (Index pi=0; pi<diagSize; pi+=PanelWidth)
129 {
130 Index actualPanelWidth = (std::min)(PanelWidth, diagSize-pi);
131 for (Index k=0; k<actualPanelWidth; ++k)
132 {
133 Index i = pi + k;
134 Index s = IsLower ? pi : ((HasUnitDiag||HasZeroDiag) ? i+1 : i);
135 Index r = IsLower ? k+1 : actualPanelWidth-k;
136 if ((!(HasUnitDiag||HasZeroDiag)) || (--r)>0)
137 res.coeffRef(i) += alpha * (cjLhs.row(i).segment(s,r).cwiseProduct(cjRhs.segment(s,r).transpose())).sum();
138 if (HasUnitDiag)
139 res.coeffRef(i) += alpha * cjRhs.coeff(i);
140 }
141 Index r = IsLower ? pi : cols - pi - actualPanelWidth;
142 if (r>0)
143 {
144 Index s = IsLower ? 0 : pi + actualPanelWidth;
147 LhsMapper(&lhs.coeffRef(pi,s), lhsStride),
148 RhsMapper(&rhs.coeffRef(s), rhsIncr),
149 &res.coeffRef(pi), resIncr, alpha);
150 }
151 }
152 if(IsLower && rows>diagSize)
153 {
155 rows-diagSize, cols,
156 LhsMapper(&lhs.coeffRef(diagSize,0), lhsStride),
157 RhsMapper(&rhs.coeffRef(0), rhsIncr),
158 &res.coeffRef(diagSize), resIncr, alpha);
159 }
160 }
161
162/***************************************************************************
163* Wrapper to product_triangular_vector
164***************************************************************************/
165
166template<int Mode,int StorageOrder>
168
169} // end namespace internal
170
171namespace internal {
172
173template<int Mode, typename Lhs, typename Rhs>
175{
176 template<typename Dest> static void run(Dest& dst, const Lhs &lhs, const Rhs &rhs, const typename Dest::Scalar& alpha)
177 {
178 eigen_assert(dst.rows()==lhs.rows() && dst.cols()==rhs.cols());
179
181 }
182};
183
184template<int Mode, typename Lhs, typename Rhs>
186{
187 template<typename Dest> static void run(Dest& dst, const Lhs &lhs, const Rhs &rhs, const typename Dest::Scalar& alpha)
188 {
189 eigen_assert(dst.rows()==lhs.rows() && dst.cols()==rhs.cols());
190
192 internal::trmv_selector<(Mode & (UnitDiag|ZeroDiag)) | ((Mode & Lower) ? Upper : Lower),
194 ::run(rhs.transpose(),lhs.transpose(), dstT, alpha);
195 }
196};
197
198} // end namespace internal
199
200namespace internal {
201
202// TODO: find a way to factorize this piece of code with gemv_selector since the logic is exactly the same.
203
204template<int Mode> struct trmv_selector<Mode,ColMajor>
205{
206 template<typename Lhs, typename Rhs, typename Dest>
207 static void run(const Lhs &lhs, const Rhs &rhs, Dest& dest, const typename Dest::Scalar& alpha)
208 {
209 typedef typename Lhs::Scalar LhsScalar;
210 typedef typename Rhs::Scalar RhsScalar;
211 typedef typename Dest::Scalar ResScalar;
212 typedef typename Dest::RealScalar RealScalar;
213
214 typedef internal::blas_traits<Lhs> LhsBlasTraits;
215 typedef typename LhsBlasTraits::DirectLinearAccessType ActualLhsType;
216 typedef internal::blas_traits<Rhs> RhsBlasTraits;
217 typedef typename RhsBlasTraits::DirectLinearAccessType ActualRhsType;
218
220
221 typename internal::add_const_on_value_type<ActualLhsType>::type actualLhs = LhsBlasTraits::extract(lhs);
222 typename internal::add_const_on_value_type<ActualRhsType>::type actualRhs = RhsBlasTraits::extract(rhs);
223
224 LhsScalar lhs_alpha = LhsBlasTraits::extractScalarFactor(lhs);
225 RhsScalar rhs_alpha = RhsBlasTraits::extractScalarFactor(rhs);
226 ResScalar actualAlpha = alpha * lhs_alpha * rhs_alpha;
227
228 enum {
229 // FIXME find a way to allow an inner stride on the result if packet_traits<Scalar>::size==1
230 // on, the other hand it is good for the cache to pack the vector anyways...
231 EvalToDestAtCompileTime = Dest::InnerStrideAtCompileTime==1,
233 MightCannotUseDest = (Dest::InnerStrideAtCompileTime!=1) || ComplexByReal
234 };
235
237
238 bool alphaIsCompatible = (!ComplexByReal) || (numext::imag(actualAlpha)==RealScalar(0));
240
242
244 evalToDest ? dest.data() : static_dest.data());
245
246 if(!evalToDest)
247 {
248 #ifdef EIGEN_DENSE_STORAGE_CTOR_PLUGIN
249 Index size = dest.size();
251 #endif
253 {
254 MappedDest(actualDestPtr, dest.size()).setZero();
255 compatibleAlpha = RhsScalar(1);
256 }
257 else
259 }
260
262 <Index,Mode,
263 LhsScalar, LhsBlasTraits::NeedToConjugate,
264 RhsScalar, RhsBlasTraits::NeedToConjugate,
265 ColMajor>
266 ::run(actualLhs.rows(),actualLhs.cols(),
267 actualLhs.data(),actualLhs.outerStride(),
268 actualRhs.data(),actualRhs.innerStride(),
270
271 if (!evalToDest)
272 {
275 else
277 }
278
279 if ( ((Mode&UnitDiag)==UnitDiag) && (lhs_alpha!=LhsScalar(1)) )
280 {
281 Index diagSize = (std::min)(lhs.rows(),lhs.cols());
282 dest.head(diagSize) -= (lhs_alpha-LhsScalar(1))*rhs.head(diagSize);
283 }
284 }
285};
286
287template<int Mode> struct trmv_selector<Mode,RowMajor>
288{
289 template<typename Lhs, typename Rhs, typename Dest>
290 static void run(const Lhs &lhs, const Rhs &rhs, Dest& dest, const typename Dest::Scalar& alpha)
291 {
292 typedef typename Lhs::Scalar LhsScalar;
293 typedef typename Rhs::Scalar RhsScalar;
294 typedef typename Dest::Scalar ResScalar;
295
296 typedef internal::blas_traits<Lhs> LhsBlasTraits;
297 typedef typename LhsBlasTraits::DirectLinearAccessType ActualLhsType;
298 typedef internal::blas_traits<Rhs> RhsBlasTraits;
299 typedef typename RhsBlasTraits::DirectLinearAccessType ActualRhsType;
300 typedef typename internal::remove_all<ActualRhsType>::type ActualRhsTypeCleaned;
301
302 typename add_const<ActualLhsType>::type actualLhs = LhsBlasTraits::extract(lhs);
303 typename add_const<ActualRhsType>::type actualRhs = RhsBlasTraits::extract(rhs);
304
305 LhsScalar lhs_alpha = LhsBlasTraits::extractScalarFactor(lhs);
306 RhsScalar rhs_alpha = RhsBlasTraits::extractScalarFactor(rhs);
307 ResScalar actualAlpha = alpha * lhs_alpha * rhs_alpha;
308
309 enum {
310 DirectlyUseRhs = ActualRhsTypeCleaned::InnerStrideAtCompileTime==1
311 };
312
314
316 DirectlyUseRhs ? const_cast<RhsScalar*>(actualRhs.data()) : static_rhs.data());
317
318 if(!DirectlyUseRhs)
319 {
320 #ifdef EIGEN_DENSE_STORAGE_CTOR_PLUGIN
321 Index size = actualRhs.size();
323 #endif
325 }
326
328 <Index,Mode,
329 LhsScalar, LhsBlasTraits::NeedToConjugate,
330 RhsScalar, RhsBlasTraits::NeedToConjugate,
331 RowMajor>
332 ::run(actualLhs.rows(),actualLhs.cols(),
333 actualLhs.data(),actualLhs.outerStride(),
334 actualRhsPtr,1,
335 dest.data(),dest.innerStride(),
337
338 if ( ((Mode&UnitDiag)==UnitDiag) && (lhs_alpha!=LhsScalar(1)) )
339 {
340 Index diagSize = (std::min)(lhs.rows(),lhs.cols());
341 dest.head(diagSize) -= (lhs_alpha-LhsScalar(1))*rhs.head(diagSize);
342 }
343 }
344};
345
346} // end namespace internal
347
348} // end namespace Eigen
349
350#endif // EIGEN_TRIANGULARMATRIXVECTOR_H
int i
Definition BiCGSTAB_step_by_step.cpp:9
#define EIGEN_PLAIN_ENUM_MIN(a, b)
Definition Macros.h:1288
#define EIGEN_DONT_INLINE
Definition Macros.h:940
#define eigen_assert(x)
Definition Macros.h:1037
#define ei_declare_aligned_stack_constructed_variable(TYPE, NAME, SIZE, BUFFER)
Definition Memory.h:768
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
Definition PartialRedux_count.cpp:3
#define EIGEN_TUNE_TRIANGULAR_PANEL_WIDTH
Definition Settings.h:38
int rows
Definition Tutorial_commainit_02.cpp:1
int cols
Definition Tutorial_commainit_02.cpp:1
Scalar Scalar int size
Definition benchVecAdd.cpp:17
NumTraits< Scalar >::Real RealScalar
Definition bench_gemm.cpp:47
@ UnitDiag
Definition Constants.h:213
@ ZeroDiag
Definition Constants.h:215
@ Lower
Definition Constants.h:209
@ Upper
Definition Constants.h:211
@ AlignedMax
Definition Constants.h:252
@ ColMajor
Definition Constants.h:319
@ RowMajor
Definition Constants.h:321
const unsigned int RowMajorBit
Definition Constants.h:66
RealScalar s
Definition level1_cplx_impl.h:126
return int(ret)+1
RealScalar alpha
Definition level1_cplx_impl.h:147
@ Lhs
Definition TensorContractionMapper.h:19
@ Rhs
Definition TensorContractionMapper.h:18
Namespace containing all symbols from the Eigen library.
Definition bench_norm.cpp:85
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:74
Definition BandTriangularSolver.h:13
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition NumTraits.h:233
Determines whether the given binary operation of two numeric types is allowed and what the scalar ret...
Definition XprHelper.h:806
static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE To run(const From &x)
Definition BlasUtil.h:43
Definition GenericPacketMath.h:107
Definition ForwardDeclarations.h:17
ScalarBinaryOpTraits< LhsScalar, RhsScalar >::ReturnType ResScalar
Definition TriangularMatrixVector.h:94
ScalarBinaryOpTraits< LhsScalar, RhsScalar >::ReturnType ResScalar
Definition TriangularMatrixVector.h:23
Definition TriangularMatrixVector.h:18
static void run(Dest &dst, const Lhs &lhs, const Rhs &rhs, const typename Dest::Scalar &alpha)
Definition TriangularMatrixVector.h:187
static void run(Dest &dst, const Lhs &lhs, const Rhs &rhs, const typename Dest::Scalar &alpha)
Definition TriangularMatrixVector.h:176
Definition ProductEvaluators.h:758
static void run(const Lhs &lhs, const Rhs &rhs, Dest &dest, const typename Dest::Scalar &alpha)
Definition TriangularMatrixVector.h:207
static void run(const Lhs &lhs, const Rhs &rhs, Dest &dest, const typename Dest::Scalar &alpha)
Definition TriangularMatrixVector.h:290
Definition TriangularMatrixVector.h:167