AccelerateSupport.h
Go to the documentation of this file.
1 #ifndef EIGEN_ACCELERATESUPPORT_H
2 #define EIGEN_ACCELERATESUPPORT_H
3 
4 #include <Accelerate/Accelerate.h>
5 
6 #include <Eigen/Sparse>
7 
8 namespace Eigen {
9 
10 template <typename MatrixType_, int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
11 class AccelerateImpl;
12 
24 template <typename MatrixType, int UpLo = Lower>
26 
38 template <typename MatrixType, int UpLo = Lower>
40 
52 template <typename MatrixType, int UpLo = Lower>
54 
67 template <typename MatrixType, int UpLo = Lower>
69 
81 template <typename MatrixType, int UpLo = Lower>
83 
94 template <typename MatrixType>
96 
107 template <typename MatrixType>
109 
110 namespace internal {
111 template <typename T>
113  void operator()(T* sym) {
114  if (sym) {
115  SparseCleanup(*sym);
116  delete sym;
117  sym = nullptr;
118  }
119  }
120 };
121 
122 template <typename DenseVecT, typename DenseMatT, typename SparseMatT, typename NumFactT>
124  typedef DenseVecT AccelDenseVector;
125  typedef DenseMatT AccelDenseMatrix;
126  typedef SparseMatT AccelSparseMatrix;
127 
128  typedef SparseOpaqueSymbolicFactorization SymbolicFactorization;
129  typedef NumFactT NumericFactorization;
130 
133 };
134 
135 template <typename Scalar>
137 
138 template <>
139 struct SparseTypesTrait<double> : SparseTypesTraitBase<DenseVector_Double, DenseMatrix_Double, SparseMatrix_Double,
140  SparseOpaqueFactorization_Double> {};
141 
142 template <>
143 struct SparseTypesTrait<float>
144  : SparseTypesTraitBase<DenseVector_Float, DenseMatrix_Float, SparseMatrix_Float, SparseOpaqueFactorization_Float> {
145 };
146 
147 } // end namespace internal
148 
149 template <typename MatrixType_, int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
150 class AccelerateImpl : public SparseSolverBase<AccelerateImpl<MatrixType_, UpLo_, Solver_, EnforceSquare_> > {
151  protected:
153  using Base::derived;
154  using Base::m_isInitialized;
155 
156  public:
157  using Base::_solve_impl;
158 
159  typedef MatrixType_ MatrixType;
160  typedef typename MatrixType::Scalar Scalar;
161  typedef typename MatrixType::StorageIndex StorageIndex;
163  enum { UpLo = UpLo_ };
164 
172 
174  m_isInitialized = false;
175 
176  auto check_flag_set = [](int value, int flag) { return ((value & flag) == flag); };
177 
178  if (check_flag_set(UpLo_, Symmetric)) {
179  m_sparseKind = SparseSymmetric;
180  m_triType = (UpLo_ & Lower) ? SparseLowerTriangle : SparseUpperTriangle;
181  } else if (check_flag_set(UpLo_, UnitLower)) {
182  m_sparseKind = SparseUnitTriangular;
183  m_triType = SparseLowerTriangle;
184  } else if (check_flag_set(UpLo_, UnitUpper)) {
185  m_sparseKind = SparseUnitTriangular;
186  m_triType = SparseUpperTriangle;
187  } else if (check_flag_set(UpLo_, StrictlyLower)) {
188  m_sparseKind = SparseTriangular;
189  m_triType = SparseLowerTriangle;
190  } else if (check_flag_set(UpLo_, StrictlyUpper)) {
191  m_sparseKind = SparseTriangular;
192  m_triType = SparseUpperTriangle;
193  } else if (check_flag_set(UpLo_, Lower)) {
194  m_sparseKind = SparseTriangular;
195  m_triType = SparseLowerTriangle;
196  } else if (check_flag_set(UpLo_, Upper)) {
197  m_sparseKind = SparseTriangular;
198  m_triType = SparseUpperTriangle;
199  } else {
200  m_sparseKind = SparseOrdinary;
201  m_triType = (UpLo_ & Lower) ? SparseLowerTriangle : SparseUpperTriangle;
202  }
203 
204  m_order = SparseOrderDefault;
205  }
206 
208 
210 
211  inline Index cols() const { return m_nCols; }
212  inline Index rows() const { return m_nRows; }
213 
215  eigen_assert(m_isInitialized && "Decomposition is not initialized.");
216  return m_info;
217  }
218 
219  void analyzePattern(const MatrixType& matrix);
220 
221  void factorize(const MatrixType& matrix);
222 
223  void compute(const MatrixType& matrix);
224 
225  template <typename Rhs, typename Dest>
226  void _solve_impl(const MatrixBase<Rhs>& b, MatrixBase<Dest>& dest) const;
227 
229  void setOrder(SparseOrder_t order) { m_order = order; }
230 
231  private:
232  template <typename T>
233  void buildAccelSparseMatrix(const SparseMatrix<T>& a, AccelSparseMatrix& A, std::vector<long>& columnStarts) {
234  const Index nColumnsStarts = a.cols() + 1;
235 
236  columnStarts.resize(nColumnsStarts);
237 
238  for (Index i = 0; i < nColumnsStarts; i++) columnStarts[i] = a.outerIndexPtr()[i];
239 
240  SparseAttributes_t attributes{};
241  attributes.transpose = false;
242  attributes.triangle = m_triType;
243  attributes.kind = m_sparseKind;
244 
245  SparseMatrixStructure structure{};
246  structure.attributes = attributes;
247  structure.rowCount = static_cast<int>(a.rows());
248  structure.columnCount = static_cast<int>(a.cols());
249  structure.blockSize = 1;
250  structure.columnStarts = columnStarts.data();
251  structure.rowIndices = const_cast<int*>(a.innerIndexPtr());
252 
253  A.structure = structure;
254  A.data = const_cast<T*>(a.valuePtr());
255  }
256 
258  m_numericFactorization.reset(nullptr);
259 
260  SparseSymbolicFactorOptions opts{};
261  opts.control = SparseDefaultControl;
262  opts.orderMethod = m_order;
263  opts.order = nullptr;
264  opts.ignoreRowsAndColumns = nullptr;
265  opts.malloc = malloc;
266  opts.free = free;
267  opts.reportError = nullptr;
268 
269  m_symbolicFactorization.reset(new SymbolicFactorization(SparseFactor(Solver_, A.structure, opts)));
270 
271  SparseStatus_t status = m_symbolicFactorization->status;
272 
273  updateInfoStatus(status);
274 
275  if (status != SparseStatusOK) m_symbolicFactorization.reset(nullptr);
276  }
277 
279  SparseStatus_t status = SparseStatusReleased;
280 
283 
284  status = m_numericFactorization->status;
285 
286  if (status != SparseStatusOK) m_numericFactorization.reset(nullptr);
287  }
288 
289  updateInfoStatus(status);
290  }
291 
292  protected:
293  void updateInfoStatus(SparseStatus_t status) const {
294  switch (status) {
295  case SparseStatusOK:
296  m_info = Success;
297  break;
298  case SparseFactorizationFailed:
299  case SparseMatrixIsSingular:
301  break;
302  case SparseInternalError:
303  case SparseParameterError:
304  case SparseStatusReleased:
305  default:
307  break;
308  }
309  }
310 
313  std::unique_ptr<SymbolicFactorization, SymbolicFactorizationDeleter> m_symbolicFactorization;
314  std::unique_ptr<NumericFactorization, NumericFactorizationDeleter> m_numericFactorization;
315  SparseKind_t m_sparseKind;
316  SparseTriangle_t m_triType;
317  SparseOrder_t m_order;
318 };
319 
321 template <typename MatrixType_, int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
323  if (EnforceSquare_) eigen_assert(a.rows() == a.cols());
324 
325  m_nRows = a.rows();
326  m_nCols = a.cols();
327 
329  std::vector<long> columnStarts;
330 
331  buildAccelSparseMatrix(a, A, columnStarts);
332 
333  doAnalysis(A);
334 
335  if (m_symbolicFactorization) doFactorization(A);
336 
337  m_isInitialized = true;
338 }
339 
346 template <typename MatrixType_, int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
348  if (EnforceSquare_) eigen_assert(a.rows() == a.cols());
349 
350  m_nRows = a.rows();
351  m_nCols = a.cols();
352 
354  std::vector<long> columnStarts;
355 
356  buildAccelSparseMatrix(a, A, columnStarts);
357 
358  doAnalysis(A);
359 
360  m_isInitialized = true;
361 }
362 
370 template <typename MatrixType_, int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
372  eigen_assert(m_symbolicFactorization && "You must first call analyzePattern()");
373  eigen_assert(m_nRows == a.rows() && m_nCols == a.cols());
374 
375  if (EnforceSquare_) eigen_assert(a.rows() == a.cols());
376 
378  std::vector<long> columnStarts;
379 
380  buildAccelSparseMatrix(a, A, columnStarts);
381 
382  doFactorization(A);
383 }
384 
385 template <typename MatrixType_, int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
386 template <typename Rhs, typename Dest>
388  MatrixBase<Dest>& x) const {
389  if (!m_numericFactorization) {
390  m_info = InvalidInput;
391  return;
392  }
393 
394  eigen_assert(m_nRows == b.rows());
395  eigen_assert(((b.cols() == 1) || b.outerStride() == b.rows()));
396 
397  SparseStatus_t status = SparseStatusOK;
398 
399  Scalar* b_ptr = const_cast<Scalar*>(b.derived().data());
400  Scalar* x_ptr = const_cast<Scalar*>(x.derived().data());
401 
402  AccelDenseMatrix xmat{};
403  xmat.attributes = SparseAttributes_t();
404  xmat.columnCount = static_cast<int>(x.cols());
405  xmat.rowCount = static_cast<int>(x.rows());
406  xmat.columnStride = xmat.rowCount;
407  xmat.data = x_ptr;
408 
409  AccelDenseMatrix bmat{};
410  bmat.attributes = SparseAttributes_t();
411  bmat.columnCount = static_cast<int>(b.cols());
412  bmat.rowCount = static_cast<int>(b.rows());
413  bmat.columnStride = bmat.rowCount;
414  bmat.data = b_ptr;
415 
416  SparseSolve(*m_numericFactorization, bmat, xmat);
417 
418  updateInfoStatus(status);
419 }
420 
421 } // end namespace Eigen
422 
423 #endif // EIGEN_ACCELERATESUPPORT_H
int i
Definition: BiCGSTAB_step_by_step.cpp:9
#define eigen_assert(x)
Definition: Macros.h:910
Scalar * b
Definition: benchVecAdd.cpp:17
SCALAR Scalar
Definition: bench_gemm.cpp:45
MatrixXf MatrixType
Definition: benchmark-blocking-sizes.cpp:52
Definition: AccelerateSupport.h:150
typename internal::SparseTypesTrait< Scalar >::SymbolicFactorizationDeleter SymbolicFactorizationDeleter
Definition: AccelerateSupport.h:170
MatrixType::StorageIndex StorageIndex
Definition: AccelerateSupport.h:161
@ UpLo
Definition: AccelerateSupport.h:163
typename internal::SparseTypesTrait< Scalar >::SymbolicFactorization SymbolicFactorization
Definition: AccelerateSupport.h:168
Index cols() const
Definition: AccelerateSupport.h:211
Index m_nRows
Definition: AccelerateSupport.h:312
std::unique_ptr< SymbolicFactorization, SymbolicFactorizationDeleter > m_symbolicFactorization
Definition: AccelerateSupport.h:313
SparseOrder_t m_order
Definition: AccelerateSupport.h:317
@ ColsAtCompileTime
Definition: AccelerateSupport.h:162
@ MaxColsAtCompileTime
Definition: AccelerateSupport.h:162
Index rows() const
Definition: AccelerateSupport.h:212
SparseKind_t m_sparseKind
Definition: AccelerateSupport.h:315
void analyzePattern(const MatrixType &matrix)
Definition: AccelerateSupport.h:347
ComputationInfo info() const
Definition: AccelerateSupport.h:214
std::unique_ptr< NumericFactorization, NumericFactorizationDeleter > m_numericFactorization
Definition: AccelerateSupport.h:314
void doFactorization(AccelSparseMatrix &A)
Definition: AccelerateSupport.h:278
~AccelerateImpl()
Definition: AccelerateSupport.h:209
MatrixType_ MatrixType
Definition: AccelerateSupport.h:159
typename internal::SparseTypesTrait< Scalar >::AccelDenseVector AccelDenseVector
Definition: AccelerateSupport.h:165
void factorize(const MatrixType &matrix)
Definition: AccelerateSupport.h:371
typename internal::SparseTypesTrait< Scalar >::AccelSparseMatrix AccelSparseMatrix
Definition: AccelerateSupport.h:167
typename internal::SparseTypesTrait< Scalar >::NumericFactorizationDeleter NumericFactorizationDeleter
Definition: AccelerateSupport.h:171
void doAnalysis(AccelSparseMatrix &A)
Definition: AccelerateSupport.h:257
AccelerateImpl()
Definition: AccelerateSupport.h:173
typename internal::SparseTypesTrait< Scalar >::AccelDenseMatrix AccelDenseMatrix
Definition: AccelerateSupport.h:166
ComputationInfo m_info
Definition: AccelerateSupport.h:311
MatrixType::Scalar Scalar
Definition: AccelerateSupport.h:160
SparseTriangle_t m_triType
Definition: AccelerateSupport.h:316
AccelerateImpl(const MatrixType &matrix)
Definition: AccelerateSupport.h:207
void setOrder(SparseOrder_t order)
Definition: AccelerateSupport.h:229
bool m_isInitialized
Definition: SparseSolverBase.h:110
typename internal::SparseTypesTrait< Scalar >::NumericFactorization NumericFactorization
Definition: AccelerateSupport.h:169
void buildAccelSparseMatrix(const SparseMatrix< T > &a, AccelSparseMatrix &A, std::vector< long > &columnStarts)
Definition: AccelerateSupport.h:233
void compute(const MatrixType &matrix)
Definition: AccelerateSupport.h:322
Index m_nCols
Definition: AccelerateSupport.h:312
void _solve_impl(const MatrixBase< Rhs > &b, MatrixBase< Dest > &dest) const
Definition: AccelerateSupport.h:387
void updateInfoStatus(SparseStatus_t status) const
Definition: AccelerateSupport.h:293
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
The matrix class, also used for vectors and row-vectors.
Definition: Eigen/Eigen/src/Core/Matrix.h:186
constexpr EIGEN_DEVICE_FUNC const Scalar * data() const
Definition: PlainObjectBase.h:273
A versatible sparse matrix representation.
Definition: SparseMatrix.h:121
A base class for sparse solvers.
Definition: SparseSolverBase.h:67
void _solve_impl(const SparseMatrixBase< Rhs > &b, SparseMatrixBase< Dest > &dest) const
Definition: SparseSolverBase.h:104
bool m_isInitialized
Definition: SparseSolverBase.h:110
Derived & derived()
Definition: SparseSolverBase.h:76
Eigen::Map< Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic, Eigen::ColMajor >, 0, Eigen::OuterStride<> > matrix(T *data, int rows, int cols, int stride)
Definition: common.h:85
ComputationInfo
Definition: Constants.h:438
@ StrictlyLower
Definition: Constants.h:223
@ StrictlyUpper
Definition: Constants.h:225
@ UnitLower
Definition: Constants.h:219
@ Symmetric
Definition: Constants.h:229
@ UnitUpper
Definition: Constants.h:221
@ Lower
Definition: Constants.h:211
@ Upper
Definition: Constants.h:213
@ NumericalIssue
Definition: Constants.h:442
@ InvalidInput
Definition: Constants.h:447
@ Success
Definition: Constants.h:440
const Scalar * a
Definition: level2_cplx_impl.h:32
Namespace containing all symbols from the Eigen library.
Definition: bench_norm.cpp:70
squared absolute value
Definition: GlobalFunctions.h:87
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:83
const int Dynamic
Definition: Constants.h:25
Definition: Eigen_Colamd.h:49
list x
Definition: plotDoE.py:28
Definition: AccelerateSupport.h:112
void operator()(T *sym)
Definition: AccelerateSupport.h:113
Definition: AccelerateSupport.h:123
DenseVecT AccelDenseVector
Definition: AccelerateSupport.h:124
SparseOpaqueSymbolicFactorization SymbolicFactorization
Definition: AccelerateSupport.h:128
DenseMatT AccelDenseMatrix
Definition: AccelerateSupport.h:125
AccelFactorizationDeleter< NumericFactorization > NumericFactorizationDeleter
Definition: AccelerateSupport.h:132
NumFactT NumericFactorization
Definition: AccelerateSupport.h:129
AccelFactorizationDeleter< SymbolicFactorization > SymbolicFactorizationDeleter
Definition: AccelerateSupport.h:131
SparseMatT AccelSparseMatrix
Definition: AccelerateSupport.h:126
Definition: AccelerateSupport.h:136