#include <AccelerateSupport.h>
|
| 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 |
|
| 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 |
|
◆ AccelDenseMatrix
template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
◆ AccelDenseVector
template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
◆ AccelSparseMatrix
template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
◆ Base
template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
◆ MatrixType
template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
◆ NumericFactorization
template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
◆ NumericFactorizationDeleter
template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
◆ Scalar
template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
◆ StorageIndex
template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
◆ SymbolicFactorization
template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
◆ SymbolicFactorizationDeleter
template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
◆ anonymous enum
template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
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_>
@ UpLo
Definition: AccelerateSupport.h:163
◆ AccelerateImpl() [1/2]
template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
176 auto check_flag_set = [](
int value,
int flag) {
return ((
value & flag) == flag); };
180 m_triType = (UpLo_ &
Lower) ? SparseLowerTriangle : SparseUpperTriangle;
181 }
else if (check_flag_set(UpLo_,
UnitLower)) {
184 }
else if (check_flag_set(UpLo_,
UnitUpper)) {
193 }
else if (check_flag_set(UpLo_,
Lower)) {
196 }
else if (check_flag_set(UpLo_,
Upper)) {
201 m_triType = (UpLo_ &
Lower) ? SparseLowerTriangle : SparseUpperTriangle;
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_>
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_>
◆ _solve_impl() [1/2]
template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
template<typename Rhs , typename Dest >
397 SparseStatus_t status = SparseStatusOK;
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;
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;
#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 >
default implementation of solving with a sparse rhs
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_>
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()
354 std::vector<long> columnStarts;
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 >
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;
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());
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()
template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
◆ compute()
template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
◆ derived() [1/2]
template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
76 {
return *
static_cast<Derived*
>(
this); }
◆ derived() [2/2]
template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
77 {
return *
static_cast<const Derived*
>(
this); }
◆ doAnalysis()
template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
260 SparseSymbolicFactorOptions opts{};
261 opts.control = SparseDefaultControl;
263 opts.order =
nullptr;
264 opts.ignoreRowsAndColumns =
nullptr;
265 opts.malloc = malloc;
267 opts.reportError =
nullptr;
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_>
279 SparseStatus_t status = SparseStatusReleased;
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_>
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()
378 std::vector<long> columnStarts;
References a, and eigen_assert.
◆ info()
template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
◆ rows()
template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
◆ setOrder()
template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
void Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::setOrder |
( |
SparseOrder_t |
order | ) |
|
|
inline |
◆ updateInfoStatus()
template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
void Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::updateInfoStatus |
( |
SparseStatus_t |
status | ) |
const |
|
inlineprotected |
298 case SparseFactorizationFailed:
299 case SparseMatrixIsSingular:
302 case SparseInternalError:
303 case SparseParameterError:
304 case SparseStatusReleased:
@ 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().
◆ m_info
template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
◆ m_isInitialized
template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
◆ m_nCols
template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
◆ m_nRows
template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
◆ m_numericFactorization
template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
◆ m_order
template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
Referenced by Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::AccelerateImpl(), Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::doAnalysis(), and Eigen::AccelerateImpl< MatrixType_, UpLo_, Solver_, EnforceSquare_ >::setOrder().
◆ m_sparseKind
template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
◆ m_symbolicFactorization
template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
◆ m_triType
template<typename MatrixType_ , int UpLo_, SparseFactorization_t Solver_, bool EnforceSquare_>
The documentation for this class was generated from the following file: