SparseAssign.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) 2008-2014 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_SPARSEASSIGN_H
11 #define EIGEN_SPARSEASSIGN_H
12 
13 #include "./InternalHeaderCheck.h"
14 
15 namespace Eigen {
16 
17 template<typename Derived>
18 template<typename OtherDerived>
20 {
21  internal::call_assignment_no_alias(derived(), other.derived());
22  return derived();
23 }
24 
25 template<typename Derived>
26 template<typename OtherDerived>
28 {
29  // TODO use the evaluator mechanism
30  other.evalTo(derived());
31  return derived();
32 }
33 
34 template<typename Derived>
35 template<typename OtherDerived>
37 {
38  // by default sparse evaluation do not alias, so we can safely bypass the generic call_assignment routine
39  internal::Assignment<Derived,OtherDerived,internal::assign_op<Scalar,typename OtherDerived::Scalar> >
40  ::run(derived(), other.derived(), internal::assign_op<Scalar,typename OtherDerived::Scalar>());
41  return derived();
42 }
43 
44 template<typename Derived>
45 inline Derived& SparseMatrixBase<Derived>::operator=(const Derived& other)
46 {
47  internal::call_assignment_no_alias(derived(), other.derived());
48  return derived();
49 }
50 
51 namespace internal {
52 
53 template<>
54 struct storage_kind_to_evaluator_kind<Sparse> {
55  typedef IteratorBased Kind;
56 };
57 
58 template<>
59 struct storage_kind_to_shape<Sparse> {
60  typedef SparseShape Shape;
61 };
62 
63 struct Sparse2Sparse {};
64 struct Sparse2Dense {};
65 
66 template<> struct AssignmentKind<SparseShape, SparseShape> { typedef Sparse2Sparse Kind; };
67 template<> struct AssignmentKind<SparseShape, SparseTriangularShape> { typedef Sparse2Sparse Kind; };
68 template<> struct AssignmentKind<DenseShape, SparseShape> { typedef Sparse2Dense Kind; };
69 template<> struct AssignmentKind<DenseShape, SparseTriangularShape> { typedef Sparse2Dense Kind; };
70 
71 
72 template<typename DstXprType, typename SrcXprType>
73 void assign_sparse_to_sparse(DstXprType &dst, const SrcXprType &src)
74 {
75  typedef typename DstXprType::Scalar Scalar;
76  typedef internal::evaluator<DstXprType> DstEvaluatorType;
77  typedef internal::evaluator<SrcXprType> SrcEvaluatorType;
78 
79  SrcEvaluatorType srcEvaluator(src);
80 
81  constexpr bool transpose = (DstEvaluatorType::Flags & RowMajorBit) != (SrcEvaluatorType::Flags & RowMajorBit);
82  const Index outerEvaluationSize = (SrcEvaluatorType::Flags&RowMajorBit) ? src.rows() : src.cols();
83 
84  Index reserveSize = 0;
85  for (Index j = 0; j < outerEvaluationSize; ++j)
86  for (typename SrcEvaluatorType::InnerIterator it(srcEvaluator, j); it; ++it)
87  reserveSize++;
88 
89  if ((!transpose) && src.isRValue())
90  {
91  // eval without temporary
92  dst.resize(src.rows(), src.cols());
93  dst.setZero();
94  dst.reserve(reserveSize);
95  for (Index j=0; j<outerEvaluationSize; ++j)
96  {
97  dst.startVec(j);
98  for (typename SrcEvaluatorType::InnerIterator it(srcEvaluator, j); it; ++it)
99  {
100  Scalar v = it.value();
101  dst.insertBackByOuterInner(j,it.index()) = v;
102  }
103  }
104  dst.finalize();
105  }
106  else
107  {
108  // eval through a temporary
109  eigen_assert(( ((internal::traits<DstXprType>::SupportedAccessPatterns & OuterRandomAccessPattern)==OuterRandomAccessPattern) ||
110  (!((DstEvaluatorType::Flags & RowMajorBit) != (SrcEvaluatorType::Flags & RowMajorBit)))) &&
111  "the transpose operation is supposed to be handled in SparseMatrix::operator=");
112 
113  enum { Flip = (DstEvaluatorType::Flags & RowMajorBit) != (SrcEvaluatorType::Flags & RowMajorBit) };
114 
115 
116  DstXprType temp(src.rows(), src.cols());
117 
118  temp.reserve(reserveSize);
119  for (Index j=0; j<outerEvaluationSize; ++j)
120  {
121  temp.startVec(j);
122  for (typename SrcEvaluatorType::InnerIterator it(srcEvaluator, j); it; ++it)
123  {
124  Scalar v = it.value();
125  temp.insertBackByOuterInner(Flip?it.index():j,Flip?j:it.index()) = v;
126  }
127  }
128  temp.finalize();
129 
130  dst = temp.markAsRValue();
131  }
132 }
133 
134 // Generic Sparse to Sparse assignment
135 template< typename DstXprType, typename SrcXprType, typename Functor>
136 struct Assignment<DstXprType, SrcXprType, Functor, Sparse2Sparse>
137 {
138  static void run(DstXprType &dst, const SrcXprType &src, const internal::assign_op<typename DstXprType::Scalar,typename SrcXprType::Scalar> &/*func*/)
139  {
140  assign_sparse_to_sparse(dst.derived(), src.derived());
141  }
142 };
143 
144 // Generic Sparse to Dense assignment
145 template< typename DstXprType, typename SrcXprType, typename Functor, typename Weak>
146 struct Assignment<DstXprType, SrcXprType, Functor, Sparse2Dense, Weak>
147 {
148  static void run(DstXprType &dst, const SrcXprType &src, const Functor &func)
149  {
150  if(internal::is_same<Functor,internal::assign_op<typename DstXprType::Scalar,typename SrcXprType::Scalar> >::value)
151  dst.setZero();
152 
153  internal::evaluator<SrcXprType> srcEval(src);
154  resize_if_allowed(dst, src, func);
155  internal::evaluator<DstXprType> dstEval(dst);
156 
157  const Index outerEvaluationSize = (internal::evaluator<SrcXprType>::Flags&RowMajorBit) ? src.rows() : src.cols();
158  for (Index j=0; j<outerEvaluationSize; ++j)
159  for (typename internal::evaluator<SrcXprType>::InnerIterator i(srcEval,j); i; ++i)
160  func.assignCoeff(dstEval.coeffRef(i.row(),i.col()), i.value());
161  }
162 };
163 
164 // Specialization for dense ?= dense +/- sparse and dense ?= sparse +/- dense
165 template<typename DstXprType, typename Func1, typename Func2>
166 struct assignment_from_dense_op_sparse
167 {
168  template<typename SrcXprType, typename InitialFunc>
169  static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE
170  void run(DstXprType &dst, const SrcXprType &src, const InitialFunc& /*func*/)
171  {
172  #ifdef EIGEN_SPARSE_ASSIGNMENT_FROM_DENSE_OP_SPARSE_PLUGIN
173  EIGEN_SPARSE_ASSIGNMENT_FROM_DENSE_OP_SPARSE_PLUGIN
174  #endif
175 
176  call_assignment_no_alias(dst, src.lhs(), Func1());
177  call_assignment_no_alias(dst, src.rhs(), Func2());
178  }
179 
180  // Specialization for dense1 = sparse + dense2; -> dense1 = dense2; dense1 += sparse;
181  template<typename Lhs, typename Rhs, typename Scalar>
182  static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE
183  std::enable_if_t<internal::is_same<typename internal::evaluator_traits<Rhs>::Shape,DenseShape>::value>
184  run(DstXprType &dst, const CwiseBinaryOp<internal::scalar_sum_op<Scalar,Scalar>, const Lhs, const Rhs> &src,
185  const internal::assign_op<typename DstXprType::Scalar,Scalar>& /*func*/)
186  {
187  #ifdef EIGEN_SPARSE_ASSIGNMENT_FROM_SPARSE_ADD_DENSE_PLUGIN
188  EIGEN_SPARSE_ASSIGNMENT_FROM_SPARSE_ADD_DENSE_PLUGIN
189  #endif
190 
191  // Apply the dense matrix first, then the sparse one.
192  call_assignment_no_alias(dst, src.rhs(), Func1());
193  call_assignment_no_alias(dst, src.lhs(), Func2());
194  }
195 
196  // Specialization for dense1 = sparse - dense2; -> dense1 = -dense2; dense1 += sparse;
197  template<typename Lhs, typename Rhs, typename Scalar>
198  static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE
199  std::enable_if_t<internal::is_same<typename internal::evaluator_traits<Rhs>::Shape,DenseShape>::value>
200  run(DstXprType &dst, const CwiseBinaryOp<internal::scalar_difference_op<Scalar,Scalar>, const Lhs, const Rhs> &src,
201  const internal::assign_op<typename DstXprType::Scalar,Scalar>& /*func*/)
202  {
203  #ifdef EIGEN_SPARSE_ASSIGNMENT_FROM_SPARSE_SUB_DENSE_PLUGIN
204  EIGEN_SPARSE_ASSIGNMENT_FROM_SPARSE_SUB_DENSE_PLUGIN
205  #endif
206 
207  // Apply the dense matrix first, then the sparse one.
208  call_assignment_no_alias(dst, -src.rhs(), Func1());
209  call_assignment_no_alias(dst, src.lhs(), add_assign_op<typename DstXprType::Scalar,typename Lhs::Scalar>());
210  }
211 };
212 
213 #define EIGEN_CATCH_ASSIGN_DENSE_OP_SPARSE(ASSIGN_OP,BINOP,ASSIGN_OP2) \
214  template< typename DstXprType, typename Lhs, typename Rhs, typename Scalar> \
215  struct Assignment<DstXprType, CwiseBinaryOp<internal::BINOP<Scalar,Scalar>, const Lhs, const Rhs>, internal::ASSIGN_OP<typename DstXprType::Scalar,Scalar>, \
216  Sparse2Dense, \
217  std::enable_if_t< internal::is_same<typename internal::evaluator_traits<Lhs>::Shape,DenseShape>::value \
218  || internal::is_same<typename internal::evaluator_traits<Rhs>::Shape,DenseShape>::value>> \
219  : assignment_from_dense_op_sparse<DstXprType, internal::ASSIGN_OP<typename DstXprType::Scalar,typename Lhs::Scalar>, internal::ASSIGN_OP2<typename DstXprType::Scalar,typename Rhs::Scalar> > \
220  {}
221 
222 EIGEN_CATCH_ASSIGN_DENSE_OP_SPARSE(assign_op, scalar_sum_op,add_assign_op);
223 EIGEN_CATCH_ASSIGN_DENSE_OP_SPARSE(add_assign_op,scalar_sum_op,add_assign_op);
224 EIGEN_CATCH_ASSIGN_DENSE_OP_SPARSE(sub_assign_op,scalar_sum_op,sub_assign_op);
225 
226 EIGEN_CATCH_ASSIGN_DENSE_OP_SPARSE(assign_op, scalar_difference_op,sub_assign_op);
227 EIGEN_CATCH_ASSIGN_DENSE_OP_SPARSE(add_assign_op,scalar_difference_op,sub_assign_op);
228 EIGEN_CATCH_ASSIGN_DENSE_OP_SPARSE(sub_assign_op,scalar_difference_op,add_assign_op);
229 
230 
231 // Specialization for "dst = dec.solve(rhs)"
232 // NOTE we need to specialize it for Sparse2Sparse to avoid ambiguous specialization error
233 template<typename DstXprType, typename DecType, typename RhsType, typename Scalar>
234 struct Assignment<DstXprType, Solve<DecType,RhsType>, internal::assign_op<Scalar,Scalar>, Sparse2Sparse>
235 {
236  typedef Solve<DecType,RhsType> SrcXprType;
237  static void run(DstXprType &dst, const SrcXprType &src, const internal::assign_op<Scalar,Scalar> &)
238  {
239  Index dstRows = src.rows();
240  Index dstCols = src.cols();
241  if((dst.rows()!=dstRows) || (dst.cols()!=dstCols))
242  dst.resize(dstRows, dstCols);
243 
244  src.dec()._solve_impl(src.rhs(), dst);
245  }
246 };
247 
248 struct Diagonal2Sparse {};
249 
250 template<> struct AssignmentKind<SparseShape,DiagonalShape> { typedef Diagonal2Sparse Kind; };
251 
252 template< typename DstXprType, typename SrcXprType, typename Functor>
253 struct Assignment<DstXprType, SrcXprType, Functor, Diagonal2Sparse>
254 {
255  typedef typename DstXprType::StorageIndex StorageIndex;
256  typedef typename DstXprType::Scalar Scalar;
257 
258  template<int Options, typename AssignFunc>
259  static void run(SparseMatrix<Scalar,Options,StorageIndex> &dst, const SrcXprType &src, const AssignFunc &func)
260  { dst.assignDiagonal(src.diagonal(), func); }
261 
262  template<typename DstDerived>
263  static void run(SparseMatrixBase<DstDerived> &dst, const SrcXprType &src, const internal::assign_op<typename DstXprType::Scalar,typename SrcXprType::Scalar> &/*func*/)
264  { dst.derived().diagonal() = src.diagonal(); }
265 
266  template<typename DstDerived>
267  static void run(SparseMatrixBase<DstDerived> &dst, const SrcXprType &src, const internal::add_assign_op<typename DstXprType::Scalar,typename SrcXprType::Scalar> &/*func*/)
268  { dst.derived().diagonal() += src.diagonal(); }
269 
270  template<typename DstDerived>
271  static void run(SparseMatrixBase<DstDerived> &dst, const SrcXprType &src, const internal::sub_assign_op<typename DstXprType::Scalar,typename SrcXprType::Scalar> &/*func*/)
272  { dst.derived().diagonal() -= src.diagonal(); }
273 };
274 } // end namespace internal
275 
276 } // end namespace Eigen
277 
278 #endif // EIGEN_SPARSEASSIGN_H
Array< int, Dynamic, 1 > v
#define EIGEN_DEVICE_FUNC
Definition: Macros.h:883
#define eigen_assert(x)
Definition: Macros.h:902
void evalTo(Dest &dst) const
Definition: ReturnByValue.h:63
Pseudo expression representing a solving operation.
Definition: Solve.h:65
EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition: Solve.h:74
Base class of any sparse matrices or sparse expressions.
Derived & operator=(const EigenBase< OtherDerived > &other)
Definition: SparseAssign.h:19
const unsigned int RowMajorBit
Definition: Constants.h:68
void resize_if_allowed(DstXprType &dst, const SrcXprType &src, const Functor &)
void assign_sparse_to_sparse(DstXprType &dst, const SrcXprType &src)
Definition: SparseAssign.h:73
EIGEN_CONSTEXPR void call_assignment_no_alias(Dst &dst, const Src &src, const Func &func)
EIGEN_CATCH_ASSIGN_DENSE_OP_SPARSE(assign_op, scalar_sum_op, add_assign_op)
: InteropHeaders
Definition: Core:139
const int OuterRandomAccessPattern
Definition: SparseUtil.h:51
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
Derived & derived()
Definition: EigenBase.h:48
std::ptrdiff_t j