1 #ifndef EIGEN_ACCELERATESUPPORT_H
2 #define EIGEN_ACCELERATESUPPORT_H
4 #include <Accelerate/Accelerate.h>
6 #include <Eigen/Sparse>
10 template <
typename MatrixType_,
int UpLo_, SparseFactorization_t Solver_,
bool EnforceSquare_>
24 template <
typename MatrixType,
int UpLo = Lower>
25 using AccelerateLLT = AccelerateImpl<MatrixType, UpLo | Symmetric, SparseFactorizationCholesky, true>;
38 template <
typename MatrixType,
int UpLo = Lower>
39 using AccelerateLDLT = AccelerateImpl<MatrixType, UpLo | Symmetric, SparseFactorizationLDLT, true>;
52 template <
typename MatrixType,
int UpLo = Lower>
53 using AccelerateLDLTUnpivoted = AccelerateImpl<MatrixType, UpLo | Symmetric, SparseFactorizationLDLTUnpivoted, true>;
66 template <
typename MatrixType,
int UpLo = Lower>
67 using AccelerateLDLTSBK = AccelerateImpl<MatrixType, UpLo | Symmetric, SparseFactorizationLDLTSBK, true>;
80 template <
typename MatrixType,
int UpLo = Lower>
81 using AccelerateLDLTTPP = AccelerateImpl<MatrixType, UpLo | Symmetric, SparseFactorizationLDLTTPP, true>;
93 template <
typename MatrixType>
94 using AccelerateQR = AccelerateImpl<MatrixType, 0, SparseFactorizationQR, false>;
106 template <
typename MatrixType>
110 template <
typename T>
111 struct AccelFactorizationDeleter {
112 void operator()(T* sym) {
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;
127 typedef SparseOpaqueSymbolicFactorization SymbolicFactorization;
128 typedef NumFactT NumericFactorization;
130 typedef AccelFactorizationDeleter<SymbolicFactorization> SymbolicFactorizationDeleter;
131 typedef AccelFactorizationDeleter<NumericFactorization> NumericFactorizationDeleter;
134 template <
typename Scalar>
135 struct SparseTypesTrait {};
138 struct SparseTypesTrait<double> : SparseTypesTraitBase<DenseVector_Double, DenseMatrix_Double, SparseMatrix_Double,
139 SparseOpaqueFactorization_Double> {};
142 struct SparseTypesTrait<float>
143 : SparseTypesTraitBase<DenseVector_Float, DenseMatrix_Float, SparseMatrix_Float, SparseOpaqueFactorization_Float> {
148 template <
typename MatrixType_,
int UpLo_, SparseFactorization_t Solver_,
bool EnforceSquare_>
149 class AccelerateImpl :
public SparseSolverBase<AccelerateImpl<MatrixType_, UpLo_, Solver_, EnforceSquare_> > {
151 using Base = SparseSolverBase<AccelerateImpl>;
153 using Base::m_isInitialized;
156 using Base::_solve_impl;
158 typedef MatrixType_ MatrixType;
159 typedef typename MatrixType::Scalar Scalar;
160 typedef typename MatrixType::StorageIndex StorageIndex;
161 enum { ColsAtCompileTime =
Dynamic, MaxColsAtCompileTime =
Dynamic };
162 enum { UpLo = UpLo_ };
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;
173 m_isInitialized =
false;
175 auto check_flag_set = [](
int value,
int flag) {
return ((value & flag) == flag); };
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;
187 m_sparseKind = SparseTriangular;
188 m_triType = SparseLowerTriangle;
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;
199 m_sparseKind = SparseOrdinary;
200 m_triType = (UpLo_ &
Lower) ? SparseLowerTriangle : SparseUpperTriangle;
203 m_order = SparseOrderDefault;
206 explicit AccelerateImpl(
const MatrixType& matrix) : AccelerateImpl() { compute(matrix); }
210 inline Index cols()
const {
return m_nCols; }
211 inline Index rows()
const {
return m_nRows; }
214 eigen_assert(m_isInitialized &&
"Decomposition is not initialized.");
218 void analyzePattern(
const MatrixType& matrix);
220 void factorize(
const MatrixType& matrix);
222 void compute(
const MatrixType& matrix);
224 template <
typename Rhs,
typename Dest>
225 void _solve_impl(
const MatrixBase<Rhs>& b, MatrixBase<Dest>& dest)
const;
228 void setOrder(SparseOrder_t order) { m_order = order; }
231 template <
typename T>
232 void buildAccelSparseMatrix(
const SparseMatrix<T>& a, AccelSparseMatrix& A, std::vector<long>& columnStarts) {
233 const Index nColumnsStarts = a.cols() + 1;
235 columnStarts.resize(nColumnsStarts);
237 for (
Index i = 0; i < nColumnsStarts; i++) columnStarts[i] = a.outerIndexPtr()[i];
239 SparseAttributes_t attributes{};
240 attributes.transpose =
false;
241 attributes.triangle = m_triType;
242 attributes.kind = m_sparseKind;
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());
252 A.structure = structure;
253 A.data =
const_cast<T*
>(a.valuePtr());
256 void doAnalysis(AccelSparseMatrix& A) {
257 m_numericFactorization.reset(
nullptr);
259 SparseSymbolicFactorOptions opts{};
260 opts.control = SparseDefaultControl;
261 opts.orderMethod = m_order;
262 opts.order =
nullptr;
263 opts.ignoreRowsAndColumns =
nullptr;
264 opts.malloc = malloc;
266 opts.reportError =
nullptr;
268 m_symbolicFactorization.reset(
new SymbolicFactorization(SparseFactor(Solver_, A.structure, opts)));
270 SparseStatus_t status = m_symbolicFactorization->status;
272 updateInfoStatus(status);
274 if (status != SparseStatusOK) m_symbolicFactorization.reset(
nullptr);
277 void doFactorization(AccelSparseMatrix& A) {
278 SparseStatus_t status = SparseStatusReleased;
280 if (m_symbolicFactorization) {
281 m_numericFactorization.reset(
new NumericFactorization(SparseFactor(*m_symbolicFactorization, A)));
283 status = m_numericFactorization->status;
285 if (status != SparseStatusOK) m_numericFactorization.reset(
nullptr);
288 updateInfoStatus(status);
292 void updateInfoStatus(SparseStatus_t status)
const {
297 case SparseFactorizationFailed:
298 case SparseMatrixIsSingular:
301 case SparseInternalError:
302 case SparseParameterError:
303 case SparseStatusReleased:
311 Index m_nRows, m_nCols;
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;
320 template <
typename MatrixType_,
int UpLo_, SparseFactorization_t Solver_,
bool EnforceSquare_>
321 void AccelerateImpl<MatrixType_, UpLo_, Solver_, EnforceSquare_>::compute(
const MatrixType& a) {
322 if (EnforceSquare_) eigen_assert(a.rows() == a.cols());
327 AccelSparseMatrix A{};
328 std::vector<long> columnStarts;
330 buildAccelSparseMatrix(a, A, columnStarts);
334 if (m_symbolicFactorization) doFactorization(A);
336 m_isInitialized =
true;
345 template <
typename MatrixType_,
int UpLo_, SparseFactorization_t Solver_,
bool EnforceSquare_>
346 void AccelerateImpl<MatrixType_, UpLo_, Solver_, EnforceSquare_>::analyzePattern(
const MatrixType& a) {
347 if (EnforceSquare_) eigen_assert(a.rows() == a.cols());
352 AccelSparseMatrix A{};
353 std::vector<long> columnStarts;
355 buildAccelSparseMatrix(a, A, columnStarts);
359 m_isInitialized =
true;
368 template <
typename MatrixType_,
int UpLo_, SparseFactorization_t Solver_,
bool EnforceSquare_>
369 void AccelerateImpl<MatrixType_, UpLo_, Solver_, EnforceSquare_>::factorize(
const MatrixType& a) {
370 eigen_assert(m_symbolicFactorization &&
"You must first call analyzePattern()");
371 eigen_assert(m_nRows == a.rows() && m_nCols == a.cols());
373 if (EnforceSquare_) eigen_assert(a.rows() == a.cols());
375 AccelSparseMatrix A{};
376 std::vector<long> columnStarts;
378 buildAccelSparseMatrix(a, A, columnStarts);
383 template <
typename MatrixType_,
int UpLo_, SparseFactorization_t Solver_,
bool EnforceSquare_>
384 template <
typename Rhs,
typename Dest>
385 void AccelerateImpl<MatrixType_, UpLo_, Solver_, EnforceSquare_>::_solve_impl(
const MatrixBase<Rhs>& b,
386 MatrixBase<Dest>& x)
const {
387 if (!m_numericFactorization) {
392 eigen_assert(m_nRows == b.rows());
393 eigen_assert(((b.cols() == 1) || b.outerStride() == b.rows()));
395 SparseStatus_t status = SparseStatusOK;
397 Scalar* b_ptr =
const_cast<Scalar*
>(b.derived().data());
398 Scalar* x_ptr =
const_cast<Scalar*
>(x.derived().data());
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;
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;
414 SparseSolve(*m_numericFactorization, bmat, xmat);
416 updateInfoStatus(status);
A QR factorization and solver based on Accelerate without storing Q (equivalent to A^TA = R^T R)
A direct Cholesky (LDLT) factorization and solver based on Accelerate with Supernode Bunch-Kaufman an...
A direct Cholesky (LDLT) factorization and solver based on Accelerate with full threshold partial piv...
A direct Cholesky-like LDL^T factorization and solver based on Accelerate with only 1x1 pivots and no...
The default Cholesky (LDLT) factorization and solver based on Accelerate.
A direct Cholesky (LLT) factorization and solver based on Accelerate.
Definition: AccelerateSupport.h:11
A QR factorization and solver based on Accelerate.
ComputationInfo
Definition: Constants.h:442
@ 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:446
@ InvalidInput
Definition: Constants.h:451
@ Success
Definition: Constants.h:444
Namespace containing all symbols from the Eigen library.
Definition: Core:139
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:59
const int Dynamic
Definition: Constants.h:24