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 
66 template <typename MatrixType, int UpLo = Lower>
68 
80 template <typename MatrixType, int UpLo = Lower>
82 
93 template <typename MatrixType>
95 
106 template <typename MatrixType>
108 
109 namespace internal {
110 template <typename T>
111 struct AccelFactorizationDeleter {
112  void operator()(T* sym) {
113  if (sym) {
114  SparseCleanup(*sym);
115  delete sym;
116  sym = nullptr;
117  }
118  }
119 };
120 
121 template <typename DenseVecT, typename DenseMatT, typename SparseMatT, typename NumFactT>
122 struct SparseTypesTraitBase {
123  typedef DenseVecT AccelDenseVector;
124  typedef DenseMatT AccelDenseMatrix;
125  typedef SparseMatT AccelSparseMatrix;
126 
127  typedef SparseOpaqueSymbolicFactorization SymbolicFactorization;
128  typedef NumFactT NumericFactorization;
129 
130  typedef AccelFactorizationDeleter<SymbolicFactorization> SymbolicFactorizationDeleter;
131  typedef AccelFactorizationDeleter<NumericFactorization> NumericFactorizationDeleter;
132 };
133 
134 template <typename Scalar>
135 struct SparseTypesTrait {};
136 
137 template <>
138 struct SparseTypesTrait<double> : SparseTypesTraitBase<DenseVector_Double, DenseMatrix_Double, SparseMatrix_Double,
139  SparseOpaqueFactorization_Double> {};
140 
141 template <>
142 struct SparseTypesTrait<float>
143  : SparseTypesTraitBase<DenseVector_Float, DenseMatrix_Float, SparseMatrix_Float, SparseOpaqueFactorization_Float> {
144 };
145 
146 } // end namespace internal
147 
148 template <typename MatrixType_, int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
149 class AccelerateImpl : public SparseSolverBase<AccelerateImpl<MatrixType_, UpLo_, Solver_, EnforceSquare_> > {
150  protected:
152  using Base::derived;
153  using Base::m_isInitialized;
154 
155  public:
156  using Base::_solve_impl;
157 
158  typedef MatrixType_ MatrixType;
159  typedef typename MatrixType::Scalar Scalar;
160  typedef typename MatrixType::StorageIndex StorageIndex;
162  enum { UpLo = UpLo_ };
163 
164  using AccelDenseVector = typename internal::SparseTypesTrait<Scalar>::AccelDenseVector;
165  using AccelDenseMatrix = typename internal::SparseTypesTrait<Scalar>::AccelDenseMatrix;
166  using AccelSparseMatrix = typename internal::SparseTypesTrait<Scalar>::AccelSparseMatrix;
167  using SymbolicFactorization = typename internal::SparseTypesTrait<Scalar>::SymbolicFactorization;
168  using NumericFactorization = typename internal::SparseTypesTrait<Scalar>::NumericFactorization;
169  using SymbolicFactorizationDeleter = typename internal::SparseTypesTrait<Scalar>::SymbolicFactorizationDeleter;
170  using NumericFactorizationDeleter = typename internal::SparseTypesTrait<Scalar>::NumericFactorizationDeleter;
171 
173  m_isInitialized = false;
174 
175  auto check_flag_set = [](int value, int flag) { return ((value & flag) == flag); };
176 
177  if (check_flag_set(UpLo_, Symmetric)) {
178  m_sparseKind = SparseSymmetric;
179  m_triType = (UpLo_ & Lower) ? SparseLowerTriangle : SparseUpperTriangle;
180  } else if (check_flag_set(UpLo_, UnitLower)) {
181  m_sparseKind = SparseUnitTriangular;
182  m_triType = SparseLowerTriangle;
183  } else if (check_flag_set(UpLo_, UnitUpper)) {
184  m_sparseKind = SparseUnitTriangular;
185  m_triType = SparseUpperTriangle;
186  } else if (check_flag_set(UpLo_, StrictlyLower)) {
187  m_sparseKind = SparseTriangular;
188  m_triType = SparseLowerTriangle;
189  } else if (check_flag_set(UpLo_, StrictlyUpper)) {
190  m_sparseKind = SparseTriangular;
191  m_triType = SparseUpperTriangle;
192  } else if (check_flag_set(UpLo_, Lower)) {
193  m_sparseKind = SparseTriangular;
194  m_triType = SparseLowerTriangle;
195  } else if (check_flag_set(UpLo_, Upper)) {
196  m_sparseKind = SparseTriangular;
197  m_triType = SparseUpperTriangle;
198  } else {
199  m_sparseKind = SparseOrdinary;
200  m_triType = (UpLo_ & Lower) ? SparseLowerTriangle : SparseUpperTriangle;
201  }
202 
203  m_order = SparseOrderDefault;
204  }
205 
206  explicit AccelerateImpl(const MatrixType& matrix) : AccelerateImpl() { compute(matrix); }
207 
209 
210  inline Index cols() const { return m_nCols; }
211  inline Index rows() const { return m_nRows; }
212 
214  eigen_assert(m_isInitialized && "Decomposition is not initialized.");
215  return m_info;
216  }
217 
218  void analyzePattern(const MatrixType& matrix);
219 
220  void factorize(const MatrixType& matrix);
221 
222  void compute(const MatrixType& matrix);
223 
224  template <typename Rhs, typename Dest>
225  void _solve_impl(const MatrixBase<Rhs>& b, MatrixBase<Dest>& dest) const;
226 
228  void setOrder(SparseOrder_t order) { m_order = order; }
229 
230  private:
231  template <typename T>
232  void buildAccelSparseMatrix(const SparseMatrix<T>& a, AccelSparseMatrix& A, std::vector<long>& columnStarts) {
233  const Index nColumnsStarts = a.cols() + 1;
234 
235  columnStarts.resize(nColumnsStarts);
236 
237  for (Index i = 0; i < nColumnsStarts; i++) columnStarts[i] = a.outerIndexPtr()[i];
238 
239  SparseAttributes_t attributes{};
240  attributes.transpose = false;
241  attributes.triangle = m_triType;
242  attributes.kind = m_sparseKind;
243 
244  SparseMatrixStructure structure{};
245  structure.attributes = attributes;
246  structure.rowCount = static_cast<int>(a.rows());
247  structure.columnCount = static_cast<int>(a.cols());
248  structure.blockSize = 1;
249  structure.columnStarts = columnStarts.data();
250  structure.rowIndices = const_cast<int*>(a.innerIndexPtr());
251 
252  A.structure = structure;
253  A.data = const_cast<T*>(a.valuePtr());
254  }
255 
257  m_numericFactorization.reset(nullptr);
258 
259  SparseSymbolicFactorOptions opts{};
260  opts.control = SparseDefaultControl;
261  opts.orderMethod = m_order;
262  opts.order = nullptr;
263  opts.ignoreRowsAndColumns = nullptr;
264  opts.malloc = malloc;
265  opts.free = free;
266  opts.reportError = nullptr;
267 
268  m_symbolicFactorization.reset(new SymbolicFactorization(SparseFactor(Solver_, A.structure, opts)));
269 
270  SparseStatus_t status = m_symbolicFactorization->status;
271 
272  updateInfoStatus(status);
273 
274  if (status != SparseStatusOK) m_symbolicFactorization.reset(nullptr);
275  }
276 
278  SparseStatus_t status = SparseStatusReleased;
279 
282 
283  status = m_numericFactorization->status;
284 
285  if (status != SparseStatusOK) m_numericFactorization.reset(nullptr);
286  }
287 
288  updateInfoStatus(status);
289  }
290 
291  protected:
292  void updateInfoStatus(SparseStatus_t status) const {
293  switch (status) {
294  case SparseStatusOK:
295  m_info = Success;
296  break;
297  case SparseFactorizationFailed:
298  case SparseMatrixIsSingular:
300  break;
301  case SparseInternalError:
302  case SparseParameterError:
303  case SparseStatusReleased:
304  default:
306  break;
307  }
308  }
309 
312  std::unique_ptr<SymbolicFactorization, SymbolicFactorizationDeleter> m_symbolicFactorization;
313  std::unique_ptr<NumericFactorization, NumericFactorizationDeleter> m_numericFactorization;
314  SparseKind_t m_sparseKind;
315  SparseTriangle_t m_triType;
316  SparseOrder_t m_order;
317 };
318 
320 template <typename MatrixType_, int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
322  if (EnforceSquare_) eigen_assert(a.rows() == a.cols());
323 
324  m_nRows = a.rows();
325  m_nCols = a.cols();
326 
328  std::vector<long> columnStarts;
329 
330  buildAccelSparseMatrix(a, A, columnStarts);
331 
332  doAnalysis(A);
333 
334  if (m_symbolicFactorization) doFactorization(A);
335 
336  m_isInitialized = true;
337 }
338 
345 template <typename MatrixType_, int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
347  if (EnforceSquare_) eigen_assert(a.rows() == a.cols());
348 
349  m_nRows = a.rows();
350  m_nCols = a.cols();
351 
353  std::vector<long> columnStarts;
354 
355  buildAccelSparseMatrix(a, A, columnStarts);
356 
357  doAnalysis(A);
358 
359  m_isInitialized = true;
360 }
361 
368 template <typename MatrixType_, int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
370  eigen_assert(m_symbolicFactorization && "You must first call analyzePattern()");
371  eigen_assert(m_nRows == a.rows() && m_nCols == a.cols());
372 
373  if (EnforceSquare_) eigen_assert(a.rows() == a.cols());
374 
376  std::vector<long> columnStarts;
377 
378  buildAccelSparseMatrix(a, A, columnStarts);
379 
380  doFactorization(A);
381 }
382 
383 template <typename MatrixType_, int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
384 template <typename Rhs, typename Dest>
386  MatrixBase<Dest>& x) const {
387  if (!m_numericFactorization) {
388  m_info = InvalidInput;
389  return;
390  }
391 
392  eigen_assert(m_nRows == b.rows());
393  eigen_assert(((b.cols() == 1) || b.outerStride() == b.rows()));
394 
395  SparseStatus_t status = SparseStatusOK;
396 
397  Scalar* b_ptr = const_cast<Scalar*>(b.derived().data());
398  Scalar* x_ptr = const_cast<Scalar*>(x.derived().data());
399 
400  AccelDenseMatrix xmat{};
401  xmat.attributes = SparseAttributes_t();
402  xmat.columnCount = static_cast<int>(x.cols());
403  xmat.rowCount = static_cast<int>(x.rows());
404  xmat.columnStride = xmat.rowCount;
405  xmat.data = x_ptr;
406 
407  AccelDenseMatrix bmat{};
408  bmat.attributes = SparseAttributes_t();
409  bmat.columnCount = static_cast<int>(b.cols());
410  bmat.rowCount = static_cast<int>(b.rows());
411  bmat.columnStride = bmat.rowCount;
412  bmat.data = b_ptr;
413 
414  SparseSolve(*m_numericFactorization, bmat, xmat);
415 
416  updateInfoStatus(status);
417 }
418 
419 } // end namespace Eigen
420 
421 #endif // EIGEN_ACCELERATESUPPORT_H
Array< int, 3, 1 > b
MatrixXcf A
IndexedView_or_Block operator()(const RowIndices &rowIndices, const ColIndices &colIndices)
#define eigen_assert(x)
Definition: Macros.h:902
Matrix< float, 1, Dynamic > MatrixType
typename internal::SparseTypesTrait< Scalar >::SymbolicFactorizationDeleter SymbolicFactorizationDeleter
MatrixType::StorageIndex StorageIndex
typename internal::SparseTypesTrait< Scalar >::SymbolicFactorization SymbolicFactorization
std::unique_ptr< SymbolicFactorization, SymbolicFactorizationDeleter > m_symbolicFactorization
void analyzePattern(const MatrixType &matrix)
ComputationInfo info() const
std::unique_ptr< NumericFactorization, NumericFactorizationDeleter > m_numericFactorization
void doFactorization(AccelSparseMatrix &A)
typename internal::SparseTypesTrait< Scalar >::AccelDenseVector AccelDenseVector
void factorize(const MatrixType &matrix)
typename internal::SparseTypesTrait< Scalar >::AccelSparseMatrix AccelSparseMatrix
typename internal::SparseTypesTrait< Scalar >::NumericFactorizationDeleter NumericFactorizationDeleter
void doAnalysis(AccelSparseMatrix &A)
typename internal::SparseTypesTrait< Scalar >::AccelDenseMatrix AccelDenseMatrix
MatrixType::Scalar Scalar
SparseTriangle_t m_triType
AccelerateImpl(const MatrixType &matrix)
void setOrder(SparseOrder_t order)
typename internal::SparseTypesTrait< Scalar >::NumericFactorization NumericFactorization
void buildAccelSparseMatrix(const SparseMatrix< T > &a, AccelSparseMatrix &A, std::vector< long > &columnStarts)
void compute(const MatrixType &matrix)
void _solve_impl(const MatrixBase< Rhs > &b, MatrixBase< Dest > &dest) const
void updateInfoStatus(SparseStatus_t status) const
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
A versatible sparse matrix representation.
Definition: SparseMatrix.h:125
A base class for sparse solvers.
ComputationInfo
Definition: Constants.h:444
@ 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:448
@ InvalidInput
Definition: Constants.h:453
@ Success
Definition: Constants.h:446
: InteropHeaders
Definition: Core:139
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
const int Dynamic
Definition: Constants.h:24