11 #ifndef EIGEN_SPARSELU_SUPERNODAL_MATRIX_H
12 #define EIGEN_SPARSELU_SUPERNODAL_MATRIX_H
35 template <
typename Scalar_,
typename StorageIndex_>
47 setInfos(
m,
n, nzval, nzval_colptr, rowind, rowind_colptr, col_to_sup, sup_to_col);
128 template <
typename Dest>
130 template <
bool Conjugate,
typename Dest>
151 template <
typename Scalar,
typename StorageIndex>
159 m_startidval(m_idval),
168 inline Scalar value()
const {
return m_matrix.valuePtr()[m_idval]; }
172 inline Index index()
const {
return m_matrix.rowIndex()[m_idrow]; }
179 return ((m_idval < m_endidval) && (m_idval >= m_startidval) && (m_idrow < m_endidrow));
197 template <
typename Scalar,
typename Index_>
198 template <
typename Dest>
205 const Scalar* Lval = valuePtr();
208 for (
Index k = 0;
k <= nsuper();
k++) {
209 Index fsupc = supToCol()[
k];
210 Index istart = rowIndexPtr()[fsupc];
211 Index nsupr = rowIndexPtr()[fsupc + 1] - istart;
212 Index nsupc = supToCol()[
k + 1] - fsupc;
213 Index nrow = nsupr - nsupc;
227 Index luptr = colIndexPtr()[fsupc];
228 Index lda = colIndexPtr()[fsupc + 1] - luptr;
233 typename Dest::RowsBlockXpr
U =
X.derived().middleRows(fsupc, nsupc);
234 U =
A.template triangularView<UnitLower>().solve(
U);
238 work.topRows(nrow).noalias() =
A *
U;
242 Index iptr = istart + nsupc;
244 irow = rowIndex()[iptr];
245 X(irow,
j) -= work(
i,
j);
254 template <
typename Scalar,
typename Index_>
255 template <
bool Conjugate,
typename Dest>
260 const Scalar* Lval = valuePtr();
263 for (
Index k = nsuper();
k >= 0;
k--) {
264 Index fsupc = supToCol()[
k];
265 Index istart = rowIndexPtr()[fsupc];
266 Index nsupr = rowIndexPtr()[fsupc + 1] - istart;
267 Index nsupc = supToCol()[
k + 1] - fsupc;
268 Index nrow = nsupr - nsupc;
282 Index luptr = colIndexPtr()[fsupc];
283 Index lda = colIndexPtr()[fsupc + 1] - luptr;
287 Index iptr = istart + nsupc;
289 irow = rowIndex()[iptr];
290 work.topRows(nrow)(
i,
j) =
X(irow,
j);
298 typename Dest::RowsBlockXpr
U =
X.derived().middleRows(fsupc, nsupc);
300 U =
U -
A.adjoint() * work.topRows(nrow);
302 U =
U -
A.transpose() * work.topRows(nrow);
308 U =
A.adjoint().template triangularView<UnitUpper>().solve(
U);
310 U =
A.transpose().template triangularView<UnitUpper>().solve(
U);
AnnoyingScalar conj(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:133
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
SCALAR Scalar
Definition: bench_gemm.cpp:45
Matrix< SCALARA, Dynamic, Dynamic, opt_A > A
Definition: bench_gemm.cpp:47
Definition: ForwardDeclarations.h:102
An InnerIterator allows to loop over the element of any matrix expression.
Definition: CoreIterators.h:37
A matrix or vector expression mapping an existing array of data.
Definition: Map.h:96
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
Convenience specialization of Stride to specify only an outer stride See class Map for some examples.
Definition: Stride.h:104
constexpr EIGEN_DEVICE_FUNC const Scalar * data() const
Definition: PlainObjectBase.h:273
EIGEN_DEVICE_FUNC Derived & setZero(Index size)
Definition: CwiseNullaryOp.h:569
InnerIterator class to iterate over nonzero values of the current column in the supernodal matrix L.
Definition: SparseLU_SupernodalMatrix.h:152
const MappedSuperNodalMatrix & m_matrix
Definition: SparseLU_SupernodalMatrix.h:183
Index m_idrow
Definition: SparseLU_SupernodalMatrix.h:189
Index row() const
Definition: SparseLU_SupernodalMatrix.h:173
Scalar value() const
Definition: SparseLU_SupernodalMatrix.h:168
const Index m_outer
Definition: SparseLU_SupernodalMatrix.h:184
InnerIterator & operator++()
Definition: SparseLU_SupernodalMatrix.h:163
Index index() const
Definition: SparseLU_SupernodalMatrix.h:172
Scalar & valueRef()
Definition: SparseLU_SupernodalMatrix.h:170
Index m_idval
Definition: SparseLU_SupernodalMatrix.h:186
Index supIndex() const
Definition: SparseLU_SupernodalMatrix.h:176
Index m_endidrow
Definition: SparseLU_SupernodalMatrix.h:190
const Index m_startidval
Definition: SparseLU_SupernodalMatrix.h:187
const Index m_supno
Definition: SparseLU_SupernodalMatrix.h:185
Index col() const
Definition: SparseLU_SupernodalMatrix.h:174
const Index m_endidval
Definition: SparseLU_SupernodalMatrix.h:188
InnerIterator(const MappedSuperNodalMatrix &mat, Index outer)
Definition: SparseLU_SupernodalMatrix.h:154
a class to manipulate the L supernodal factor from the SparseLU factorization
Definition: SparseLU_SupernodalMatrix.h:36
Matrix< Scalar, Dynamic, 1 > ScalarVector
Definition: SparseLU_SupernodalMatrix.h:41
StorageIndex * m_col_to_sup
Definition: SparseLU_SupernodalMatrix.h:141
StorageIndex * colIndexPtr()
Definition: SparseLU_SupernodalMatrix.h:91
StorageIndex * rowIndex()
Definition: SparseLU_SupernodalMatrix.h:98
const StorageIndex * colIndexPtr() const
Definition: SparseLU_SupernodalMatrix.h:93
Scalar * m_nzval
Definition: SparseLU_SupernodalMatrix.h:137
StorageIndex * supToCol()
Definition: SparseLU_SupernodalMatrix.h:118
void setInfos(Index m, Index n, ScalarVector &nzval, IndexVector &nzval_colptr, IndexVector &rowind, IndexVector &rowind_colptr, IndexVector &col_to_sup, IndexVector &sup_to_col)
Definition: SparseLU_SupernodalMatrix.h:57
Scalar_ Scalar
Definition: SparseLU_SupernodalMatrix.h:38
const StorageIndex * colToSup() const
Definition: SparseLU_SupernodalMatrix.h:114
StorageIndex * colToSup()
Definition: SparseLU_SupernodalMatrix.h:112
const StorageIndex * rowIndexPtr() const
Definition: SparseLU_SupernodalMatrix.h:107
Index nsuper() const
Definition: SparseLU_SupernodalMatrix.h:125
StorageIndex * m_sup_to_col
Definition: SparseLU_SupernodalMatrix.h:142
MappedSuperNodalMatrix(Index m, Index n, ScalarVector &nzval, IndexVector &nzval_colptr, IndexVector &rowind, IndexVector &rowind_colptr, IndexVector &col_to_sup, IndexVector &sup_to_col)
Definition: SparseLU_SupernodalMatrix.h:45
Index m_row
Definition: SparseLU_SupernodalMatrix.h:134
const StorageIndex * supToCol() const
Definition: SparseLU_SupernodalMatrix.h:120
StorageIndex * m_rowind_colptr
Definition: SparseLU_SupernodalMatrix.h:140
void solveTransposedInPlace(MatrixBase< Dest > &X) const
Definition: SparseLU_SupernodalMatrix.h:256
Index m_col
Definition: SparseLU_SupernodalMatrix.h:135
const Scalar * valuePtr() const
Definition: SparseLU_SupernodalMatrix.h:87
Scalar * valuePtr()
Definition: SparseLU_SupernodalMatrix.h:85
MappedSuperNodalMatrix()
Definition: SparseLU_SupernodalMatrix.h:44
const StorageIndex * rowIndex() const
Definition: SparseLU_SupernodalMatrix.h:100
StorageIndex * rowIndexPtr()
Definition: SparseLU_SupernodalMatrix.h:105
void solveInPlace(MatrixBase< Dest > &X) const
Solve with the supernode triangular matrix.
Definition: SparseLU_SupernodalMatrix.h:199
Matrix< StorageIndex, Dynamic, 1 > IndexVector
Definition: SparseLU_SupernodalMatrix.h:40
StorageIndex * m_nzval_colptr
Definition: SparseLU_SupernodalMatrix.h:138
StorageIndex_ StorageIndex
Definition: SparseLU_SupernodalMatrix.h:39
Index rows() const
Definition: SparseLU_SupernodalMatrix.h:73
StorageIndex * m_rowind
Definition: SparseLU_SupernodalMatrix.h:139
Index m_nsuper
Definition: SparseLU_SupernodalMatrix.h:136
Index cols() const
Definition: SparseLU_SupernodalMatrix.h:78
~MappedSuperNodalMatrix()
Definition: SparseLU_SupernodalMatrix.h:50
#define X
Definition: icosphere.cpp:20
const char const int const RealScalar const RealScalar const int * lda
Definition: level2_cplx_impl.h:20
int * m
Definition: level2_cplx_impl.h:294
char char char int int * k
Definition: level2_impl.h:374
Namespace containing all symbols from the Eigen library.
Definition: bench_norm.cpp:70
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:83
const AutoDiffScalar< DerType > & conj(const AutoDiffScalar< DerType > &x)
Definition: AutoDiffScalar.h:482
double U
Swimming speed.
Definition: two_d_variable_diff_adapt.cc:53
Definition: Eigen_Colamd.h:49
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2