Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ > Class Template Reference

#include <AccelerateSupport.h>

+ Inheritance diagram for Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >:

Public Types

enum  { ColsAtCompileTime = Dynamic , MaxColsAtCompileTime = Dynamic }
 
enum  { UpLo = UpLo_ }
 
typedef MatrixType_ MatrixType
 
typedef MatrixType::Scalar Scalar
 
typedef MatrixType::StorageIndex StorageIndex
 
using AccelDenseVector = typename internal::SparseTypesTrait< Scalar >::AccelDenseVector
 
using AccelDenseMatrix = typename internal::SparseTypesTrait< Scalar >::AccelDenseMatrix
 
using AccelSparseMatrix = typename internal::SparseTypesTrait< Scalar >::AccelSparseMatrix
 
using SymbolicFactorization = typename internal::SparseTypesTrait< Scalar >::SymbolicFactorization
 
using NumericFactorization = typename internal::SparseTypesTrait< Scalar >::NumericFactorization
 
using SymbolicFactorizationDeleter = typename internal::SparseTypesTrait< Scalar >::SymbolicFactorizationDeleter
 
using NumericFactorizationDeleter = typename internal::SparseTypesTrait< Scalar >::NumericFactorizationDeleter
 

Public Member Functions

 AccelerateImpl ()
 
 AccelerateImpl (const MatrixType &matrix)
 
 ~AccelerateImpl ()
 
Index cols () const
 
Index rows () const
 
ComputationInfo info () const
 
void analyzePattern (const MatrixType &matrix)
 
void factorize (const MatrixType &matrix)
 
void compute (const MatrixType &matrix)
 
template<typename Rhs , typename Dest >
void _solve_impl (const MatrixBase< Rhs > &b, MatrixBase< Dest > &dest) const
 
void setOrder (SparseOrder_t order)
 
template<typename Rhs , typename Dest >
void _solve_impl (const SparseMatrixBase< Rhs > &b, SparseMatrixBase< Dest > &dest) const
 
- Public Member Functions inherited from Eigen::SparseSolverBase< AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ > >
 SparseSolverBase ()
 
 SparseSolverBase (SparseSolverBase &&other)
 
 ~SparseSolverBase ()
 
AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ > & derived ()
 
const AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ > & derived () const
 
const Solve< AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >, Rhs > solve (const MatrixBase< Rhs > &b) const
 
const Solve< AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >, Rhs > solve (const SparseMatrixBase< Rhs > &b) const
 
void _solve_impl (const SparseMatrixBase< Rhs > &b, SparseMatrixBase< Dest > &dest) const
 

Protected Types

using Base = SparseSolverBase< AccelerateImpl >
 

Protected Member Functions

void updateInfoStatus (SparseStatus_t status) const
 
Derived & derived ()
 
const Derived & derived () const
 

Protected Attributes

ComputationInfo m_info
 
Index m_nRows
 
Index m_nCols
 
std::unique_ptr< SymbolicFactorization, SymbolicFactorizationDeleterm_symbolicFactorization
 
std::unique_ptr< NumericFactorization, NumericFactorizationDeleterm_numericFactorization
 
SparseKind_t m_sparseKind
 
SparseTriangle_t m_triType
 
SparseOrder_t m_order
 
bool m_isInitialized
 
- Protected Attributes inherited from Eigen::SparseSolverBase< AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ > >
bool m_isInitialized
 

Private Member Functions

template<typename T >
void buildAccelSparseMatrix (const SparseMatrix< T > &a, AccelSparseMatrix &A, std::vector< long > &columnStarts)
 
void doAnalysis (AccelSparseMatrix &A)
 
void doFactorization (AccelSparseMatrix &A)
 

Member Typedef Documentation

◆ AccelDenseMatrix

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
using Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::AccelDenseMatrix = typename internal::SparseTypesTrait<Scalar>::AccelDenseMatrix

◆ AccelDenseVector

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
using Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::AccelDenseVector = typename internal::SparseTypesTrait<Scalar>::AccelDenseVector

◆ AccelSparseMatrix

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
using Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::AccelSparseMatrix = typename internal::SparseTypesTrait<Scalar>::AccelSparseMatrix

◆ Base

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
using Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::Base = SparseSolverBase<AccelerateImpl>
protected

◆ MatrixType

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
typedef MatrixType_ Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::MatrixType

◆ NumericFactorization

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
using Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::NumericFactorization = typename internal::SparseTypesTrait<Scalar>::NumericFactorization

◆ NumericFactorizationDeleter

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
using Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::NumericFactorizationDeleter = typename internal::SparseTypesTrait<Scalar>::NumericFactorizationDeleter

◆ Scalar

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
typedef MatrixType::Scalar Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::Scalar

◆ StorageIndex

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
typedef MatrixType::StorageIndex Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::StorageIndex

◆ SymbolicFactorization

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
using Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::SymbolicFactorization = typename internal::SparseTypesTrait<Scalar>::SymbolicFactorization

◆ SymbolicFactorizationDeleter

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
using Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::SymbolicFactorizationDeleter = typename internal::SparseTypesTrait<Scalar>::SymbolicFactorizationDeleter

Member Enumeration Documentation

◆ anonymous enum

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
anonymous enum
Enumerator
ColsAtCompileTime 
MaxColsAtCompileTime 
@ ColsAtCompileTime
Definition: AccelerateSupport.h:162
@ MaxColsAtCompileTime
Definition: AccelerateSupport.h:162
const int Dynamic
Definition: Constants.h:25

◆ anonymous enum

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
anonymous enum
Enumerator
UpLo 
163 { UpLo = UpLo_ };
@ UpLo
Definition: AccelerateSupport.h:163

Constructor & Destructor Documentation

◆ AccelerateImpl() [1/2]

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::AccelerateImpl ( )
inline
173  {
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  }
SparseOrder_t m_order
Definition: AccelerateSupport.h:317
SparseKind_t m_sparseKind
Definition: AccelerateSupport.h:315
SparseTriangle_t m_triType
Definition: AccelerateSupport.h:316
bool m_isInitialized
Definition: SparseSolverBase.h:110
@ 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
squared absolute value
Definition: GlobalFunctions.h:87

References Eigen::Lower, Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::m_isInitialized, Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::m_order, Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::m_sparseKind, Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::m_triType, Eigen::StrictlyLower, Eigen::StrictlyUpper, Eigen::Symmetric, Eigen::UnitLower, Eigen::UnitUpper, Eigen::Upper, and Eigen::value.

◆ AccelerateImpl() [2/2]

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::AccelerateImpl ( const MatrixType matrix)
inlineexplicit
207 : AccelerateImpl() { compute(matrix); }
AccelerateImpl()
Definition: AccelerateSupport.h:173
void compute(const MatrixType &matrix)
Definition: AccelerateSupport.h:322
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

References Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::compute(), and matrix().

◆ ~AccelerateImpl()

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::~AccelerateImpl ( )
inline
209 {}

Member Function Documentation

◆ _solve_impl() [1/2]

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
template<typename Rhs , typename Dest >
void Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::_solve_impl ( const MatrixBase< Rhs > &  b,
MatrixBase< Dest > &  dest 
) const
388  {
389  if (!m_numericFactorization) {
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 }
#define eigen_assert(x)
Definition: Macros.h:910
Scalar * b
Definition: benchVecAdd.cpp:17
SCALAR Scalar
Definition: bench_gemm.cpp:45
Index m_nRows
Definition: AccelerateSupport.h:312
std::unique_ptr< NumericFactorization, NumericFactorizationDeleter > m_numericFactorization
Definition: AccelerateSupport.h:314
typename internal::SparseTypesTrait< Scalar >::AccelDenseMatrix AccelDenseMatrix
Definition: AccelerateSupport.h:166
ComputationInfo m_info
Definition: AccelerateSupport.h:311
void updateInfoStatus(SparseStatus_t status) const
Definition: AccelerateSupport.h:293
@ InvalidInput
Definition: Constants.h:447
list x
Definition: plotDoE.py:28

References b, eigen_assert, Eigen::InvalidInput, and plotDoE::x.

◆ _solve_impl() [2/2]

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
template<typename Rhs , typename Dest >
void Eigen::SparseSolverBase< Derived >::_solve_impl ( typename Rhs  ,
typename Dest   
)
inline

default implementation of solving with a sparse rhs

104  {
105  internal::solve_sparse_through_dense_panels(derived(), b.derived(), dest.derived());
106  }
Derived & derived()
Definition: SparseSolverBase.h:76
std::enable_if_t< Rhs::ColsAtCompileTime !=1 &&Dest::ColsAtCompileTime !=1 > solve_sparse_through_dense_panels(const Decomposition &dec, const Rhs &rhs, Dest &dest)
Definition: SparseSolverBase.h:25

◆ analyzePattern()

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
void Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::analyzePattern ( const MatrixType a)

Performs a symbolic decomposition on the sparsity pattern of matrix a.

This function is particularly useful when solving for several problems having the same structure.

See also
factorize()
347  {
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 }
Matrix< SCALARA, Dynamic, Dynamic, opt_A > A
Definition: bench_gemm.cpp:47
typename internal::SparseTypesTrait< Scalar >::AccelSparseMatrix AccelSparseMatrix
Definition: AccelerateSupport.h:167
void doAnalysis(AccelSparseMatrix &A)
Definition: AccelerateSupport.h:257
void buildAccelSparseMatrix(const SparseMatrix< T > &a, AccelSparseMatrix &A, std::vector< long > &columnStarts)
Definition: AccelerateSupport.h:233
Index m_nCols
Definition: AccelerateSupport.h:312
const Scalar * a
Definition: level2_cplx_impl.h:32

References a, and eigen_assert.

◆ buildAccelSparseMatrix()

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
template<typename T >
void Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::buildAccelSparseMatrix ( const SparseMatrix< T > &  a,
AccelSparseMatrix A,
std::vector< long > &  columnStarts 
)
inlineprivate
233  {
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  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:83

References a, Eigen::PlainObjectBase< Derived >::data(), i, Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::m_sparseKind, and Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::m_triType.

◆ cols()

◆ compute()

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
void Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::compute ( const MatrixType a)

Computes the symbolic and numeric decomposition of matrix a

322  {
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 
336 
337  m_isInitialized = true;
338 }
std::unique_ptr< SymbolicFactorization, SymbolicFactorizationDeleter > m_symbolicFactorization
Definition: AccelerateSupport.h:313
void doFactorization(AccelSparseMatrix &A)
Definition: AccelerateSupport.h:278

References a, and eigen_assert.

Referenced by Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::AccelerateImpl().

◆ derived() [1/2]

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
Derived& Eigen::SparseSolverBase< Derived >::derived
inlineprotected
76 { return *static_cast<Derived*>(this); }

◆ derived() [2/2]

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
const Derived& Eigen::SparseSolverBase< Derived >::derived
inlineprotected
77 { return *static_cast<const Derived*>(this); }

◆ doAnalysis()

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
void Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::doAnalysis ( AccelSparseMatrix A)
inlineprivate
257  {
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  }
typename internal::SparseTypesTrait< Scalar >::SymbolicFactorization SymbolicFactorization
Definition: AccelerateSupport.h:168

References Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::m_numericFactorization, Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::m_order, Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::m_symbolicFactorization, and Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::updateInfoStatus().

◆ doFactorization()

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
void Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::doFactorization ( AccelSparseMatrix A)
inlineprivate
278  {
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  }
typename internal::SparseTypesTrait< Scalar >::NumericFactorization NumericFactorization
Definition: AccelerateSupport.h:169

References Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::m_numericFactorization, Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::m_symbolicFactorization, and Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::updateInfoStatus().

◆ factorize()

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
void Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::factorize ( const MatrixType a)

Performs a numeric decomposition of matrix a.

The given matrix must have the same sparsity pattern as the matrix on which the symbolic decomposition has been performed.

See also
analyzePattern()
371  {
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 
383 }

References a, and eigen_assert.

◆ info()

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
ComputationInfo Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::info ( ) const
inline

◆ rows()

◆ setOrder()

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
void Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::setOrder ( SparseOrder_t  order)
inline

Sets the ordering algorithm to use.

229 { m_order = order; }

References Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::m_order.

◆ updateInfoStatus()

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
void Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::updateInfoStatus ( SparseStatus_t  status) const
inlineprotected
293  {
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  }
@ NumericalIssue
Definition: Constants.h:442
@ Success
Definition: Constants.h:440

References Eigen::InvalidInput, Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::m_info, Eigen::NumericalIssue, and Eigen::Success.

Referenced by Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::doAnalysis(), and Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::doFactorization().

Member Data Documentation

◆ m_info

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
ComputationInfo Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::m_info
mutableprotected

◆ m_isInitialized

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
bool Eigen::SparseSolverBase< Derived >::m_isInitialized
mutableprotected

◆ m_nCols

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
Index Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::m_nCols
protected

◆ m_nRows

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
Index Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::m_nRows
protected

◆ m_numericFactorization

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
std::unique_ptr<NumericFactorization, NumericFactorizationDeleter> Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::m_numericFactorization
protected

◆ m_order

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
SparseOrder_t Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::m_order
protected

◆ m_sparseKind

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
SparseKind_t Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::m_sparseKind
protected

◆ m_symbolicFactorization

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
std::unique_ptr<SymbolicFactorization, SymbolicFactorizationDeleter> Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::m_symbolicFactorization
protected

◆ m_triType

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
SparseTriangle_t Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::m_triType
protected

The documentation for this class was generated from the following file: