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

Public Types

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

Public Member Functions

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

Protected Types

using Base = SparseSolverBase< AccelerateImpl >
 

Protected Member Functions

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

Protected Attributes

ComputationInfo m_info
 
bool m_isInitialized
 
Index m_nCols
 
Index m_nRows
 
std::unique_ptr< NumericFactorization, NumericFactorizationDeleterm_numericFactorization
 
SparseOrder_t m_order
 
SparseKind_t m_sparseKind
 
std::unique_ptr< SymbolicFactorization, SymbolicFactorizationDeleterm_symbolicFactorization
 
SparseTriangle_t m_triType
 
- 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)
 

Detailed Description

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

Definition at line 149 of file AccelerateSupport.h.

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

Definition at line 165 of file AccelerateSupport.h.

◆ AccelDenseVector

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

Definition at line 164 of file AccelerateSupport.h.

◆ AccelSparseMatrix

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

Definition at line 166 of file AccelerateSupport.h.

◆ Base

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

Definition at line 151 of file AccelerateSupport.h.

◆ MatrixType

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

Definition at line 158 of file AccelerateSupport.h.

◆ NumericFactorization

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

Definition at line 168 of file AccelerateSupport.h.

◆ NumericFactorizationDeleter

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

Definition at line 170 of file AccelerateSupport.h.

◆ Scalar

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

Definition at line 159 of file AccelerateSupport.h.

◆ StorageIndex

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

Definition at line 160 of file AccelerateSupport.h.

◆ SymbolicFactorization

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

Definition at line 167 of file AccelerateSupport.h.

◆ SymbolicFactorizationDeleter

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

Definition at line 169 of file AccelerateSupport.h.

Member Enumeration Documentation

◆ anonymous enum

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
anonymous enum
Enumerator
ColsAtCompileTime 
MaxColsAtCompileTime 

Definition at line 161 of file AccelerateSupport.h.

◆ anonymous enum

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
anonymous enum
Enumerator
UpLo 

Definition at line 162 of file AccelerateSupport.h.

162 { UpLo = UpLo_ };

Constructor & Destructor Documentation

◆ AccelerateImpl() [1/2]

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

Definition at line 172 of file AccelerateSupport.h.

172  {
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  }
SparseTriangle_t m_triType
@ 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

◆ AccelerateImpl() [2/2]

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::AccelerateImpl ( const MatrixType matrix)
inlineexplicit

Definition at line 206 of file AccelerateSupport.h.

206 : AccelerateImpl() { compute(matrix); }
void compute(const MatrixType &matrix)

◆ ~AccelerateImpl()

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

Definition at line 208 of file AccelerateSupport.h.

208 {}

Member Function Documentation

◆ _solve_impl()

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

Definition at line 385 of file AccelerateSupport.h.

386  {
387  if (!m_numericFactorization) {
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 }
Array< int, 3, 1 > b
#define eigen_assert(x)
Definition: Macros.h:902
std::unique_ptr< NumericFactorization, NumericFactorizationDeleter > m_numericFactorization
typename internal::SparseTypesTrait< Scalar >::AccelDenseMatrix AccelDenseMatrix
MatrixType::Scalar Scalar
void updateInfoStatus(SparseStatus_t status) const
@ InvalidInput
Definition: Constants.h:453

◆ 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()

Definition at line 346 of file AccelerateSupport.h.

346  {
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 }
MatrixXcf A
typename internal::SparseTypesTrait< Scalar >::AccelSparseMatrix AccelSparseMatrix
void doAnalysis(AccelSparseMatrix &A)
void buildAccelSparseMatrix(const SparseMatrix< T > &a, AccelSparseMatrix &A, std::vector< long > &columnStarts)

◆ 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

Definition at line 232 of file AccelerateSupport.h.

232  {
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  }
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82

◆ cols()

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

Definition at line 210 of file AccelerateSupport.h.

210 { return m_nCols; }

◆ 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

Definition at line 321 of file AccelerateSupport.h.

321  {
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 
335 
336  m_isInitialized = true;
337 }
std::unique_ptr< SymbolicFactorization, SymbolicFactorizationDeleter > m_symbolicFactorization
void doFactorization(AccelSparseMatrix &A)

◆ derived() [1/2]

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
Derived& Eigen::SparseSolverBase< Derived >::derived
inlineprotected

Definition at line 83 of file SparseSolverBase.h.

83 { 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

Definition at line 84 of file SparseSolverBase.h.

84 { 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

Definition at line 256 of file AccelerateSupport.h.

256  {
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  }
typename internal::SparseTypesTrait< Scalar >::SymbolicFactorization SymbolicFactorization

◆ doFactorization()

template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
void Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::doFactorization ( AccelSparseMatrix A)
inlineprivate

Definition at line 277 of file AccelerateSupport.h.

277  {
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  }
typename internal::SparseTypesTrait< Scalar >::NumericFactorization NumericFactorization

◆ 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()

Definition at line 369 of file AccelerateSupport.h.

369  {
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 
381 }

◆ info()

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

Definition at line 213 of file AccelerateSupport.h.

213  {
214  eigen_assert(m_isInitialized && "Decomposition is not initialized.");
215  return m_info;
216  }

◆ rows()

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

Definition at line 211 of file AccelerateSupport.h.

211 { return m_nRows; }

◆ 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.

Definition at line 228 of file AccelerateSupport.h.

228 { m_order = 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

Definition at line 292 of file AccelerateSupport.h.

292  {
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  }
@ NumericalIssue
Definition: Constants.h:448
@ Success
Definition: Constants.h:446

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

Definition at line 310 of file AccelerateSupport.h.

◆ m_isInitialized

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

Definition at line 123 of file SparseSolverBase.h.

◆ m_nCols

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

Definition at line 311 of file AccelerateSupport.h.

◆ m_nRows

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

Definition at line 311 of file AccelerateSupport.h.

◆ 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

Definition at line 313 of file AccelerateSupport.h.

◆ m_order

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

Definition at line 316 of file AccelerateSupport.h.

◆ m_sparseKind

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

Definition at line 314 of file AccelerateSupport.h.

◆ 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

Definition at line 312 of file AccelerateSupport.h.

◆ m_triType

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

Definition at line 315 of file AccelerateSupport.h.


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