Eigen  5.0.1-dev+60122df6
 
Loading...
Searching...
No Matches
AccelerateSupport.h
1#ifndef EIGEN_ACCELERATESUPPORT_H
2#define EIGEN_ACCELERATESUPPORT_H
3
4#include <Accelerate/Accelerate.h>
5
6#include <Eigen/Sparse>
7
8namespace Eigen {
9
10template <typename MatrixType_, int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
11class AccelerateImpl;
12
24template <typename MatrixType, int UpLo = Lower>
25using AccelerateLLT = AccelerateImpl<MatrixType, UpLo | Symmetric, SparseFactorizationCholesky, true>;
26
38template <typename MatrixType, int UpLo = Lower>
39using AccelerateLDLT = AccelerateImpl<MatrixType, UpLo | Symmetric, SparseFactorizationLDLT, true>;
40
52template <typename MatrixType, int UpLo = Lower>
53using AccelerateLDLTUnpivoted = AccelerateImpl<MatrixType, UpLo | Symmetric, SparseFactorizationLDLTUnpivoted, true>;
54
67template <typename MatrixType, int UpLo = Lower>
68using AccelerateLDLTSBK = AccelerateImpl<MatrixType, UpLo | Symmetric, SparseFactorizationLDLTSBK, true>;
69
81template <typename MatrixType, int UpLo = Lower>
82using AccelerateLDLTTPP = AccelerateImpl<MatrixType, UpLo | Symmetric, SparseFactorizationLDLTTPP, true>;
83
94template <typename MatrixType>
95using AccelerateQR = AccelerateImpl<MatrixType, 0, SparseFactorizationQR, false>;
96
107template <typename MatrixType>
108using AccelerateCholeskyAtA = AccelerateImpl<MatrixType, 0, SparseFactorizationCholeskyAtA, false>;
109
110namespace internal {
111template <typename T>
112struct AccelFactorizationDeleter {
113 void operator()(T* sym) {
114 if (sym) {
115 SparseCleanup(*sym);
116 delete sym;
117 sym = nullptr;
118 }
119 }
120};
121
122template <typename DenseVecT, typename DenseMatT, typename SparseMatT, typename NumFactT>
123struct SparseTypesTraitBase {
124 typedef DenseVecT AccelDenseVector;
125 typedef DenseMatT AccelDenseMatrix;
126 typedef SparseMatT AccelSparseMatrix;
127
128 typedef SparseOpaqueSymbolicFactorization SymbolicFactorization;
129 typedef NumFactT NumericFactorization;
130
131 typedef AccelFactorizationDeleter<SymbolicFactorization> SymbolicFactorizationDeleter;
132 typedef AccelFactorizationDeleter<NumericFactorization> NumericFactorizationDeleter;
133};
134
135template <typename Scalar>
136struct SparseTypesTrait {};
137
138template <>
139struct SparseTypesTrait<double> : SparseTypesTraitBase<DenseVector_Double, DenseMatrix_Double, SparseMatrix_Double,
140 SparseOpaqueFactorization_Double> {};
141
142template <>
143struct SparseTypesTrait<float>
144 : SparseTypesTraitBase<DenseVector_Float, DenseMatrix_Float, SparseMatrix_Float, SparseOpaqueFactorization_Float> {
145};
146
147} // end namespace internal
148
149template <typename MatrixType_, int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
150class AccelerateImpl : public SparseSolverBase<AccelerateImpl<MatrixType_, UpLo_, Solver_, EnforceSquare_> > {
151 protected:
152 using Base = SparseSolverBase<AccelerateImpl>;
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;
162 enum { ColsAtCompileTime = Dynamic, MaxColsAtCompileTime = Dynamic };
163 enum { UpLo = UpLo_ };
164
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;
172
173 AccelerateImpl() {
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
207 explicit AccelerateImpl(const MatrixType& matrix) : AccelerateImpl() { compute(matrix); }
208
209 ~AccelerateImpl() {}
210
211 inline Index cols() const { return m_nCols; }
212 inline Index rows() const { return m_nRows; }
213
214 ComputationInfo info() const {
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
257 void doAnalysis(AccelSparseMatrix& A) {
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
278 void doFactorization(AccelSparseMatrix& A) {
279 SparseStatus_t status = SparseStatusReleased;
280
281 if (m_symbolicFactorization) {
282 m_numericFactorization.reset(new NumericFactorization(SparseFactor(*m_symbolicFactorization, A)));
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:
300 m_info = NumericalIssue;
301 break;
302 case SparseInternalError:
303 case SparseParameterError:
304 case SparseStatusReleased:
305 default:
306 m_info = InvalidInput;
307 break;
308 }
309 }
310
311 mutable ComputationInfo m_info;
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;
318};
319
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());
324
325 m_nRows = a.rows();
326 m_nCols = a.cols();
327
328 AccelSparseMatrix A{};
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
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());
349
350 m_nRows = a.rows();
351 m_nCols = a.cols();
352
353 AccelSparseMatrix A{};
354 std::vector<long> columnStarts;
355
356 buildAccelSparseMatrix(a, A, columnStarts);
357
358 doAnalysis(A);
359
360 m_isInitialized = true;
361}
362
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());
374
375 if (EnforceSquare_) eigen_assert(a.rows() == a.cols());
376
377 AccelSparseMatrix A{};
378 std::vector<long> columnStarts;
379
380 buildAccelSparseMatrix(a, A, columnStarts);
381
382 doFactorization(A);
383}
384
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,
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
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