1#ifndef EIGEN_ACCELERATESUPPORT_H
2#define EIGEN_ACCELERATESUPPORT_H
4#include <Accelerate/Accelerate.h>
10template <
typename MatrixType_,
int UpLo_, SparseFactorization_t Solver_,
bool EnforceSquare_>
24template <
typename MatrixType,
int UpLo = Lower>
25using AccelerateLLT = AccelerateImpl<MatrixType, UpLo | Symmetric, SparseFactorizationCholesky, true>;
38template <
typename MatrixType,
int UpLo = Lower>
39using AccelerateLDLT = AccelerateImpl<MatrixType, UpLo | Symmetric, SparseFactorizationLDLT, true>;
52template <
typename MatrixType,
int UpLo = Lower>
67template <
typename MatrixType,
int UpLo = Lower>
68using AccelerateLDLTSBK = AccelerateImpl<MatrixType, UpLo | Symmetric, SparseFactorizationLDLTSBK, true>;
81template <
typename MatrixType,
int UpLo = Lower>
82using AccelerateLDLTTPP = AccelerateImpl<MatrixType, UpLo | Symmetric, SparseFactorizationLDLTTPP, true>;
94template <
typename MatrixType>
95using AccelerateQR = AccelerateImpl<MatrixType, 0, SparseFactorizationQR, false>;
107template <
typename MatrixType>
112struct AccelFactorizationDeleter {
113 void operator()(T* sym) {
122template <
typename DenseVecT,
typename DenseMatT,
typename SparseMatT,
typename NumFactT>
123struct SparseTypesTraitBase {
124 typedef DenseVecT AccelDenseVector;
125 typedef DenseMatT AccelDenseMatrix;
126 typedef SparseMatT AccelSparseMatrix;
128 typedef SparseOpaqueSymbolicFactorization SymbolicFactorization;
129 typedef NumFactT NumericFactorization;
131 typedef AccelFactorizationDeleter<SymbolicFactorization> SymbolicFactorizationDeleter;
132 typedef AccelFactorizationDeleter<NumericFactorization> NumericFactorizationDeleter;
135template <
typename Scalar>
136struct SparseTypesTrait {};
139struct SparseTypesTrait<double> : SparseTypesTraitBase<DenseVector_Double, DenseMatrix_Double, SparseMatrix_Double,
140 SparseOpaqueFactorization_Double> {};
143struct SparseTypesTrait<float>
144 : SparseTypesTraitBase<DenseVector_Float, DenseMatrix_Float, SparseMatrix_Float, SparseOpaqueFactorization_Float> {
149template <
typename MatrixType_,
int UpLo_, SparseFactorization_t Solver_,
bool EnforceSquare_>
150class AccelerateImpl :
public SparseSolverBase<AccelerateImpl<MatrixType_, UpLo_, Solver_, EnforceSquare_> > {
152 using Base = SparseSolverBase<AccelerateImpl>;
154 using Base::m_isInitialized;
157 using Base::_solve_impl;
159 typedef MatrixType_ MatrixType;
160 typedef typename MatrixType::Scalar Scalar;
161 typedef typename MatrixType::StorageIndex StorageIndex;
162 enum { ColsAtCompileTime =
Dynamic, MaxColsAtCompileTime =
Dynamic };
163 enum { UpLo = UpLo_ };
165 using AccelDenseVector =
typename internal::SparseTypesTrait<Scalar>::AccelDenseVector;
166 using AccelDenseMatrix =
typename internal::SparseTypesTrait<Scalar>::AccelDenseMatrix;
167 using AccelSparseMatrix =
typename internal::SparseTypesTrait<Scalar>::AccelSparseMatrix;
168 using SymbolicFactorization =
typename internal::SparseTypesTrait<Scalar>::SymbolicFactorization;
169 using NumericFactorization =
typename internal::SparseTypesTrait<Scalar>::NumericFactorization;
170 using SymbolicFactorizationDeleter =
typename internal::SparseTypesTrait<Scalar>::SymbolicFactorizationDeleter;
171 using NumericFactorizationDeleter =
typename internal::SparseTypesTrait<Scalar>::NumericFactorizationDeleter;
174 m_isInitialized =
false;
176 auto check_flag_set = [](
int value,
int flag) {
return ((value & flag) == flag); };
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;
188 m_sparseKind = SparseTriangular;
189 m_triType = SparseLowerTriangle;
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;
200 m_sparseKind = SparseOrdinary;
201 m_triType = (UpLo_ &
Lower) ? SparseLowerTriangle : SparseUpperTriangle;
204 m_order = SparseOrderDefault;
207 explicit AccelerateImpl(
const MatrixType& matrix) : AccelerateImpl() { compute(matrix); }
211 inline Index cols()
const {
return m_nCols; }
212 inline Index rows()
const {
return m_nRows; }
215 eigen_assert(m_isInitialized &&
"Decomposition is not initialized.");
219 void analyzePattern(
const MatrixType& matrix);
221 void factorize(
const MatrixType& matrix);
223 void compute(
const MatrixType& matrix);
225 template <
typename Rhs,
typename Dest>
226 void _solve_impl(
const MatrixBase<Rhs>& b, MatrixBase<Dest>& dest)
const;
229 void setOrder(SparseOrder_t order) { m_order = order; }
232 template <
typename T>
233 void buildAccelSparseMatrix(
const SparseMatrix<T>& a, AccelSparseMatrix& A, std::vector<long>& columnStarts) {
234 const Index nColumnsStarts = a.cols() + 1;
236 columnStarts.resize(nColumnsStarts);
238 for (
Index i = 0; i < nColumnsStarts; i++) columnStarts[i] = a.outerIndexPtr()[i];
240 SparseAttributes_t attributes{};
241 attributes.transpose =
false;
242 attributes.triangle = m_triType;
243 attributes.kind = m_sparseKind;
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());
253 A.structure = structure;
254 A.data =
const_cast<T*
>(a.valuePtr());
257 void doAnalysis(AccelSparseMatrix& A) {
258 m_numericFactorization.reset(
nullptr);
260 SparseSymbolicFactorOptions opts{};
261 opts.control = SparseDefaultControl;
262 opts.orderMethod = m_order;
263 opts.order =
nullptr;
264 opts.ignoreRowsAndColumns =
nullptr;
265 opts.malloc = malloc;
267 opts.reportError =
nullptr;
269 m_symbolicFactorization.reset(
new SymbolicFactorization(SparseFactor(Solver_, A.structure, opts)));
271 SparseStatus_t status = m_symbolicFactorization->status;
273 updateInfoStatus(status);
275 if (status != SparseStatusOK) m_symbolicFactorization.reset(
nullptr);
278 void doFactorization(AccelSparseMatrix& A) {
279 SparseStatus_t status = SparseStatusReleased;
281 if (m_symbolicFactorization) {
282 m_numericFactorization.reset(
new NumericFactorization(SparseFactor(*m_symbolicFactorization, A)));
284 status = m_numericFactorization->status;
286 if (status != SparseStatusOK) m_numericFactorization.reset(
nullptr);
289 updateInfoStatus(status);
293 void updateInfoStatus(SparseStatus_t status)
const {
298 case SparseFactorizationFailed:
299 case SparseMatrixIsSingular:
302 case SparseInternalError:
303 case SparseParameterError:
304 case SparseStatusReleased:
312 Index m_nRows, m_nCols;
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;
321template <
typename MatrixType_,
int UpLo_, SparseFactorization_t Solver_,
bool EnforceSquare_>
322void AccelerateImpl<MatrixType_, UpLo_, Solver_, EnforceSquare_>::compute(
const MatrixType& a) {
323 if (EnforceSquare_) eigen_assert(a.rows() == a.cols());
328 AccelSparseMatrix A{};
329 std::vector<long> columnStarts;
331 buildAccelSparseMatrix(a, A, columnStarts);
335 if (m_symbolicFactorization) doFactorization(A);
337 m_isInitialized =
true;
346template <
typename MatrixType_,
int UpLo_, SparseFactorization_t Solver_,
bool EnforceSquare_>
347void AccelerateImpl<MatrixType_, UpLo_, Solver_, EnforceSquare_>::analyzePattern(
const MatrixType& a) {
348 if (EnforceSquare_) eigen_assert(a.rows() == a.cols());
353 AccelSparseMatrix A{};
354 std::vector<long> columnStarts;
356 buildAccelSparseMatrix(a, A, columnStarts);
360 m_isInitialized =
true;
370template <
typename MatrixType_,
int UpLo_, SparseFactorization_t Solver_,
bool EnforceSquare_>
371void AccelerateImpl<MatrixType_, UpLo_, Solver_, EnforceSquare_>::factorize(
const MatrixType& a) {
372 eigen_assert(m_symbolicFactorization &&
"You must first call analyzePattern()");
373 eigen_assert(m_nRows == a.rows() && m_nCols == a.cols());
375 if (EnforceSquare_) eigen_assert(a.rows() == a.cols());
377 AccelSparseMatrix A{};
378 std::vector<long> columnStarts;
380 buildAccelSparseMatrix(a, A, columnStarts);
385template <
typename MatrixType_,
int UpLo_, SparseFactorization_t Solver_,
bool EnforceSquare_>
386template <
typename Rhs,
typename Dest>
387void AccelerateImpl<MatrixType_, UpLo_, Solver_, EnforceSquare_>::_solve_impl(
const MatrixBase<Rhs>& b,
389 if (!m_numericFactorization) {
394 eigen_assert(m_nRows == b.rows());
395 eigen_assert(((b.cols() == 1) || b.outerStride() == b.rows()));
397 SparseStatus_t status = SparseStatusOK;
399 Scalar* b_ptr =
const_cast<Scalar*
>(b.derived().data());
400 Scalar* x_ptr =
const_cast<Scalar*
>(x.derived().data());
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;
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;
416 SparseSolve(*m_numericFactorization, bmat, xmat);
418 updateInfoStatus(status);
Base class for all dense matrices, vectors, and expressions.
Definition MatrixBase.h:52
AccelerateImpl< MatrixType, UpLo|Symmetric, SparseFactorizationCholesky, true > AccelerateLLT
A direct Cholesky (LLT) factorization and solver based on Accelerate.
Definition AccelerateSupport.h:25
AccelerateImpl< MatrixType, UpLo|Symmetric, SparseFactorizationLDLTUnpivoted, true > AccelerateLDLTUnpivoted
A direct Cholesky-like LDL^T factorization and solver based on Accelerate with only 1x1 pivots and no...
Definition AccelerateSupport.h:53
AccelerateImpl< MatrixType, UpLo|Symmetric, SparseFactorizationLDLT, true > AccelerateLDLT
The default Cholesky (LDLT) factorization and solver based on Accelerate.
Definition AccelerateSupport.h:39
AccelerateImpl< MatrixType, UpLo|Symmetric, SparseFactorizationLDLTTPP, true > AccelerateLDLTTPP
A direct Cholesky (LDLT) factorization and solver based on Accelerate with full threshold partial piv...
Definition AccelerateSupport.h:82
AccelerateImpl< MatrixType, 0, SparseFactorizationCholeskyAtA, false > AccelerateCholeskyAtA
A QR factorization and solver based on Accelerate without storing Q (equivalent to A^TA = R^T R)
Definition AccelerateSupport.h:108
AccelerateImpl< MatrixType, UpLo|Symmetric, SparseFactorizationLDLTSBK, true > AccelerateLDLTSBK
A direct Cholesky (LDLT) factorization and solver based on Accelerate with Supernode Bunch-Kaufman an...
Definition AccelerateSupport.h:68
AccelerateImpl< MatrixType, 0, SparseFactorizationQR, false > AccelerateQR
A QR factorization and solver based on Accelerate.
Definition AccelerateSupport.h:95
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
Namespace containing all symbols from the Eigen library.
Definition B01_Experimental.dox:1
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:82
const int Dynamic
Definition Constants.h:25