12#ifndef EIGEN_SPARSE_LU_H
13#define EIGEN_SPARSE_LU_H
17template <
typename _MatrixType,
typename _OrderingType = COLAMDOrdering<
typename _MatrixType::StorageIndex> >
class SparseLU;
18template <
typename MappedSparseMatrixType>
struct SparseLUMatrixLReturnType;
19template <
typename MatrixLType,
typename MatrixUType>
struct SparseLUMatrixUReturnType;
21template <
bool Conjugate,
class SparseLUType>
22class SparseLUTransposeView :
public SparseSolverBase<SparseLUTransposeView<Conjugate,SparseLUType> >
26 using APIBase::m_isInitialized;
28 typedef typename SparseLUType::Scalar Scalar;
29 typedef typename SparseLUType::StorageIndex StorageIndex;
30 typedef typename SparseLUType::MatrixType MatrixType;
31 typedef typename SparseLUType::OrderingType OrderingType;
34 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
35 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
38 SparseLUTransposeView() : APIBase(), m_sparseLU(NULL) {}
39 SparseLUTransposeView(
const SparseLUTransposeView& view) : APIBase() {
40 this->m_sparseLU = view.m_sparseLU;
41 this->m_isInitialized = view.m_isInitialized;
43 void setIsInitialized(
const bool isInitialized) {this->m_isInitialized = isInitialized;}
44 void setSparseLU(SparseLUType* sparseLU) {m_sparseLU = sparseLU;}
45 using APIBase::_solve_impl;
46 template<
typename Rhs,
typename Dest>
47 bool _solve_impl(
const MatrixBase<Rhs> &B, MatrixBase<Dest> &X_base)
const
49 Dest& X(X_base.derived());
50 eigen_assert(m_sparseLU->info() ==
Success &&
"The matrix should be factorized first");
52 THIS_METHOD_IS_ONLY_FOR_COLUMN_MAJOR_MATRICES);
56 for(
Index j = 0; j < B.cols(); ++j){
57 X.col(j) = m_sparseLU->colsPermutation() * B.const_cast_derived().col(j);
60 m_sparseLU->matrixU().template solveTransposedInPlace<Conjugate>(X);
63 m_sparseLU->matrixL().template solveTransposedInPlace<Conjugate>(X);
66 for (
Index j = 0; j < B.cols(); ++j)
67 X.col(j) = m_sparseLU->rowsPermutation().transpose() * X.col(j);
70 inline Index rows()
const {
return m_sparseLU->rows(); }
71 inline Index cols()
const {
return m_sparseLU->cols(); }
74 SparseLUType *m_sparseLU;
75 SparseLUTransposeView& operator=(
const SparseLUTransposeView&);
131template <
typename _MatrixType,
typename _OrderingType>
136 using APIBase::m_isInitialized;
138 using APIBase::_solve_impl;
140 typedef _MatrixType MatrixType;
141 typedef _OrderingType OrderingType;
142 typedef typename MatrixType::Scalar Scalar;
143 typedef typename MatrixType::RealScalar RealScalar;
144 typedef typename MatrixType::StorageIndex StorageIndex;
153 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
154 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
159 SparseLU():m_lastError(
""),m_Ustore(0,0,0,0,0,0),m_symmetricmode(
false),m_diagpivotthresh(1.0),m_detPermR(1)
163 explicit SparseLU(
const MatrixType& matrix)
164 : m_lastError(
""),m_Ustore(0,0,0,0,0,0),m_symmetricmode(
false),m_diagpivotthresh(1.0),m_detPermR(1)
176 void factorize (
const MatrixType& matrix);
177 void simplicialfactorize(
const MatrixType& matrix);
201 const SparseLUTransposeView<false,SparseLU<_MatrixType,_OrderingType> >
transpose()
203 SparseLUTransposeView<false, SparseLU<_MatrixType,_OrderingType> > transposeView;
204 transposeView.setSparseLU(
this);
205 transposeView.setIsInitialized(this->m_isInitialized);
206 return transposeView;
222 const SparseLUTransposeView<true, SparseLU<_MatrixType,_OrderingType> >
adjoint()
224 SparseLUTransposeView<true, SparseLU<_MatrixType,_OrderingType> > adjointView;
225 adjointView.setSparseLU(
this);
226 adjointView.setIsInitialized(this->m_isInitialized);
230 inline Index rows()
const {
return m_mat.
rows(); }
231 inline Index cols()
const {
return m_mat.
cols(); }
235 m_symmetricmode = sym;
244 SparseLUMatrixLReturnType<SCMatrix>
matrixL()
const
246 return SparseLUMatrixLReturnType<SCMatrix>(m_Lstore);
254 SparseLUMatrixUReturnType<SCMatrix,MappedSparseMatrix<Scalar,ColMajor,StorageIndex> >
matrixU()
const
256 return SparseLUMatrixUReturnType<SCMatrix, MappedSparseMatrix<Scalar,ColMajor,StorageIndex> >(m_Lstore, m_Ustore);
278 m_diagpivotthresh = thresh;
281#ifdef EIGEN_PARSED_BY_DOXYGEN
288 template<
typename Rhs>
302 eigen_assert(m_isInitialized &&
"Decomposition is not initialized.");
314 template<
typename Rhs,
typename Dest>
317 Dest& X(X_base.derived());
318 eigen_assert(m_factorizationIsOk &&
"The matrix should be factorized first");
320 THIS_METHOD_IS_ONLY_FOR_COLUMN_MAJOR_MATRICES);
324 X.resize(B.rows(),B.cols());
327 for(
Index j = 0; j < B.cols(); ++j)
331 this->
matrixL().solveInPlace(X);
332 this->
matrixU().solveInPlace(X);
335 for (
Index j = 0; j < B.cols(); ++j)
354 eigen_assert(m_factorizationIsOk &&
"The matrix should be factorized first.");
356 Scalar det = Scalar(1.);
359 for (
Index j = 0; j < this->cols(); ++j)
365 det *=
abs(it.value());
386 eigen_assert(m_factorizationIsOk &&
"The matrix should be factorized first.");
387 Scalar det = Scalar(0.);
388 for (
Index j = 0; j < this->cols(); ++j)
392 if(it.row() < j)
continue;
395 det +=
log(
abs(it.value()));
409 eigen_assert(m_factorizationIsOk &&
"The matrix should be factorized first.");
414 for (
Index j = 0; j < this->cols(); ++j)
422 else if(it.value()==0)
428 return det * m_detPermR * m_detPermC;
437 eigen_assert(m_factorizationIsOk &&
"The matrix should be factorized first.");
439 Scalar det = Scalar(1.);
442 for (
Index j = 0; j < this->cols(); ++j)
453 return (m_detPermR * m_detPermC) > 0 ? det : -det;
456 Index nnzL()
const {
return m_nnzL; };
457 Index nnzU()
const {
return m_nnzU; };
461 void initperfvalues()
463 m_perfv.panel_size = 16;
465 m_perfv.maxsuper = 128;
468 m_perfv.fillfactor = 20;
473 bool m_factorizationIsOk;
475 std::string m_lastError;
478 MappedSparseMatrix<Scalar,ColMajor,StorageIndex> m_Ustore;
479 PermutationType m_perm_c;
480 PermutationType m_perm_r ;
483 typename Base::GlobalLU_t m_glu;
486 bool m_symmetricmode;
488 internal::perfvalues m_perfv;
489 RealScalar m_diagpivotthresh;
490 Index m_nnzL, m_nnzU;
491 Index m_detPermR, m_detPermC;
494 SparseLU (
const SparseLU& );
510template <
typename MatrixType,
typename OrderingType>
528 ei_declare_aligned_stack_constructed_variable(StorageIndex,outerIndexPtr,mat.cols()+1,mat.isCompressed()?
const_cast<StorageIndex*
>(mat.outerIndexPtr()):0);
531 if(!mat.isCompressed())
535 for (
Index i = 0; i < mat.cols(); i++)
537 m_mat.outerIndexPtr()[m_perm_c.indices()(i)] = outerIndexPtr[i];
538 m_mat.innerNonZeroPtr()[m_perm_c.indices()(i)] = outerIndexPtr[i+1] - outerIndexPtr[i];
543 IndexVector firstRowElt;
544 internal::coletree(m_mat, m_etree,firstRowElt);
547 if (!m_symmetricmode) {
548 IndexVector post, iwork;
550 internal::treePostorder(StorageIndex(m_mat.cols()), m_etree, post);
554 Index m = m_mat.cols();
556 for (
Index i = 0; i < m; ++i) iwork(post(i)) = post(m_etree(i));
560 PermutationType post_perm(m);
561 for (
Index i = 0; i < m; i++)
562 post_perm.
indices()(i) = post(i);
565 if(m_perm_c.size()) {
566 m_perm_c = post_perm * m_perm_c;
571 m_analysisIsOk =
true;
595template <
typename MatrixType,
typename OrderingType>
598 using internal::emptyIdxLU;
599 eigen_assert(m_analysisIsOk &&
"analyzePattern() should be called first");
600 eigen_assert((matrix.rows() == matrix.cols()) &&
"Only for squared matrices");
602 m_isInitialized =
true;
611 const StorageIndex * outerIndexPtr;
612 if (matrix.isCompressed()) outerIndexPtr = matrix.outerIndexPtr();
615 StorageIndex* outerIndexPtr_t =
new StorageIndex[matrix.cols()+1];
616 for(
Index i = 0; i <= matrix.cols(); i++) outerIndexPtr_t[i] = m_mat.outerIndexPtr()[i];
617 outerIndexPtr = outerIndexPtr_t;
619 for (
Index i = 0; i < matrix.cols(); i++)
621 m_mat.outerIndexPtr()[m_perm_c.indices()(i)] = outerIndexPtr[i];
622 m_mat.innerNonZeroPtr()[m_perm_c.indices()(i)] = outerIndexPtr[i+1] - outerIndexPtr[i];
624 if(!matrix.isCompressed())
delete[] outerIndexPtr;
628 m_perm_c.resize(matrix.cols());
629 for(StorageIndex i = 0; i < matrix.cols(); ++i) m_perm_c.indices()(i) = i;
632 Index m = m_mat.rows();
633 Index n = m_mat.cols();
634 Index nnz = m_mat.nonZeros();
635 Index maxpanel = m_perfv.panel_size * m;
641 m_lastError =
"UNABLE TO ALLOCATE WORKING MEMORY\n\n" ;
642 m_factorizationIsOk =
false;
647 IndexVector segrep(m); segrep.
setZero();
648 IndexVector parent(m); parent.
setZero();
649 IndexVector xplore(m); xplore.
setZero();
650 IndexVector repfnz(maxpanel);
651 IndexVector panel_lsub(maxpanel);
652 IndexVector xprune(n); xprune.
setZero();
653 IndexVector marker(m*internal::LUNoMarker); marker.
setZero();
662 tempv.
setZero(internal::LUnumTempV(m, m_perfv.panel_size, m_perfv.maxsuper, m) );
665 PermutationType iperm_c(m_perm_c.inverse());
668 IndexVector relax_end(n);
669 if ( m_symmetricmode ==
true )
676 m_perm_r.indices().setConstant(-1);
680 m_glu.supno(0) = emptyIdxLU; m_glu.xsup.setConstant(0);
681 m_glu.xsup(0) = m_glu.xlsub(0) = m_glu.xusub(0) = m_glu.xlusup(0) =
Index(0);
692 for (jcol = 0; jcol < n; )
695 Index panel_size = m_perfv.panel_size;
696 for (k = jcol + 1; k < (std::min)(jcol+panel_size, n); k++)
698 if (relax_end(k) != emptyIdxLU)
700 panel_size = k - jcol;
705 panel_size = n - jcol;
708 Base::panel_dfs(m, panel_size, jcol, m_mat, m_perm_r.indices(), nseg1, dense, panel_lsub, segrep, repfnz, xprune, marker, parent, xplore, m_glu);
711 Base::panel_bmod(m, panel_size, jcol, nseg1, dense, tempv, segrep, repfnz, m_glu);
714 for ( jj = jcol; jj< jcol + panel_size; jj++)
722 info =
Base::column_dfs(m, jj, m_perm_r.indices(), m_perfv.maxsuper, nseg, panel_lsubk, segrep, repfnz_k, xprune, marker, parent, xplore, m_glu);
725 m_lastError =
"UNABLE TO EXPAND MEMORY IN COLUMN_DFS() ";
727 m_factorizationIsOk =
false;
736 m_lastError =
"UNABLE TO EXPAND MEMORY IN COLUMN_BMOD() ";
738 m_factorizationIsOk =
false;
746 m_lastError =
"UNABLE TO EXPAND MEMORY IN COPY_TO_UCOL() ";
748 m_factorizationIsOk =
false;
756 m_lastError =
"THE MATRIX IS STRUCTURALLY SINGULAR";
758 std::ostringstream returnInfo;
759 returnInfo <<
" ... ZERO COLUMN AT ";
761 m_lastError += returnInfo.str();
764 m_factorizationIsOk =
false;
770 if (pivrow != jj) m_detPermR = -m_detPermR;
773 Base::pruneL(jj, m_perm_r.indices(), pivrow, nseg, segrep, repfnz_k, xprune, m_glu);
776 for (i = 0; i < nseg; i++)
779 repfnz_k(irep) = emptyIdxLU;
785 m_detPermR = m_perm_r.determinant();
786 m_detPermC = m_perm_c.determinant();
794 m_Lstore.setInfos(m, n, m_glu.lusup, m_glu.xlusup, m_glu.lsub, m_glu.xlsub, m_glu.supno, m_glu.xsup);
799 m_factorizationIsOk =
true;
802template<
typename MappedSupernodalType>
803struct SparseLUMatrixLReturnType : internal::no_assignment_operator
805 typedef typename MappedSupernodalType::Scalar
Scalar;
806 explicit SparseLUMatrixLReturnType(
const MappedSupernodalType& mapL) : m_mapL(mapL)
808 Index rows()
const {
return m_mapL.rows(); }
809 Index cols()
const {
return m_mapL.cols(); }
810 template<
typename Dest>
811 void solveInPlace( MatrixBase<Dest> &X)
const
813 m_mapL.solveInPlace(X);
815 template<
bool Conjugate,
typename Dest>
816 void solveTransposedInPlace( MatrixBase<Dest> &X)
const
818 m_mapL.template solveTransposedInPlace<Conjugate>(X);
821 const MappedSupernodalType& m_mapL;
824template<
typename MatrixLType,
typename MatrixUType>
825struct SparseLUMatrixUReturnType : internal::no_assignment_operator
827 typedef typename MatrixLType::Scalar Scalar;
828 SparseLUMatrixUReturnType(
const MatrixLType& mapL,
const MatrixUType& mapU)
829 : m_mapL(mapL),m_mapU(mapU)
831 Index rows()
const {
return m_mapL.rows(); }
832 Index cols()
const {
return m_mapL.cols(); }
834 template<
typename Dest>
void solveInPlace(MatrixBase<Dest> &X)
const
836 Index nrhs = X.cols();
838 for (
Index k = m_mapL.nsuper(); k >= 0; k--)
840 Index fsupc = m_mapL.supToCol()[k];
841 Index lda = m_mapL.colIndexPtr()[fsupc+1] - m_mapL.colIndexPtr()[fsupc];
842 Index nsupc = m_mapL.supToCol()[k+1] - fsupc;
843 Index luptr = m_mapL.colIndexPtr()[fsupc];
847 for (
Index j = 0; j < nrhs; j++)
849 X(fsupc, j) /= m_mapL.valuePtr()[luptr];
855 Map<const Matrix<Scalar,Dynamic,Dynamic, ColMajor>, 0, OuterStride<> > A( &(m_mapL.valuePtr()[luptr]), nsupc, nsupc, OuterStride<>(lda) );
856 typename Dest::RowsBlockXpr U = X.derived().middleRows(fsupc, nsupc);
857 U = A.template triangularView<Upper>().solve(U);
860 for (
Index j = 0; j < nrhs; ++j)
862 for (
Index jcol = fsupc; jcol < fsupc + nsupc; jcol++)
864 typename MatrixUType::InnerIterator it(m_mapU, jcol);
867 Index irow = it.index();
868 X(irow, j) -= X(jcol, j) * it.value();
875 template<
bool Conjugate,
typename Dest>
void solveTransposedInPlace(MatrixBase<Dest> &X)
const
878 Index nrhs = X.cols();
880 for (
Index k = 0; k <= m_mapL.nsuper(); k++)
882 Index fsupc = m_mapL.supToCol()[k];
883 Index lda = m_mapL.colIndexPtr()[fsupc+1] - m_mapL.colIndexPtr()[fsupc];
884 Index nsupc = m_mapL.supToCol()[k+1] - fsupc;
885 Index luptr = m_mapL.colIndexPtr()[fsupc];
887 for (
Index j = 0; j < nrhs; ++j)
889 for (
Index jcol = fsupc; jcol < fsupc + nsupc; jcol++)
891 typename MatrixUType::InnerIterator it(m_mapU, jcol);
894 Index irow = it.index();
895 X(jcol, j) -= X(irow, j) * (Conjugate?
conj(it.value()): it.value());
901 for (
Index j = 0; j < nrhs; j++)
903 X(fsupc, j) /= (Conjugate?
conj(m_mapL.valuePtr()[luptr]) : m_mapL.valuePtr()[luptr]);
908 Map<const Matrix<Scalar,Dynamic,Dynamic, ColMajor>, 0, OuterStride<> > A( &(m_mapL.valuePtr()[luptr]), nsupc, nsupc, OuterStride<>(lda) );
909 typename Dest::RowsBlockXpr U = X.derived().middleRows(fsupc, nsupc);
911 U = A.adjoint().template triangularView<Lower>().solve(U);
913 U = A.transpose().template triangularView<Lower>().solve(U);
919 const MatrixLType& m_mapL;
920 const MatrixUType& m_mapU;
ColXpr col(Index i)
Definition DenseBase.h:1098
Sparse matrix.
Definition MappedSparseMatrix.h:34
Base class for all dense matrices, vectors, and expressions.
Definition MatrixBase.h:50
The matrix class, also used for vectors and row-vectors.
Definition Matrix.h:180
InverseReturnType inverse() const
Definition PermutationMatrix.h:185
Permutation matrix.
Definition PermutationMatrix.h:298
const IndicesType & indices() const
Definition PermutationMatrix.h:360
friend class Eigen::Map
Definition PlainObjectBase.h:987
void resize(Index rows, Index cols)
Definition PlainObjectBase.h:271
Derived & setConstant(Index size, const Scalar &val)
Definition CwiseNullaryOp.h:361
Derived & setZero(Index size)
Definition CwiseNullaryOp.h:562
Pseudo expression representing a solving operation.
Definition Solve.h:63
Sparse supernodal LU factorization for general matrices.
Definition SparseLU.h:133
Scalar determinant()
Definition SparseLU.h:435
Scalar absDeterminant()
Definition SparseLU.h:351
const PermutationType & colsPermutation() const
Definition SparseLU.h:271
const SparseLUTransposeView< false, SparseLU< _MatrixType, _OrderingType > > transpose()
Definition SparseLU.h:201
void factorize(const MatrixType &matrix)
Definition SparseLU.h:596
SparseLUMatrixLReturnType< SCMatrix > matrixL() const
Definition SparseLU.h:244
std::string lastErrorMessage() const
Definition SparseLU.h:309
Scalar signDeterminant()
Definition SparseLU.h:407
const Solve< SparseLU, Rhs > solve(const MatrixBase< Rhs > &B) const
Scalar logAbsDeterminant() const
Definition SparseLU.h:381
void setPivotThreshold(const RealScalar &thresh)
Definition SparseLU.h:276
void compute(const MatrixType &matrix)
Definition SparseLU.h:183
void analyzePattern(const MatrixType &matrix)
Definition SparseLU.h:511
const SparseLUTransposeView< true, SparseLU< _MatrixType, _OrderingType > > adjoint()
Definition SparseLU.h:222
ComputationInfo info() const
Reports whether previous computation was successful.
Definition SparseLU.h:300
const PermutationType & rowsPermutation() const
Definition SparseLU.h:263
SparseLUMatrixUReturnType< SCMatrix, MappedSparseMatrix< Scalar, ColMajor, StorageIndex > > matrixU() const
Definition SparseLU.h:254
void isSymmetric(bool sym)
Definition SparseLU.h:233
A versatible sparse matrix representation.
Definition SparseMatrix.h:98
Index rows() const
Definition SparseMatrix.h:138
Index cols() const
Definition SparseMatrix.h:140
A base class for sparse solvers.
Definition SparseSolverBase.h:68
SparseSolverBase()
Definition SparseSolverBase.h:72
Expression of a fixed-size or dynamic-size sub-vector.
Definition VectorBlock.h:60
a class to manipulate the L supernodal factor from the SparseLU factorization
Definition SparseLU_SupernodalMatrix.h:34
Definition SparseLUImpl.h:21
void panel_bmod(const Index m, const Index w, const Index jcol, const Index nseg, ScalarVector &dense, ScalarVector &tempv, IndexVector &segrep, IndexVector &repfnz, GlobalLU_t &glu)
Performs numeric block updates (sup-panel) in topological order.
Definition SparseLU_panel_bmod.h:56
void relax_snode(const Index n, IndexVector &et, const Index relax_columns, IndexVector &descendants, IndexVector &relax_end)
Identify the initial relaxed supernodes.
Definition SparseLU_relax_snode.h:47
void pruneL(const Index jcol, const IndexVector &perm_r, const Index pivrow, const Index nseg, const IndexVector &segrep, BlockIndexVector repfnz, IndexVector &xprune, GlobalLU_t &glu)
Prunes the L-structure.
Definition SparseLU_pruneL.h:53
Index column_dfs(const Index m, const Index jcol, IndexVector &perm_r, Index maxsuper, Index &nseg, BlockIndexVector lsub_col, IndexVector &segrep, BlockIndexVector repfnz, IndexVector &xprune, IndexVector &marker, IndexVector &parent, IndexVector &xplore, GlobalLU_t &glu)
Performs a symbolic factorization on column jcol and decide the supernode boundary.
Definition SparseLU_column_dfs.h:93
void heap_relax_snode(const Index n, IndexVector &et, const Index relax_columns, IndexVector &descendants, IndexVector &relax_end)
Identify the initial relaxed supernodes.
Definition SparseLU_heap_relax_snode.h:46
Index pivotL(const Index jcol, const RealScalar &diagpivotthresh, IndexVector &perm_r, IndexVector &iperm_c, Index &pivrow, GlobalLU_t &glu)
Performs the numerical pivotin on the current column of L, and the CDIV operation.
Definition SparseLU_pivotL.h:60
Index memInit(Index m, Index n, Index annz, Index lwork, Index fillratio, Index panel_size, GlobalLU_t &glu)
Allocate various working space for the numerical factorization phase.
Definition SparseLU_Memory.h:151
void panel_dfs(const Index m, const Index w, const Index jcol, MatrixType &A, IndexVector &perm_r, Index &nseg, ScalarVector &dense, IndexVector &panel_lsub, IndexVector &segrep, IndexVector &repfnz, IndexVector &xprune, IndexVector &marker, IndexVector &parent, IndexVector &xplore, GlobalLU_t &glu)
Performs a symbolic factorization on a panel of columns [jcol, jcol+w)
Definition SparseLU_panel_dfs.h:219
void countnz(const Index n, Index &nnzL, Index &nnzU, GlobalLU_t &glu)
Count Nonzero elements in the factors.
Definition SparseLU_Utils.h:21
void fixupL(const Index n, const IndexVector &perm_r, GlobalLU_t &glu)
Fix up the data storage lsub for L-subscripts.
Definition SparseLU_Utils.h:52
Index column_bmod(const Index jcol, const Index nseg, BlockScalarVector dense, ScalarVector &tempv, BlockIndexVector segrep, BlockIndexVector repfnz, Index fpanelc, GlobalLU_t &glu)
Performs numeric block updates (sup-col) in topological order.
Definition SparseLU_column_bmod.h:53
Index copy_to_ucol(const Index jcol, const Index nseg, IndexVector &segrep, BlockIndexVector repfnz, IndexVector &perm_r, BlockScalarVector dense, GlobalLU_t &glu)
Performs numeric block updates (sup-col) in topological order.
Definition SparseLU_copy_to_ucol.h:50
ComputationInfo
Definition Constants.h:440
@ NumericalIssue
Definition Constants.h:444
@ Success
Definition Constants.h:442
const unsigned int RowMajorBit
Definition Constants.h:66
Namespace containing all symbols from the Eigen library.
Definition B01_Experimental.dox:1
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_conjugate_op< typename Derived::Scalar >, const Derived > conj(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:74
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_log_op< typename Derived::Scalar >, const Derived > log(const Eigen::ArrayBase< Derived > &x)