10 #ifndef EIGEN_SPARSE_SELFADJOINTVIEW_H 11 #define EIGEN_SPARSE_SELFADJOINTVIEW_H 31 template<
typename MatrixType,
unsigned int Mode>
35 template<
int SrcMode,
int DstMode,
typename MatrixType,
int DestOrder>
38 template<
int Mode,
typename MatrixType,
int DestOrder>
44 :
public EigenBase<SparseSelfAdjointView<MatrixType,_Mode> >
55 typedef typename MatrixType::Scalar Scalar;
56 typedef typename MatrixType::StorageIndex StorageIndex;
63 eigen_assert(rows()==cols() &&
"SelfAdjointView is only for squared matrices");
66 inline Index rows()
const {
return m_matrix.rows(); }
67 inline Index cols()
const {
return m_matrix.cols(); }
70 const _MatrixTypeNested& matrix()
const {
return m_matrix; }
78 template<
typename OtherDerived>
90 template<
typename OtherDerived>
friend 98 template<
typename OtherDerived>
106 template<
typename OtherDerived>
friend 121 template<
typename DerivedU>
131 template<
typename SrcMatrixType,
int SrcMode>
134 internal::call_assignment_no_alias_no_transpose(*
this, permutedMatrix);
144 template<
typename SrcMatrixType,
unsigned int SrcMode>
153 EIGEN_ONLY_USED_FOR_DEBUG(rows);
154 EIGEN_ONLY_USED_FOR_DEBUG(cols);
155 eigen_assert(rows == this->rows() && cols == this->cols()
156 &&
"SparseSelfadjointView::resize() does not actually allow to resize.");
161 MatrixTypeNested m_matrix;
165 template<
typename Dest>
void evalTo(Dest &)
const;
172 template<
typename Derived>
173 template<
unsigned int UpLo>
179 template<
typename Derived>
180 template<
unsigned int UpLo>
190 template<
typename MatrixType,
unsigned int Mode>
191 template<
typename DerivedU>
197 m_matrix = tmp.template triangularView<Mode>();
199 m_matrix += alpha * tmp.template triangularView<Mode>();
209 template<
typename MatrixType,
unsigned int Mode>
221 template<
typename DstXprType,
typename SrcXprType,
typename Functor>
224 typedef typename DstXprType::StorageIndex StorageIndex;
225 template<
typename DestScalar,
int StorageOrder>
228 internal::permute_symm_to_fullsymm<SrcXprType::Mode>(src.matrix(), dst);
231 template<
typename DestScalar>
236 internal::permute_symm_to_fullsymm<SrcXprType::Mode>(src.matrix(), tmp);
249 template<
int Mode,
typename SparseLhsType,
typename DenseRhsType,
typename DenseResType,
typename AlphaType>
250 inline void sparse_selfadjoint_time_dense_product(
const SparseLhsType& lhs,
const DenseRhsType& rhs, DenseResType& res,
const AlphaType& alpha)
252 EIGEN_ONLY_USED_FOR_DEBUG(alpha);
257 typedef typename LhsEval::InnerIterator LhsIterator;
258 typedef typename SparseLhsType::Scalar LhsScalar;
264 || ( (Mode&
Upper) && !LhsIsRowMajor)
265 || ( (Mode&
Lower) && LhsIsRowMajor),
266 ProcessSecondHalf = !ProcessFirstHalf
269 SparseLhsTypeNested lhs_nested(lhs);
270 LhsEval lhsEval(lhs_nested);
273 for (
Index k=0; k<rhs.cols(); ++k)
275 for (
Index j=0; j<lhs.outerSize(); ++j)
277 LhsIterator i(lhsEval,j);
279 if (ProcessSecondHalf)
281 while (i && i.index()<j) ++i;
282 if(i && i.index()==j)
284 res(j,k) += alpha * i.value() * rhs(j,k);
292 typename DenseResType::Scalar res_j(0);
293 for(; (ProcessFirstHalf ? i && i.index() < j : i) ; ++i)
295 LhsScalar lhs_ij = i.value();
296 if(!LhsIsRowMajor) lhs_ij = numext::conj(lhs_ij);
297 res_j += lhs_ij * rhs(i.index(),k);
298 res(i.index(),k) += numext::conj(lhs_ij) * rhs_j;
300 res(j,k) += alpha * res_j;
303 if (ProcessFirstHalf && i && (i.index()==j))
304 res(j,k) += alpha * i.value() * rhs(j,k);
310 template<
typename LhsView,
typename Rhs,
int ProductType>
312 :
generic_product_impl_base<LhsView, Rhs, generic_product_impl<LhsView, Rhs, SparseSelfAdjointShape, DenseShape, ProductType> >
314 template<
typename Dest>
315 static void scaleAndAddTo(Dest& dst,
const LhsView& lhsView,
const Rhs& rhs,
const typename Dest::Scalar& alpha)
317 typedef typename LhsView::_MatrixTypeNested Lhs;
320 LhsNested lhsNested(lhsView.matrix());
321 RhsNested rhsNested(rhs);
323 internal::sparse_selfadjoint_time_dense_product<LhsView::Mode>(lhsNested, rhsNested, dst, alpha);
327 template<
typename Lhs,
typename RhsView,
int ProductType>
329 :
generic_product_impl_base<Lhs, RhsView, generic_product_impl<Lhs, RhsView, DenseShape, SparseSelfAdjointShape, ProductType> >
331 template<
typename Dest>
332 static void scaleAndAddTo(Dest& dst,
const Lhs& lhs,
const RhsView& rhsView,
const typename Dest::Scalar& alpha)
334 typedef typename RhsView::_MatrixTypeNested Rhs;
337 LhsNested lhsNested(lhs);
338 RhsNested rhsNested(rhsView.matrix());
342 internal::sparse_selfadjoint_time_dense_product<RhsView::Mode>(rhsNested.transpose(), lhsNested.transpose(), dstT, alpha);
349 template<
typename LhsView,
typename Rhs,
int ProductTag>
351 :
public evaluator<typename Product<typename Rhs::PlainObject, Rhs, DefaultProduct>::PlainObject>
354 typedef typename XprType::PlainObject PlainObject;
358 : m_lhs(xpr.lhs()), m_result(xpr.rows(), xpr.cols())
360 ::new (static_cast<Base*>(
this)) Base(m_result);
365 typename Rhs::PlainObject m_lhs;
366 PlainObject m_result;
369 template<
typename Lhs,
typename RhsView,
int ProductTag>
371 :
public evaluator<typename Product<Lhs, typename Lhs::PlainObject, DefaultProduct>::PlainObject>
374 typedef typename XprType::PlainObject PlainObject;
378 : m_rhs(xpr.rhs()), m_result(xpr.rows(), xpr.cols())
380 ::new (static_cast<Base*>(
this)) Base(m_result);
385 typename Lhs::PlainObject m_rhs;
386 PlainObject m_result;
396 template<
int Mode,
typename MatrixType,
int DestOrder>
399 typedef typename MatrixType::StorageIndex StorageIndex;
400 typedef typename MatrixType::Scalar Scalar;
406 MatEval matEval(mat);
407 Dest& dest(_dest.derived());
409 StorageOrderMatch = int(Dest::IsRowMajor) == int(MatrixType::IsRowMajor)
412 Index size = mat.rows();
416 dest.resize(size,size);
417 for(
Index j = 0; j<size; ++j)
419 Index jp = perm ? perm[j] : j;
420 for(MatIterator it(matEval,j); it; ++it)
422 Index i = it.index();
425 Index ip = perm ? perm[i] : i;
427 count[StorageOrderMatch ? jp : ip]++;
430 else if(( Mode==
Lower && r>c) || ( Mode==
Upper && r<c))
437 Index nnz = count.sum();
440 dest.resizeNonZeros(nnz);
441 dest.outerIndexPtr()[0] = 0;
442 for(
Index j=0; j<size; ++j)
443 dest.outerIndexPtr()[j+1] = dest.outerIndexPtr()[j] + count[j];
444 for(
Index j=0; j<size; ++j)
445 count[j] = dest.outerIndexPtr()[j];
448 for(StorageIndex j = 0; j<size; ++j)
450 for(MatIterator it(matEval,j); it; ++it)
452 StorageIndex i = internal::convert_index<StorageIndex>(it.index());
456 StorageIndex jp = perm ? perm[j] : j;
457 StorageIndex ip = perm ? perm[i] : i;
461 Index k = count[StorageOrderMatch ? jp : ip]++;
462 dest.innerIndexPtr()[k] = StorageOrderMatch ? ip : jp;
463 dest.valuePtr()[k] = it.value();
467 Index k = count[ip]++;
468 dest.innerIndexPtr()[k] = ip;
469 dest.valuePtr()[k] = it.value();
471 else if(( (Mode&
Lower)==Lower && r>c) || ( (Mode&
Upper)==Upper && r<c))
473 if(!StorageOrderMatch)
475 Index k = count[jp]++;
476 dest.innerIndexPtr()[k] = ip;
477 dest.valuePtr()[k] = it.value();
479 dest.innerIndexPtr()[k] = jp;
480 dest.valuePtr()[k] = numext::conj(it.value());
486 template<
int _SrcMode,
int _DstMode,
typename MatrixType,
int DstOrder>
489 typedef typename MatrixType::StorageIndex StorageIndex;
490 typedef typename MatrixType::Scalar Scalar;
498 StorageOrderMatch = int(SrcOrder) == int(DstOrder),
503 MatEval matEval(mat);
505 Index size = mat.rows();
508 dest.resize(size,size);
509 for(StorageIndex j = 0; j<size; ++j)
511 StorageIndex jp = perm ? perm[j] : j;
512 for(MatIterator it(matEval,j); it; ++it)
514 StorageIndex i = it.index();
515 if((
int(SrcMode)==
int(
Lower) && i<j) || (
int(SrcMode)==
int(
Upper) && i>j))
518 StorageIndex ip = perm ? perm[i] : i;
519 count[int(DstMode)==int(
Lower) ? (std::min)(ip,jp) : (std::max)(ip,jp)]++;
522 dest.outerIndexPtr()[0] = 0;
523 for(
Index j=0; j<size; ++j)
524 dest.outerIndexPtr()[j+1] = dest.outerIndexPtr()[j] + count[j];
525 dest.resizeNonZeros(dest.outerIndexPtr()[size]);
526 for(
Index j=0; j<size; ++j)
527 count[j] = dest.outerIndexPtr()[j];
529 for(StorageIndex j = 0; j<size; ++j)
532 for(MatIterator it(matEval,j); it; ++it)
534 StorageIndex i = it.index();
535 if((
int(SrcMode)==
int(
Lower) && i<j) || (
int(SrcMode)==
int(
Upper) && i>j))
538 StorageIndex jp = perm ? perm[j] : j;
539 StorageIndex ip = perm? perm[i] : i;
541 Index k = count[int(DstMode)==int(
Lower) ? (std::min)(ip,jp) : (std::max)(ip,jp)]++;
542 dest.innerIndexPtr()[k] = int(DstMode)==int(
Lower) ? (std::max)(ip,jp) : (std::min)(ip,jp);
544 if(!StorageOrderMatch) std::swap(ip,jp);
545 if( ((
int(DstMode)==
int(
Lower) && ip<jp) || (
int(DstMode)==
int(
Upper) && ip>jp)))
546 dest.valuePtr()[k] = numext::conj(it.value());
548 dest.valuePtr()[k] = it.value();
559 template<
typename MatrixType,
int Mode>
565 template<
typename MatrixType,
int Mode>
567 :
public EigenBase<SparseSymmetricPermutationProduct<MatrixType,Mode> >
570 typedef typename MatrixType::Scalar Scalar;
571 typedef typename MatrixType::StorageIndex StorageIndex;
580 typedef typename MatrixType::Nested MatrixTypeNested;
584 : m_matrix(mat), m_perm(perm)
587 inline Index rows()
const {
return m_matrix.rows(); }
588 inline Index cols()
const {
return m_matrix.cols(); }
590 const NestedExpression& matrix()
const {
return m_matrix; }
591 const Perm& perm()
const {
return m_perm; }
594 MatrixTypeNested m_matrix;
601 template<
typename DstXprType,
typename MatrixType,
int Mode,
typename Scalar>
605 typedef typename DstXprType::StorageIndex DstIndex;
606 template<
int Options>
611 internal::permute_symm_to_fullsymm<Mode>(src.matrix(),tmp,src.perm().
indices().data());
615 template<
typename DestType,
unsigned int DestMode>
618 internal::permute_symm_to_symm<Mode,DestMode>(src.matrix(),dst.matrix(),src.perm().
indices().data());
626 #endif // EIGEN_SPARSE_SELFADJOINTVIEW_H Definition: NonLinearOptimization.cpp:108
Product< SparseSelfAdjointView, OtherDerived > operator*(const SparseMatrixBase< OtherDerived > &rhs) const
Definition: SparseSelfAdjointView.h:80
Storage order is column major (see TopicStorageOrders).
Definition: Constants.h:320
Definition: Constants.h:526
Expression of the product of two arbitrary matrices or vectors.
Definition: Product.h:71
Definition: ForwardDeclarations.h:162
A versatible sparse matrix representation.
Definition: SparseMatrix.h:92
Definition: SparseSelfAdjointView.h:566
Definition: SparseSelfAdjointView.h:216
Expression of the transpose of a matrix.
Definition: Transpose.h:52
friend Product< OtherDerived, SparseSelfAdjointView > operator*(const SparseMatrixBase< OtherDerived > &lhs, const SparseSelfAdjointView &rhs)
Definition: SparseSelfAdjointView.h:92
friend Product< OtherDerived, SparseSelfAdjointView > operator*(const MatrixBase< OtherDerived > &lhs, const SparseSelfAdjointView &rhs)
Efficient dense vector/matrix times sparse self-adjoint matrix product.
Definition: SparseSelfAdjointView.h:108
Definition: CoreEvaluators.h:90
Namespace containing all symbols from the Eigen library.
Definition: bench_norm.cpp:85
Pseudo expression to manipulate a triangular sparse matrix as a selfadjoint matrix.
Definition: SparseSelfAdjointView.h:43
Definition: Constants.h:521
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:37
const unsigned int RowMajorBit
for a matrix, this means that the storage order is row-major.
Definition: Constants.h:61
Definition: Constants.h:512
Definition: ProductEvaluators.h:343
Definition: AssignmentFunctors.h:21
Definition: AssignEvaluator.h:753
Common base class for all classes T such that MatrixBase has an operator=(T) and a constructor Matrix...
Definition: EigenBase.h:28
Base class of any sparse matrices or sparse expressions.
Definition: ForwardDeclarations.h:281
View matrix as a lower triangular matrix.
Definition: Constants.h:204
Definition: ProductEvaluators.h:86
const IndicesType & indices() const
const version of indices().
Definition: PermutationMatrix.h:388
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:33
A sparse matrix class designed for matrix assembly purpose.
Definition: SparseUtil.h:53
View matrix as an upper triangular matrix.
Definition: Constants.h:206
Product< SparseSelfAdjointView, OtherDerived > operator*(const MatrixBase< OtherDerived > &rhs) const
Efficient sparse self-adjoint matrix times dense vector/matrix product.
Definition: SparseSelfAdjointView.h:100
Definition: SparseUtil.h:138
Definition: BandTriangularSolver.h:13
SparseSelfAdjointView & rankUpdate(const SparseMatrixBase< DerivedU > &u, const Scalar &alpha=Scalar(1))
Perform a symmetric rank K update of the selfadjoint matrix *this: where u is a vector or matrix...
Definition: CoreEvaluators.h:79
SparseSymmetricPermutationProduct< _MatrixTypeNested, Mode > twistedBy(const PermutationMatrix< Dynamic, Dynamic, StorageIndex > &perm) const
Definition: SparseSelfAdjointView.h:126
Definition: AssignEvaluator.h:740
Storage order is row major (see TopicStorageOrders).
Definition: Constants.h:322
Determines whether the given binary operation of two numeric types is allowed and what the scalar ret...
Definition: XprHelper.h:757
Generic expression where a coefficient-wise unary operator is applied to an expression.
Definition: CwiseUnaryOp.h:55
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
Definition: ForwardDeclarations.h:17
Definition: SparseAssign.h:61