aboutsummaryrefslogtreecommitdiffhomepage
path: root/Eigen/src/Core
diff options
context:
space:
mode:
authorGravatar Jitse Niesen <jitse@maths.leeds.ac.uk>2011-05-03 17:08:14 +0100
committerGravatar Jitse Niesen <jitse@maths.leeds.ac.uk>2011-05-03 17:08:14 +0100
commita96c849c20cd787c40c2cbf30a496fed7bf2bf1e (patch)
treeafb8853eb08908695e45b3d9ed237557ae2ac571 /Eigen/src/Core
parent1947da39ab18813ac5611f81c034eaca3ddc98a5 (diff)
Document enums in Constants.h (bug #248).
To get the links to work, I also had to document the Eigen namespace. Unfortunately, this means that the word Eigen is linked whenever it appears in the docs.
Diffstat (limited to 'Eigen/src/Core')
-rw-r--r--Eigen/src/Core/BandMatrix.h2
-rw-r--r--Eigen/src/Core/DenseCoeffsBase.h14
-rw-r--r--Eigen/src/Core/GenericPacketMath.h4
-rw-r--r--Eigen/src/Core/Map.h8
-rw-r--r--Eigen/src/Core/Matrix.h4
-rw-r--r--Eigen/src/Core/SelfAdjointView.h2
-rw-r--r--Eigen/src/Core/TriangularMatrix.h14
-rw-r--r--Eigen/src/Core/VectorwiseOp.h8
-rw-r--r--Eigen/src/Core/util/Constants.h167
9 files changed, 170 insertions, 53 deletions
diff --git a/Eigen/src/Core/BandMatrix.h b/Eigen/src/Core/BandMatrix.h
index a1f71d5f6..dda8efba3 100644
--- a/Eigen/src/Core/BandMatrix.h
+++ b/Eigen/src/Core/BandMatrix.h
@@ -180,7 +180,7 @@ class BandMatrixBase : public EigenBase<Derived>
* \param Cols Number of columns, or \b Dynamic
* \param Supers Number of super diagonal
* \param Subs Number of sub diagonal
- * \param _Options A combination of either \b RowMajor or \b ColMajor, and of \b SelfAdjoint
+ * \param _Options A combination of either \b #RowMajor or \b #ColMajor, and of \b #SelfAdjoint
* The former controls \ref TopicStorageOrders "storage order", and defaults to
* column-major. The latter controls whether the matrix represents a selfadjoint
* matrix in which case either Supers of Subs have to be null.
diff --git a/Eigen/src/Core/DenseCoeffsBase.h b/Eigen/src/Core/DenseCoeffsBase.h
index 7838a1cfb..e45238fb5 100644
--- a/Eigen/src/Core/DenseCoeffsBase.h
+++ b/Eigen/src/Core/DenseCoeffsBase.h
@@ -35,7 +35,7 @@ template<typename T> struct add_const_on_value_type_if_arithmetic
/** \brief Base class providing read-only coefficient access to matrices and arrays.
* \ingroup Core_Module
* \tparam Derived Type of the derived class
- * \tparam ReadOnlyAccessors Constant indicating read-only access
+ * \tparam #ReadOnlyAccessors Constant indicating read-only access
*
* This class defines the \c operator() \c const function and friends, which can be used to read specific
* entries of a matrix or array.
@@ -212,7 +212,7 @@ class DenseCoeffsBase<Derived,ReadOnlyAccessors> : public EigenBase<Derived>
* to ensure that a packet really starts there. This method is only available on expressions having the
* PacketAccessBit.
*
- * The \a LoadMode parameter may have the value \a Aligned or \a Unaligned. Its effect is to select
+ * The \a LoadMode parameter may have the value \a #Aligned or \a #Unaligned. Its effect is to select
* the appropriate vectorization instruction. Aligned access is faster, but is only possible for packets
* starting at an address which is a multiple of the packet size.
*/
@@ -239,7 +239,7 @@ class DenseCoeffsBase<Derived,ReadOnlyAccessors> : public EigenBase<Derived>
* to ensure that a packet really starts there. This method is only available on expressions having the
* PacketAccessBit and the LinearAccessBit.
*
- * The \a LoadMode parameter may have the value \a Aligned or \a Unaligned. Its effect is to select
+ * The \a LoadMode parameter may have the value \a #Aligned or \a #Unaligned. Its effect is to select
* the appropriate vectorization instruction. Aligned access is faster, but is only possible for packets
* starting at an address which is a multiple of the packet size.
*/
@@ -275,7 +275,7 @@ class DenseCoeffsBase<Derived,ReadOnlyAccessors> : public EigenBase<Derived>
/** \brief Base class providing read/write coefficient access to matrices and arrays.
* \ingroup Core_Module
* \tparam Derived Type of the derived class
- * \tparam WriteAccessors Constant indicating read/write access
+ * \tparam #WriteAccessors Constant indicating read/write access
*
* This class defines the non-const \c operator() function and friends, which can be used to write specific
* entries of a matrix or array. This class inherits DenseCoeffsBase<Derived, ReadOnlyAccessors> which
@@ -433,7 +433,7 @@ class DenseCoeffsBase<Derived, WriteAccessors> : public DenseCoeffsBase<Derived,
* to ensure that a packet really starts there. This method is only available on expressions having the
* PacketAccessBit.
*
- * The \a LoadMode parameter may have the value \a Aligned or \a Unaligned. Its effect is to select
+ * The \a LoadMode parameter may have the value \a #Aligned or \a #Unaligned. Its effect is to select
* the appropriate vectorization instruction. Aligned access is faster, but is only possible for packets
* starting at an address which is a multiple of the packet size.
*/
@@ -567,7 +567,7 @@ class DenseCoeffsBase<Derived, WriteAccessors> : public DenseCoeffsBase<Derived,
/** \brief Base class providing direct read-only coefficient access to matrices and arrays.
* \ingroup Core_Module
* \tparam Derived Type of the derived class
- * \tparam DirectAccessors Constant indicating direct access
+ * \tparam #DirectAccessors Constant indicating direct access
*
* This class defines functions to work with strides which can be used to access entries directly. This class
* inherits DenseCoeffsBase<Derived, ReadOnlyAccessors> which defines functions to access entries read-only using
@@ -637,7 +637,7 @@ class DenseCoeffsBase<Derived, DirectAccessors> : public DenseCoeffsBase<Derived
/** \brief Base class providing direct read/write coefficient access to matrices and arrays.
* \ingroup Core_Module
* \tparam Derived Type of the derived class
- * \tparam DirectAccessors Constant indicating direct access
+ * \tparam #DirectWriteAccessors Constant indicating direct access
*
* This class defines functions to work with strides which can be used to access entries directly. This class
* inherits DenseCoeffsBase<Derived, WriteAccessors> which defines functions to access entries read/write using
diff --git a/Eigen/src/Core/GenericPacketMath.h b/Eigen/src/Core/GenericPacketMath.h
index 4ba322a32..2cedab04a 100644
--- a/Eigen/src/Core/GenericPacketMath.h
+++ b/Eigen/src/Core/GenericPacketMath.h
@@ -286,7 +286,7 @@ pmadd(const Packet& a,
{ return padd(pmul(a, b),c); }
/** \internal \returns a packet version of \a *from.
- * \If LoadMode equals Aligned, \a from must be 16 bytes aligned */
+ * If LoadMode equals #Aligned, \a from must be 16 bytes aligned */
template<typename Packet, int LoadMode>
inline Packet ploadt(const typename unpacket_traits<Packet>::type* from)
{
@@ -297,7 +297,7 @@ inline Packet ploadt(const typename unpacket_traits<Packet>::type* from)
}
/** \internal copy the packet \a from to \a *to.
- * If StoreMode equals Aligned, \a to must be 16 bytes aligned */
+ * If StoreMode equals #Aligned, \a to must be 16 bytes aligned */
template<typename Scalar, typename Packet, int LoadMode>
inline void pstoret(Scalar* to, const Packet& from)
{
diff --git a/Eigen/src/Core/Map.h b/Eigen/src/Core/Map.h
index 692d0a179..81e3979f3 100644
--- a/Eigen/src/Core/Map.h
+++ b/Eigen/src/Core/Map.h
@@ -31,10 +31,10 @@
*
* \brief A matrix or vector expression mapping an existing array of data.
*
- * \param PlainObjectType the equivalent matrix type of the mapped data
- * \param MapOptions specifies whether the pointer is \c Aligned, or \c Unaligned.
- * The default is \c Unaligned.
- * \param StrideType optionnally specifies strides. By default, Map assumes the memory layout
+ * \tparam PlainObjectType the equivalent matrix type of the mapped data
+ * \tparam MapOptions specifies whether the pointer is \c #Aligned, or \c #Unaligned.
+ * The default is \c #Unaligned.
+ * \tparam StrideType optionnally specifies strides. By default, Map assumes the memory layout
* of an ordinary, contiguous array. This can be overridden by specifying strides.
* The type passed here must be a specialization of the Stride template, see examples below.
*
diff --git a/Eigen/src/Core/Matrix.h b/Eigen/src/Core/Matrix.h
index 8ae55da6e..44de22cb4 100644
--- a/Eigen/src/Core/Matrix.h
+++ b/Eigen/src/Core/Matrix.h
@@ -43,8 +43,8 @@
* \tparam _Cols Number of columns, or \b Dynamic
*
* The remaining template parameters are optional -- in most cases you don't have to worry about them.
- * \tparam _Options \anchor matrix_tparam_options A combination of either \b RowMajor or \b ColMajor, and of either
- * \b AutoAlign or \b DontAlign.
+ * \tparam _Options \anchor matrix_tparam_options A combination of either \b #RowMajor or \b #ColMajor, and of either
+ * \b #AutoAlign or \b #DontAlign.
* The former controls \ref TopicStorageOrders "storage order", and defaults to column-major. The latter controls alignment, which is required
* for vectorization. It defaults to aligning matrices except for fixed sizes that aren't a multiple of the packet size.
* \tparam _MaxRows Maximum number of rows. Defaults to \a _Rows (\ref maxrows "note").
diff --git a/Eigen/src/Core/SelfAdjointView.h b/Eigen/src/Core/SelfAdjointView.h
index 4bb90c373..4bb68755e 100644
--- a/Eigen/src/Core/SelfAdjointView.h
+++ b/Eigen/src/Core/SelfAdjointView.h
@@ -32,7 +32,7 @@
* \brief Expression of a selfadjoint matrix from a triangular part of a dense matrix
*
* \param MatrixType the type of the dense matrix storing the coefficients
- * \param TriangularPart can be either \c Lower or \c Upper
+ * \param TriangularPart can be either \c #Lower or \c #Upper
*
* This class is an expression of a sefladjoint matrix from a triangular part of a matrix
* with given dense storage of the coefficients. It is the return type of MatrixBase::selfadjointView()
diff --git a/Eigen/src/Core/TriangularMatrix.h b/Eigen/src/Core/TriangularMatrix.h
index f9fedcb0f..fee751acd 100644
--- a/Eigen/src/Core/TriangularMatrix.h
+++ b/Eigen/src/Core/TriangularMatrix.h
@@ -134,13 +134,13 @@ template<typename Derived> class TriangularBase : public EigenBase<Derived>
* \brief Base class for triangular part in a matrix
*
* \param MatrixType the type of the object in which we are taking the triangular part
- * \param Mode the kind of triangular matrix expression to construct. Can be Upper,
- * Lower, UpperSelfadjoint, or LowerSelfadjoint. This is in fact a bit field;
- * it must have either Upper or Lower, and additionnaly it may have either
- * UnitDiag or Selfadjoint.
+ * \param Mode the kind of triangular matrix expression to construct. Can be #Upper,
+ * #Lower, #UnitUpper, #UnitLower, #StrictlyUpper, or #StrictlyLower.
+ * This is in fact a bit field; it must have either #Upper or #Lower,
+ * and additionnaly it may have #UnitDiag or #ZeroDiag or neither.
*
* This class represents a triangular part of a matrix, not necessarily square. Strictly speaking, for rectangular
- * matrices one should speak ok "trapezoid" parts. This class is the return type
+ * matrices one should speak of "trapezoid" parts. This class is the return type
* of MatrixBase::triangularView() and most of the time this is the only way it is used.
*
* \sa MatrixBase::triangularView()
@@ -756,8 +756,8 @@ typename internal::eigen2_part_return_type<Derived, Mode>::type MatrixBase<Deriv
/**
* \returns an expression of a triangular view extracted from the current matrix
*
- * The parameter \a Mode can have the following values: \c Upper, \c StrictlyUpper, \c UnitUpper,
- * \c Lower, \c StrictlyLower, \c UnitLower.
+ * The parameter \a Mode can have the following values: \c #Upper, \c #StrictlyUpper, \c #UnitUpper,
+ * \c #Lower, \c #StrictlyLower, \c #UnitLower.
*
* Example: \include MatrixBase_extract.cpp
* Output: \verbinclude MatrixBase_extract.out
diff --git a/Eigen/src/Core/VectorwiseOp.h b/Eigen/src/Core/VectorwiseOp.h
index e328d94aa..20f688157 100644
--- a/Eigen/src/Core/VectorwiseOp.h
+++ b/Eigen/src/Core/VectorwiseOp.h
@@ -31,9 +31,9 @@
*
* \brief Generic expression of a partially reduxed matrix
*
- * \param MatrixType the type of the matrix we are applying the redux operation
- * \param MemberOp type of the member functor
- * \param Direction indicates the direction of the redux (Vertical or Horizontal)
+ * \tparam MatrixType the type of the matrix we are applying the redux operation
+ * \tparam MemberOp type of the member functor
+ * \tparam Direction indicates the direction of the redux (#Vertical or #Horizontal)
*
* This class represents an expression of a partial redux operator of a matrix.
* It is the return type of some VectorwiseOp functions,
@@ -164,7 +164,7 @@ struct member_redux {
* \brief Pseudo expression providing partial reduction operations
*
* \param ExpressionType the type of the object on which to do partial reductions
- * \param Direction indicates the direction of the redux (Vertical or Horizontal)
+ * \param Direction indicates the direction of the redux (#Vertical or #Horizontal)
*
* This class represents a pseudo expression with partial reduction features.
* It is the return type of DenseBase::colwise() and DenseBase::rowwise()
diff --git a/Eigen/src/Core/util/Constants.h b/Eigen/src/Core/util/Constants.h
index 2ffeb7948..37f18f2b4 100644
--- a/Eigen/src/Core/util/Constants.h
+++ b/Eigen/src/Core/util/Constants.h
@@ -161,23 +161,72 @@ const unsigned int HereditaryBits = RowMajorBit
| EvalBeforeNestingBit
| EvalBeforeAssigningBit;
-// Possible values for the Mode parameter of triangularView()
+/** \defgroup enums Enumerations
+ * \ingroup Core_Module
+ *
+ * Various enumerations used in %Eigen. Many of these are used as template parameters.
+ */
+
+/** \ingroup enums
+ * Enum containing possible values for the \p Mode parameter of
+ * MatrixBase::selfadjointView() and MatrixBase::triangularView(). */
enum {
- Lower=0x1, Upper=0x2, UnitDiag=0x4, ZeroDiag=0x8,
- UnitLower=UnitDiag|Lower, UnitUpper=UnitDiag|Upper,
- StrictlyLower=ZeroDiag|Lower, StrictlyUpper=ZeroDiag|Upper,
- SelfAdjoint=0x10};
+ /** View matrix as a lower triangular matrix. */
+ Lower=0x1,
+ /** View matrix as an upper triangular matrix. */
+ Upper=0x2,
+ /** %Matrix has ones on the diagonal; to be used in combination with #Lower or #Upper. */
+ UnitDiag=0x4,
+ /** %Matrix has zeros on the diagonal; to be used in combination with #Lower or #Upper. */
+ ZeroDiag=0x8,
+ /** View matrix as a lower triangular matrix with ones on the diagonal. */
+ UnitLower=UnitDiag|Lower,
+ /** View matrix as an upper triangular matrix with ones on the diagonal. */
+ UnitUpper=UnitDiag|Upper,
+ /** View matrix as a lower triangular matrix with zeros on the diagonal. */
+ StrictlyLower=ZeroDiag|Lower,
+ /** View matrix as an upper triangular matrix with zeros on the diagonal. */
+ StrictlyUpper=ZeroDiag|Upper,
+ /** Used in BandMatrix and SelfAdjointView to indicate that the matrix is self-adjoint. */
+ SelfAdjoint=0x10
+};
+
+/** \ingroup enums
+ * Enum for indicating whether an object is aligned or not. */
+enum {
+ /** Object is not correctly aligned for vectorization. */
+ Unaligned=0,
+ /** Object is aligned for vectorization. */
+ Aligned=1
+};
-enum { Unaligned=0, Aligned=1 };
enum { ConditionalJumpCost = 5 };
+/** \ingroup enums
+ * Enum used by DenseBase::corner() in Eigen2 compatibility mode. */
// FIXME after the corner() API change, this was not needed anymore, except by AlignedBox
// TODO: find out what to do with that. Adapt the AlignedBox API ?
enum CornerType { TopLeft, TopRight, BottomLeft, BottomRight };
-enum DirectionType { Vertical, Horizontal, BothDirections };
+/** \ingroup enums
+ * Enum containing possible values for the \p Direction parameter of
+ * Reverse, PartialReduxExpr and VectorwiseOp. */
+enum DirectionType {
+ /** For Reverse, all columns are reversed;
+ * for PartialReduxExpr and VectorwiseOp, act on columns. */
+ Vertical,
+ /** For Reverse, all rows are reversed;
+ * for PartialReduxExpr and VectorwiseOp, act on rows. */
+ Horizontal,
+ /** For Reverse, both rows and columns are reversed;
+ * not used for PartialReduxExpr and VectorwiseOp. */
+ BothDirections
+};
+
enum ProductEvaluationMode { NormalProduct, CacheFriendlyProduct };
+/** \internal \ingroup enums
+ * Enum to specify how to traverse the entries of a matrix. */
enum {
/** \internal Default traversal, no vectorization, no index-based access */
DefaultTraversal,
@@ -196,14 +245,25 @@ enum {
InvalidTraversal
};
+/** \internal \ingroup enums
+ * Enum to specify whether to unroll loops when traversing over the entries of a matrix. */
enum {
+ /** \internal Do not unroll loops. */
NoUnrolling,
+ /** \internal Unroll only the inner loop, but not the outer loop. */
InnerUnrolling,
+ /** \internal Unroll both the inner and the outer loop. If there is only one loop,
+ * because linear traversal is used, then unroll that loop. */
CompleteUnrolling
};
+/** \ingroup enums
+ * Enum containing possible values for the \p _Options template parameter of
+ * Matrix, Array and BandMatrix. */
enum {
+ /** Storage order is column major (see \ref TopicStorageOrders). */
ColMajor = 0,
+ /** Storage order is row major (see \ref TopicStorageOrders). */
RowMajor = 0x1, // it is only a coincidence that this is equal to RowMajorBit -- don't rely on that
/** \internal Align the matrix itself if it is vectorizable fixed-size */
AutoAlign = 0,
@@ -211,11 +271,13 @@ enum {
DontAlign = 0x2
};
-/** \brief Enum for specifying whether to apply or solve on the left or right.
- */
+/** \ingroup enums
+ * Enum for specifying whether to apply or solve on the left or right. */
enum {
- OnTheLeft = 1, /**< \brief Apply transformation on the left. */
- OnTheRight = 2 /**< \brief Apply transformation on the right. */
+ /** Apply transformation on the left. */
+ OnTheLeft = 1,
+ /** Apply transformation on the right. */
+ OnTheRight = 2
};
/* the following could as well be written:
@@ -239,53 +301,104 @@ namespace {
EIGEN_UNUSED Default_t Default;
}
+/** \internal \ingroup enums
+ * Used in AmbiVector. */
enum {
IsDense = 0,
IsSparse
};
+/** \ingroup enums
+ * Used as template parameter in DenseCoeffBase and MapBase to indicate
+ * which accessors should be provided. */
enum AccessorLevels {
- ReadOnlyAccessors, WriteAccessors, DirectAccessors, DirectWriteAccessors
+ /** Read-only access via a member function. */
+ ReadOnlyAccessors,
+ /** Read/write access via member functions. */
+ WriteAccessors,
+ /** Direct read-only access to the coefficients. */
+ DirectAccessors,
+ /** Direct read/write access to the coefficients. */
+ DirectWriteAccessors
};
+/** \ingroup enums
+ * Enum with options to give to various decompositions. */
enum DecompositionOptions {
- Pivoting = 0x01, // LDLT,
- NoPivoting = 0x02, // LDLT,
- ComputeFullU = 0x04, // SVD,
- ComputeThinU = 0x08, // SVD,
- ComputeFullV = 0x10, // SVD,
- ComputeThinV = 0x20, // SVD,
- EigenvaluesOnly = 0x40, // all eigen solvers
- ComputeEigenvectors = 0x80, // all eigen solvers
+ /** \internal Not used (meant for LDLT?). */
+ Pivoting = 0x01,
+ /** \internal Not used (meant for LDLT?). */
+ NoPivoting = 0x02,
+ /** Used in JacobiSVD to indicate that the square matrix U is to be computed. */
+ ComputeFullU = 0x04,
+ /** Used in JacobiSVD to indicate that the thin matrix U is to be computed. */
+ ComputeThinU = 0x08,
+ /** Used in JacobiSVD to indicate that the square matrix V is to be computed. */
+ ComputeFullV = 0x10,
+ /** Used in JacobiSVD to indicate that the thin matrix V is to be computed. */
+ ComputeThinV = 0x20,
+ /** Used in SelfAdjointEigenSolver and GeneralizedSelfAdjointEigenSolver to specify
+ * that only the eigenvalues are to be computed and not the eigenvectors. */
+ EigenvaluesOnly = 0x40,
+ /** Used in SelfAdjointEigenSolver and GeneralizedSelfAdjointEigenSolver to specify
+ * that both the eigenvalues and the eigenvectors are to be computed. */
+ ComputeEigenvectors = 0x80,
+ /** \internal */
EigVecMask = EigenvaluesOnly | ComputeEigenvectors,
+ /** Used in GeneralizedSelfAdjointEigenSolver to indicate that it should
+ * solve the generalized eigenproblem \f$ Ax = \lambda B x \f$. */
Ax_lBx = 0x100,
+ /** Used in GeneralizedSelfAdjointEigenSolver to indicate that it should
+ * solve the generalized eigenproblem \f$ ABx = \lambda x \f$. */
ABx_lx = 0x200,
+ /** Used in GeneralizedSelfAdjointEigenSolver to indicate that it should
+ * solve the generalized eigenproblem \f$ BAx = \lambda x \f$. */
BAx_lx = 0x400,
+ /** \internal */
GenEigMask = Ax_lBx | ABx_lx | BAx_lx
};
+/** \ingroup enums
+ * Possible values for the \p QRPreconditioner template parameter of JacobiSVD. */
enum QRPreconditioners {
+ /** Do not specify what is to be done if the SVD of a non-square matrix is asked for. */
NoQRPreconditioner,
+ /** Use a QR decomposition without pivoting as the first step. */
HouseholderQRPreconditioner,
+ /** Use a QR decomposition with column pivoting as the first step. */
ColPivHouseholderQRPreconditioner,
+ /** Use a QR decomposition with full pivoting as the first step. */
FullPivHouseholderQRPreconditioner
};
-/** \brief Enum for reporting the status of a computation.
- */
+/** \ingroups enums
+ * Enum for reporting the status of a computation. */
enum ComputationInfo {
- Success = 0, /**< \brief Computation was successful. */
- NumericalIssue = 1, /**< \brief The provided data did not satisfy the prerequisites. */
- NoConvergence = 2 /**< \brief Iterative procedure did not converge. */
+ /** Computation was successful. */
+ Success = 0,
+ /** The provided data did not satisfy the prerequisites. */
+ NumericalIssue = 1,
+ /** Iterative procedure did not converge. */
+ NoConvergence = 2
};
+/** \ingroup enums
+ * Enum used to specify how a particular transformation is stored in a matrix.
+ * \sa Transform, Hyperplane::transform(). */
enum TransformTraits {
+ /** Transformation is an isometry. */
Isometry = 0x1,
+ /** Transformation is an affine transformation stored as a (Dim+1)^2 matrix whose last row is
+ * assumed to be [0 ... 0 1]. */
Affine = 0x2,
+ /** Transformation is an affine transformation stored as a (Dim) x (Dim+1) matrix. */
AffineCompact = 0x10 | Affine,
+ /** Transformation is a general projective transformation stored as a (Dim+1)^2 matrix. */
Projective = 0x20
};
+/** \internal \ingroup enums
+ * Enum used to choose between implementation depending on the computer architecture. */
namespace Architecture
{
enum Type {
@@ -302,8 +415,12 @@ namespace Architecture
};
}
+/** \internal \ingroup enums
+ * Enum used as template parameter in GeneralProduct. */
enum { CoeffBasedProductMode, LazyCoeffBasedProductMode, OuterProduct, InnerProduct, GemvProduct, GemmProduct };
+/** \internal \ingroup enums
+ * Enum used in experimental parallel implementation. */
enum Action {GetAction, SetAction};
/** The type used to identify a dense storage. */