* move memory related stuff to util/Memory.h

* clean ugly doxygen inheritence of expressions
* keep improving the documentation... slowly !
This commit is contained in:
Gael Guennebaud 2008-08-26 19:12:23 +00:00
parent 3e526dcdbd
commit 00a8d314c5
35 changed files with 445 additions and 188 deletions

View File

@ -6,16 +6,27 @@
#include <cstring>
#include <string>
#ifdef EIGEN_VECTORIZE
// it seems we cannot assume posix_memalign is defined in the stdlib header
extern "C" int posix_memalign (void **, size_t, size_t) throw ();
#endif
namespace Eigen {
/** \defgroup Core_Module Core module
* This is the main module of Eigen providing dense matrix and vector support with
* features equivalent to a BLAS library and much more...
*
* \code
* #include <Eigen/Core>
* \endcode
*
* To speedup compilation, in header files you might only include the following file:
* \code
* #include <Eigen/CoreDeclaration>
* \endcode
* which gives access to matrix and vector types as well as coefficient accessor operators.
*/
#include "src/Core/util/Memory.h"
#include "src/Core/NumTraits.h"
#include "src/Core/MathFunctions.h"
#include "src/Core/DummyPacketMath.h"
#include "src/Core/GenericPacketMath.h"
#if defined EIGEN_VECTORIZE_SSE
#include "src/Core/arch/SSE/PacketMath.h"

View File

@ -15,8 +15,8 @@ namespace Eigen {
* - 2D and 3D rotations
* - quaternions
* - \ref MatrixBase::cross() "cross product"
* - \ref MatrixBase::someOrthognal() "orthognal vector generation"
*
* - \ref MatrixBase::unitOrthogonal() "orthognal vector generation"
*
* \code
* #include <Eigen/Geometry>
* \endcode

View File

@ -27,7 +27,7 @@
#define EIGEN_PARTIAL_REDUX_H
/** \array_module \ingroup Array
*
*
* \class PartialReduxExpr
*
* \brief Generic expression of a partially reduxed matrix
@ -69,7 +69,11 @@ struct ei_traits<PartialReduxExpr<MatrixType, MemberOp, Direction> >
template< typename MatrixType, typename MemberOp, int Direction>
class PartialReduxExpr : ei_no_assignment_operator,
#ifndef EIGEN_PARSED_BY_DOXYGEN
public MatrixBase<PartialReduxExpr<MatrixType, MemberOp, Direction> >
#else
public MapBase
#endif
{
public:
@ -186,7 +190,7 @@ template<typename ExpressionType, int Direction> class PartialRedux
*
* Example: \include PartialRedux_minCoeff.cpp
* Output: \verbinclude PartialRedux_minCoeff.out
*
*
* \sa MatrixBase::minCoeff() */
const typename ReturnType<ei_member_minCoeff>::Type minCoeff() const
{ return _expression(); }
@ -196,7 +200,7 @@ template<typename ExpressionType, int Direction> class PartialRedux
*
* Example: \include PartialRedux_maxCoeff.cpp
* Output: \verbinclude PartialRedux_maxCoeff.out
*
*
* \sa MatrixBase::maxCoeff() */
const typename ReturnType<ei_member_maxCoeff>::Type maxCoeff() const
{ return _expression(); }
@ -206,7 +210,7 @@ template<typename ExpressionType, int Direction> class PartialRedux
*
* Example: \include PartialRedux_norm2.cpp
* Output: \verbinclude PartialRedux_norm2.out
*
*
* \sa MatrixBase::norm2() */
const typename ReturnType<ei_member_norm2>::Type norm2() const
{ return _expression(); }

View File

@ -88,7 +88,11 @@ struct ei_traits<Block<MatrixType, BlockRows, BlockCols, _PacketAccess, _DirectA
};
template<typename MatrixType, int BlockRows, int BlockCols, int PacketAccess, int _DirectAccessStatus> class Block
#ifndef EIGEN_PARSED_BY_DOXYGEN
: public MatrixBase<Block<MatrixType, BlockRows, BlockCols, PacketAccess, _DirectAccessStatus> >
#else
: public MatrixBase
#endif
{
public:
@ -211,7 +215,11 @@ template<typename MatrixType, int BlockRows, int BlockCols, int PacketAccess, in
/** \internal */
template<typename MatrixType, int BlockRows, int BlockCols, int PacketAccess>
class Block<MatrixType,BlockRows,BlockCols,PacketAccess,HasDirectAccess>
#ifndef EIGEN_PARSED_BY_DOXYGEN
: public MapBase<Block<MatrixType, BlockRows, BlockCols,PacketAccess,HasDirectAccess> >
#else
: public MapBase
#endif
{
public:

View File

@ -76,7 +76,11 @@ struct ei_traits<CwiseBinaryOp<BinaryOp, Lhs, Rhs> >
template<typename BinaryOp, typename Lhs, typename Rhs>
class CwiseBinaryOp : ei_no_assignment_operator,
#ifndef EIGEN_PARSED_BY_DOXYGEN
public MatrixBase<CwiseBinaryOp<BinaryOp, Lhs, Rhs> >
#else
public MatrixBase
#endif
{
public:
@ -183,7 +187,7 @@ MatrixBase<Derived>::operator+=(const MatrixBase<OtherDerived>& other)
*
* Example: \include Cwise_product.cpp
* Output: \verbinclude Cwise_product.out
*
*
* \sa class CwiseBinaryOp, operator/(), square()
*/
template<typename ExpressionType>

View File

@ -60,7 +60,11 @@ struct ei_traits<CwiseNullaryOp<NullaryOp, MatrixType> >
template<typename NullaryOp, typename MatrixType>
class CwiseNullaryOp : ei_no_assignment_operator,
#ifndef EIGEN_PARSED_BY_DOXYGEN
public MatrixBase<CwiseNullaryOp<NullaryOp, MatrixType> >
#else
public MatrixBase
#endif
{
public:

View File

@ -63,7 +63,11 @@ struct ei_traits<CwiseUnaryOp<UnaryOp, MatrixType> >
template<typename UnaryOp, typename MatrixType>
class CwiseUnaryOp : ei_no_assignment_operator,
#ifndef EIGEN_PARSED_BY_DOXYGEN
public MatrixBase<CwiseUnaryOp<UnaryOp, MatrixType> >
#else
public MatrixBase
#endif
{
public:

View File

@ -60,7 +60,11 @@ struct ei_traits<DiagonalCoeffs<MatrixType> >
};
template<typename MatrixType> class DiagonalCoeffs
: public MatrixBase<DiagonalCoeffs<MatrixType> >
#ifndef EIGEN_PARSED_BY_DOXYGEN
: public MatrixBase<DiagonalCoeffs<MatrixType> >
#else
: public MatrixBase
#endif
{
public:

View File

@ -56,7 +56,11 @@ struct ei_traits<DiagonalMatrix<CoeffsVectorType> >
template<typename CoeffsVectorType>
class DiagonalMatrix : ei_no_assignment_operator,
public MatrixBase<DiagonalMatrix<CoeffsVectorType> >
#ifndef EIGEN_PARSED_BY_DOXYGEN
public MatrixBase<DiagonalMatrix<CoeffsVectorType> >
#else
public MatrixBase
#endif
{
public:

View File

@ -61,7 +61,11 @@ struct ei_traits<Product<LhsNested, RhsNested, DiagonalProduct> >
};
template<typename LhsNested, typename RhsNested> class Product<LhsNested, RhsNested, DiagonalProduct> : ei_no_assignment_operator,
#ifndef EIGEN_PARSED_BY_DOXYGEN
public MatrixBase<Product<LhsNested, RhsNested, DiagonalProduct> >
#else
public MatrixBase
#endif
{
typedef typename ei_traits<Product>::_LhsNested _LhsNested;
typedef typename ei_traits<Product>::_RhsNested _RhsNested;

View File

@ -55,7 +55,11 @@ struct ei_traits<Flagged<ExpressionType, Added, Removed> >
};
template<typename ExpressionType, unsigned int Added, unsigned int Removed> class Flagged
#ifndef EIGEN_PARSED_BY_DOXYGEN
: public MatrixBase<Flagged<ExpressionType, Added, Removed> >
#else
: public MatrixBase
#endif
{
public:

View File

@ -2,6 +2,7 @@
// for linear algebra. Eigen itself is part of the KDE project.
//
// Copyright (C) 2008 Gael Guennebaud <g.gael@free.fr>
// Copyright (C) 2006-2008 Benoit Jacob <jacob@math.jussieu.fr>
//
// Eigen is free software; you can redistribute it and/or
// modify it under the terms of the GNU Lesser General Public
@ -22,88 +23,95 @@
// License and a copy of the GNU General Public License along with
// Eigen. If not, see <http://www.gnu.org/licenses/>.
#ifndef EIGEN_DUMMY_PACKET_MATH_H
#define EIGEN_DUMMY_PACKET_MATH_H
#ifndef EIGEN_GENERIC_PACKET_MATH_H
#define EIGEN_GENERIC_PACKET_MATH_H
/** \internal
* \file GenericPacketMath.h
*
* Default implementation for types not supported by the vectorization.
* In practice these functions are provided to make easier the writting
* of generic vectorized code.
*/
// Default implementation for types not supported by the vectorization.
// In practice these functions are provided to make easier the writting
// of generic vectorized code. However, at runtime, they should never be
// called, TODO so sould we raise an assertion or not ?
/** \internal \returns a + b (coeff-wise) */
template <typename Packet> inline Packet
template<typename Packet> inline Packet
ei_padd(const Packet& a,
const Packet& b) { return a+b; }
/** \internal \returns a - b (coeff-wise) */
template <typename Packet> inline Packet
template<typename Packet> inline Packet
ei_psub(const Packet& a,
const Packet& b) { return a-b; }
/** \internal \returns a * b (coeff-wise) */
template <typename Packet> inline Packet
template<typename Packet> inline Packet
ei_pmul(const Packet& a,
const Packet& b) { return a*b; }
/** \internal \returns a / b (coeff-wise) */
template <typename Packet> inline Packet
template<typename Packet> inline Packet
ei_pdiv(const Packet& a,
const Packet& b) { return a/b; }
/** \internal \returns the min of \a a and \a b (coeff-wise) */
template <typename Packet> inline Packet
template<typename Packet> inline Packet
ei_pmin(const Packet& a,
const Packet& b) { return std::min(a, b); }
/** \internal \returns the max of \a a and \a b (coeff-wise) */
template <typename Packet> inline Packet
template<typename Packet> inline Packet
ei_pmax(const Packet& a,
const Packet& b) { return std::max(a, b); }
/** \internal \returns a packet version of \a *from, from must be 16 bytes aligned */
template <typename Scalar> inline typename ei_packet_traits<Scalar>::type
template<typename Scalar> inline typename ei_packet_traits<Scalar>::type
ei_pload(const Scalar* from) { return *from; }
/** \internal \returns a packet version of \a *from, (un-aligned load) */
template <typename Scalar> inline typename ei_packet_traits<Scalar>::type
template<typename Scalar> inline typename ei_packet_traits<Scalar>::type
ei_ploadu(const Scalar* from) { return *from; }
/** \internal \returns a packet with constant coefficients \a a, e.g.: (a,a,a,a) */
template <typename Scalar> inline typename ei_packet_traits<Scalar>::type
template<typename Scalar> inline typename ei_packet_traits<Scalar>::type
ei_pset1(const Scalar& a) { return a; }
/** \internal copy the packet \a from to \a *to, \a to must be 16 bytes aligned */
template <typename Scalar, typename Packet> inline void ei_pstore(Scalar* to, const Packet& from)
template<typename Scalar, typename Packet> inline void ei_pstore(Scalar* to, const Packet& from)
{ (*to) = from; }
/** \internal copy the packet \a from to \a *to, (un-aligned store) */
template <typename Scalar, typename Packet> inline void ei_pstoreu(Scalar* to, const Packet& from)
template<typename Scalar, typename Packet> inline void ei_pstoreu(Scalar* to, const Packet& from)
{ (*to) = from; }
/** \internal \returns the first element of a packet */
template <typename Packet> inline typename ei_unpacket_traits<Packet>::type ei_pfirst(const Packet& a)
template<typename Packet> inline typename ei_unpacket_traits<Packet>::type ei_pfirst(const Packet& a)
{ return a; }
/** \internal \returns a packet where the element i contains the sum of the packet of \a vec[i] */
template <typename Packet> inline Packet
template<typename Packet> inline Packet
ei_preduxp(const Packet* vecs) { return vecs[0]; }
/** \internal \returns the sum of the elements of \a a*/
template <typename Packet> inline typename ei_unpacket_traits<Packet>::type ei_predux(const Packet& a)
template<typename Packet> inline typename ei_unpacket_traits<Packet>::type ei_predux(const Packet& a)
{ return a; }
////////////
/***************************************************************************
* The following functions might not have to be overwritten for vectorized types
***************************************************************************/
/** \internal \returns a * b + c (coeff-wise) */
template <typename Packet> inline Packet
template<typename Packet> inline Packet
ei_pmadd(const Packet& a,
const Packet& b,
const Packet& c)
{ return ei_padd(ei_pmul(a, b),c); }
/** \internal \returns a packet version of \a *from. If LoadMode equals Aligned, \a from must be 16 bytes aligned */
template <typename Scalar, int LoadMode> inline typename ei_packet_traits<Scalar>::type ei_ploadt(const Scalar* from)
/** \internal \returns a packet version of \a *from.
* \If LoadMode equals Aligned, \a from must be 16 bytes aligned */
template<typename Scalar, int LoadMode>
inline typename ei_packet_traits<Scalar>::type ei_ploadt(const Scalar* from)
{
if(LoadMode == Aligned)
return ei_pload(from);
@ -111,8 +119,10 @@ template <typename Scalar, int LoadMode> inline typename ei_packet_traits<Scalar
return ei_ploadu(from);
}
/** \internal copy the packet \a from to \a *to. If StoreMode equals Aligned, \a to must be 16 bytes aligned */
template <typename Scalar, typename Packet, int LoadMode> inline void ei_pstoret(Scalar* to, const Packet& from)
/** \internal copy the packet \a from to \a *to.
* If StoreMode equals Aligned, \a to must be 16 bytes aligned */
template<typename Scalar, typename Packet, int LoadMode>
inline void ei_pstoret(Scalar* to, const Packet& from)
{
if(LoadMode == Aligned)
ei_pstore(to, from);
@ -120,21 +130,7 @@ template <typename Scalar, typename Packet, int LoadMode> inline void ei_pstoret
ei_pstoreu(to, from);
}
/** \internal \returns the number of elements which have to be skipped such that data are aligned */
template<typename Scalar>
inline static int ei_alignmentOffset(const Scalar* ptr, int maxOffset)
{
typedef typename ei_packet_traits<Scalar>::type Packet;
const int PacketSize = ei_packet_traits<Scalar>::size;
const int PacketAlignedMask = PacketSize-1;
const bool Vectorized = PacketSize>1;
return Vectorized
? std::min<int>( (PacketSize - ((size_t(ptr)/sizeof(Scalar)) & PacketAlignedMask))
& PacketAlignedMask, maxOffset)
: 0;
}
/** \internal specialization of ei_palign() */
/** \internal default implementation of ei_palign() allowing partial specialization */
template<int Offset,typename PacketType>
struct ei_palign_impl
{
@ -150,5 +146,5 @@ inline void ei_palign(PacketType& first, const PacketType& second)
ei_palign_impl<Offset,PacketType>::run(first,second);
}
#endif // EIGEN_DUMMY_PACKET_MATH_H
#endif // EIGEN_GENERIC_PACKET_MATH_H

View File

@ -57,7 +57,11 @@ struct ei_traits<Map<MatrixType, _PacketAccess> > : public ei_traits<MatrixType>
};
template<typename MatrixType, int PacketAccess> class Map
#ifndef EIGEN_PARSED_BY_DOXYGEN
: public MapBase<Map<MatrixType, PacketAccess> >
#else
: public MapBase
#endif
{
public:

View File

@ -96,7 +96,12 @@ struct ei_traits<Matrix<_Scalar, _Rows, _Cols, _StorageOrder, _MaxRows, _MaxCols
};
template<typename _Scalar, int _Rows, int _Cols, int _StorageOrder, int _MaxRows, int _MaxCols>
class Matrix : public MatrixBase<Matrix<_Scalar, _Rows, _Cols, _StorageOrder, _MaxRows, _MaxCols> >
class Matrix
#ifndef EIGEN_PARSED_BY_DOXYGEN
: public MatrixBase<Matrix<_Scalar, _Rows, _Cols, _StorageOrder, _MaxRows, _MaxCols> >
#else
: public MatrixBase
#endif
{
public:
EIGEN_GENERIC_PUBLIC_INTERFACE(Matrix)
@ -369,6 +374,8 @@ class Matrix : public MatrixBase<Matrix<_Scalar, _Rows, _Cols, _StorageOrder, _M
};
/** \defgroup matrixtypedefs Global matrix typedefs
*
* \ingroup Core_Module
*
* Eigen defines several typedef shortcuts for most common matrix and vector types.
*

View File

@ -39,44 +39,6 @@
*/
template<typename T, int Size, int _Rows, int _Cols> class ei_matrix_storage;
template <typename T, int Size, bool Align> struct ei_aligned_array
{
EIGEN_ALIGN_128 T array[Size];
};
template <typename T, int Size> struct ei_aligned_array<T,Size,false>
{
T array[Size];
};
template<typename T>
inline T* ei_aligned_malloc(size_t size)
{
#ifdef EIGEN_VECTORIZE
if (ei_packet_traits<T>::size>1)
{
void* ptr;
if (posix_memalign(&ptr, 16, size*sizeof(T))==0)
return static_cast<T*>(ptr);
else
return 0;
}
else
#endif
return new T[size];
}
template<typename T>
inline void ei_aligned_free(T* ptr)
{
#ifdef EIGEN_VECTORIZE
if (ei_packet_traits<T>::size>1)
free(ptr);
else
#endif
delete[] ptr;
}
// purely fixed-size matrix
template<typename T, int Size, int _Rows, int _Cols> class ei_matrix_storage
{

View File

@ -58,7 +58,11 @@ struct ei_traits<Minor<MatrixType> >
};
template<typename MatrixType> class Minor
#ifndef EIGEN_PARSED_BY_DOXYGEN
: public MatrixBase<Minor<MatrixType> >
#else
: public MatrixBase
#endif
{
public:

View File

@ -52,7 +52,11 @@ struct ei_traits<NestByValue<ExpressionType> >
};
template<typename ExpressionType> class NestByValue
#ifndef EIGEN_PARSED_BY_DOXYGEN
: public MatrixBase<NestByValue<ExpressionType> >
#else
: public MatrixBase
#endif
{
public:

View File

@ -59,7 +59,11 @@ struct ei_traits<Part<MatrixType, Mode> >
};
template<typename MatrixType, unsigned int Mode> class Part
#ifndef EIGEN_PARSED_BY_DOXYGEN
: public MatrixBase<Part<MatrixType, Mode> >
#else
: public MatrixBase
#endif
{
public:

View File

@ -168,7 +168,11 @@ struct ei_traits<Product<LhsNested, RhsNested, ProductMode> >
};
template<typename LhsNested, typename RhsNested, int ProductMode> class Product : ei_no_assignment_operator,
#ifndef EIGEN_PARSED_BY_DOXYGEN
public MatrixBase<Product<LhsNested, RhsNested, ProductMode> >
#else
public MatrixBase
#endif
{
public:

View File

@ -46,7 +46,11 @@ struct ei_traits<SwapWrapper<ExpressionType> >
};
template<typename ExpressionType> class SwapWrapper
#ifndef EIGEN_PARSED_BY_DOXYGEN
: public MatrixBase<SwapWrapper<ExpressionType> >
#else
: public MatrixBase
#endif
{
public:

View File

@ -57,7 +57,11 @@ struct ei_traits<Transpose<MatrixType> >
};
template<typename MatrixType> class Transpose
#ifndef EIGEN_PARSED_BY_DOXYGEN
: public MatrixBase<Transpose<MatrixType> >
#else
: public MatrixBase
#endif
{
public:

View File

@ -29,6 +29,7 @@
const int Dynamic = 10000;
/** \defgroup flags
* \ingroup Core_Module
*
* These are the possible bits which can be OR'ed to constitute the flags of a matrix or
* expression.

View File

@ -153,16 +153,4 @@ _EIGEN_GENERIC_PUBLIC_INTERFACE(Derived, Eigen::MatrixBase<Derived>)
#define EIGEN_ENUM_MIN(a,b) (((int)a <= (int)b) ? (int)a : (int)b)
#define EIGEN_ENUM_MAX(a,b) (((int)a >= (int)b) ? (int)a : (int)b)
/* ei_alloc_stack(TYPE,SIZE) allocates sizeof(TYPE)*SIZE bytes on the stack if sizeof(TYPE)*SIZE is smaller
* than EIGEN_STACK_ALLOCATION_LIMIT. Otherwise the memory is allocated using the operator new.
* Data allocated with ei_alloc_stack must be freed calling ei_free_stack(PTR,TYPE,SIZE)
*/
#ifdef __linux__
# define ei_alloc_stack(TYPE,SIZE) ((sizeof(TYPE)*(SIZE)>16000000) ? new TYPE[SIZE] : (TYPE*)alloca(sizeof(TYPE)*(SIZE)))
# define ei_free_stack(PTR,TYPE,SIZE) if (sizeof(TYPE)*SIZE>16000000) delete[] PTR
#else
# define ei_alloc_stack(TYPE,SIZE) new TYPE[SIZE]
# define ei_free_stack(PTR,TYPE,SIZE) delete[] PTR
#endif
#endif // EIGEN_MACROS_H

View File

@ -0,0 +1,109 @@
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra. Eigen itself is part of the KDE project.
//
// Copyright (C) 2008 Gael Guennebaud <g.gael@free.fr>
// Copyright (C) 2006-2008 Benoit Jacob <jacob@math.jussieu.fr>
//
// Eigen is free software; you can redistribute it and/or
// modify it under the terms of the GNU Lesser General Public
// License as published by the Free Software Foundation; either
// version 3 of the License, or (at your option) any later version.
//
// Alternatively, you can redistribute it and/or
// modify it under the terms of the GNU General Public License as
// published by the Free Software Foundation; either version 2 of
// the License, or (at your option) any later version.
//
// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
// GNU General Public License for more details.
//
// You should have received a copy of the GNU Lesser General Public
// License and a copy of the GNU General Public License along with
// Eigen. If not, see <http://www.gnu.org/licenses/>.
#ifndef EIGEN_MEMORY_H
#define EIGEN_MEMORY_H
#ifdef EIGEN_VECTORIZE
// it seems we cannot assume posix_memalign is defined in the stdlib header
extern "C" int posix_memalign (void **, size_t, size_t) throw ();
#endif
/** \internal
* Static array automatically aligned if the total byte size is a multiple of 16
*/
template <typename T, int Size, bool Align> struct ei_aligned_array
{
EIGEN_ALIGN_128 T array[Size];
};
template <typename T, int Size> struct ei_aligned_array<T,Size,false>
{
T array[Size];
};
/** \internal allocates \a size * sizeof(\a T) bytes with a 16 bytes based alignement */
template<typename T>
inline T* ei_aligned_malloc(size_t size)
{
#ifdef EIGEN_VECTORIZE
if (ei_packet_traits<T>::size>1)
{
void* ptr;
if (posix_memalign(&ptr, 16, size*sizeof(T))==0)
return static_cast<T*>(ptr);
else
return 0;
}
else
#endif
return new T[size];
}
/** \internal free memory allocated with ei_aligned_malloc */
template<typename T>
inline void ei_aligned_free(T* ptr)
{
#ifdef EIGEN_VECTORIZE
if (ei_packet_traits<T>::size>1)
free(ptr);
else
#endif
delete[] ptr;
}
/** \internal \returns the number of elements which have to be skipped such that data are 16 bytes aligned */
template<typename Scalar>
inline static int ei_alignmentOffset(const Scalar* ptr, int maxOffset)
{
typedef typename ei_packet_traits<Scalar>::type Packet;
const int PacketSize = ei_packet_traits<Scalar>::size;
const int PacketAlignedMask = PacketSize-1;
const bool Vectorized = PacketSize>1;
return Vectorized
? std::min<int>( (PacketSize - ((size_t(ptr)/sizeof(Scalar)) & PacketAlignedMask))
& PacketAlignedMask, maxOffset)
: 0;
}
/** \internal
* ei_alloc_stack(TYPE,SIZE) allocates sizeof(TYPE)*SIZE bytes on the stack if sizeof(TYPE)*SIZE is
* smaller than EIGEN_STACK_ALLOCATION_LIMIT. Otherwise the memory is allocated using the operator new.
* Data allocated with ei_alloc_stack \b must be freed calling ei_free_stack(PTR,TYPE,SIZE).
* \code
* float * data = ei_alloc_stack(float,array.size());
* // ...
* ei_free_stack(data,float,array.size());
* \endcode
*/
#ifdef __linux__
# define ei_alloc_stack(TYPE,SIZE) ((sizeof(TYPE)*(SIZE)>16000000) ? new TYPE[SIZE] : (TYPE*)alloca(sizeof(TYPE)*(SIZE)))
# define ei_free_stack(PTR,TYPE,SIZE) if (sizeof(TYPE)*SIZE>16000000) delete[] PTR
#else
# define ei_alloc_stack(TYPE,SIZE) new TYPE[SIZE]
# define ei_free_stack(PTR,TYPE,SIZE) delete[] PTR
#endif
#endif // EIGEN_MEMORY_H

View File

@ -76,17 +76,17 @@
* Let's now describe precisely the parameters:
* @param numPoints the number of points
* @param points the array of pointers to the points on which to perform the linear regression
* @param retCoefficients pointer to the vector in which to store the result.
This vector must be of the same type and size as the
data points. The meaning of its coords is as follows.
For brevity, let \f$n=Size\f$,
\f$r_i=retCoefficients[i]\f$,
and \f$f=funcOfOthers\f$. Denote by
\f$x_0,\ldots,x_{n-1}\f$
the n coordinates in the n-dimensional space.
Then the result equation is:
\f[ x_f = r_0 x_0 + \cdots + r_{f-1}x_{f-1}
+ r_{f+1}x_{f+1} + \cdots + r_{n-1}x_{n-1} + r_n. \f]
* @param result pointer to the vector in which to store the result.
This vector must be of the same type and size as the
data points. The meaning of its coords is as follows.
For brevity, let \f$n=Size\f$,
\f$r_i=retCoefficients[i]\f$,
and \f$f=funcOfOthers\f$. Denote by
\f$x_0,\ldots,x_{n-1}\f$
the n coordinates in the n-dimensional space.
Then the result equation is:
\f[ x_f = r_0 x_0 + \cdots + r_{f-1}x_{f-1}
+ r_{f+1}x_{f+1} + \cdots + r_{n-1}x_{n-1} + r_n. \f]
* @param funcOfOthers Determines which coord to express as a function of the
others. Coords are numbered starting from 0, so that a
value of 0 means \f$x\f$, 1 means \f$y\f$,
@ -183,7 +183,7 @@ void fitHyperplane(int numPoints,
VectorType diff = (*(points[i]) - mean).conjugate();
covMat += diff * diff.adjoint();
}
// now we just have to pick the eigen vector with smallest eigen value
SelfAdjointEigenSolver<CovMatrixType> eig(covMat);
result->start(size) = eig.eigenvectors().col(0);

View File

@ -42,7 +42,12 @@ struct ei_traits<HashMatrix<_Scalar, _Flags> >
// TODO reimplement this class using custom linked lists
template<typename _Scalar, int _Flags>
class HashMatrix : public SparseMatrixBase<HashMatrix<_Scalar, _Flags> >
class HashMatrix
#ifndef EIGEN_PARSED_BY_DOXYGEN
: public SparseMatrixBase<HashMatrix<_Scalar, _Flags> >
#else
: public SparseMatrixBase
#endif
{
public:
EIGEN_GENERIC_PUBLIC_INTERFACE(HashMatrix)

View File

@ -52,7 +52,12 @@ struct LinkedVectorChunk
};
template<typename _Scalar, int _Flags>
class LinkedVectorMatrix : public SparseMatrixBase<LinkedVectorMatrix<_Scalar,_Flags> >
class LinkedVectorMatrix
#ifndef EIGEN_PARSED_BY_DOXYGEN
: public SparseMatrixBase<LinkedVectorMatrix<_Scalar,_Flags> >
#else
: public SparseMatrixBase
#endif
{
public:
EIGEN_GENERIC_PUBLIC_INTERFACE(LinkedVectorMatrix)

View File

@ -52,7 +52,12 @@ struct ei_traits<SparseMatrix<_Scalar, _Flags> >
template<typename _Scalar, int _Flags>
class SparseMatrix : public SparseMatrixBase<SparseMatrix<_Scalar, _Flags> >
class SparseMatrix
#ifndef EIGEN_PARSED_BY_DOXYGEN
: public SparseMatrixBase<SparseMatrix<_Scalar, _Flags> >
#else
: public SparseMatrixBase
#endif
{
public:
EIGEN_GENERIC_PUBLIC_INTERFACE(SparseMatrix)
@ -92,7 +97,7 @@ class SparseMatrix : public SparseMatrixBase<SparseMatrix<_Scalar, _Flags> >
return m_data.value(end-1);
// ^^ optimization: let's first check if it is the last coefficient
// (very common in high level algorithms)
const int* r = std::lower_bound(&m_data.index(start),&m_data.index(end),inner);
const int id = r-&m_data.index(0);
return ((*r==inner) && (id<end)) ? m_data.value(id) : Scalar(0);

View File

@ -88,7 +88,11 @@ struct ei_traits<Product<LhsNested, RhsNested, SparseProduct> >
};
template<typename LhsNested, typename RhsNested> class Product<LhsNested,RhsNested,SparseProduct> : ei_no_assignment_operator,
#ifndef EIGEN_PARSED_BY_DOXYGEN
public MatrixBase<Product<LhsNested, RhsNested, SparseProduct> >
#else
public MatrixBase
#endif
{
public:

View File

@ -4,7 +4,7 @@ namespace Eigen {
<h1>Customizing Eigen</h1>
Eigen2 can be extended in several way, for instance, by defining global methods, \link ExtendingMatrixBase by adding custom methods to MatrixBase \endlink, etc.
Eigen2 can be extended in several way, for instance, by defining global methods, \ref ExtendingMatrixBase "by adding custom methods to MatrixBase", etc.
\section ExtendingMatrixBase Extending MatrixBase

View File

@ -2,11 +2,17 @@ namespace Eigen {
o /** \mainpage Eigen
This is the API documentation for Eigen.
<div class="eimainmenu">\b Overview
| \ref TutorialCore "Core features"
| \ref TutorialGeometry "Geometry"
| \ref TutorialAdvancedLinearAlgebra "Advanced linear algebra"
</div>
This is the API documentation for Eigen.
Most of the API is available as methods in MatrixBase, so this is a good starting point for browsing. Also have a look at Matrix, as a few methods and the matrix constructors are there. Other notable classes for the Eigen API are Cwise, which contains the methods for doing certain coefficient-wise operations, and Part.
For a first contact with Eigen, the best place is to have a look at the \ref QuickStartGuide "quick start guide". Then, it is enough to look at Matrix, MatrixBase, and Cwise. In fact, except for advanced use, the only class that you'll have to explicitly name in your program, i.e. of which you'll explicitly contruct objects, is Matrix. For instance, vectors are handled as a special case of Matrix with one column. Typedefs are provided, e.g. Vector2f is a typedef for Matrix<float, 2, 1>. Finally, you might also have look at the \ref ExampleList "the list of selected examples".
For a first contact with Eigen, the best place is to have a look at the \ref TutorialCore "tutorial". Then, it is enough to look at Matrix, MatrixBase, and Cwise. In fact, except for advanced use, the only class that you'll have to explicitly name in your program, i.e. of which you'll explicitly contruct objects, is Matrix. For instance, vectors are handled as a special case of Matrix with one column. Typedefs are provided, e.g. Vector2f is a typedef for Matrix<float, 2, 1>. Finally, you might also have look at the \ref ExampleList "the list of selected examples".
Most of the other classes are just return types for MatrixBase methods.

View File

@ -1,34 +1,31 @@
namespace Eigen {
/** \page QuickStartGuide
/** \page TutorialCore Tutorial 1/3 - Core features
\ingroup Tutorial
<h1>Quick start guide</h1>
<div class="eimainmenu">\ref index "Overview"
| \b Core \b features
| \ref TutorialGeometry "Geometry"
| \ref TutorialAdvancedLinearAlgebra "Advanced linear algebra"
</div>
\b Table \b of \b contents
- Core features (Chapter I)
- \ref SimpleExampleFixedSize
- \ref SimpleExampleDynamicSize
- \ref MatrixTypes
- \ref MatrixInitialization
- \ref BasicLinearAlgebra
- \ref Reductions
- \ref SubMatrix
- \ref MatrixTransformations
- \ref TriangularMatrix
- \ref Performance
- \ref Geometry (Chapter II)
- \ref AdvancedLinearAlgebra (Chapter III)
- \ref LinearSolvers
- \ref LU
- \ref Cholesky
- \ref QR
- \ref EigenProblems
- \ref TutorialCoreSimpleExampleFixedSize
- \ref TutorialCoreSimpleExampleDynamicSize
- \ref TutorialCoreMatrixTypes
- \ref TutorialCoreMatrixInitialization
- \ref TutorialCoreBasicLinearAlgebra
- \ref TutorialCoreReductions
- \ref TutorialCoreSubMatrix
- \ref TutorialCoreMatrixTransformations
- \ref TutorialCoreTriangularMatrix
- \ref TutorialCorePerformance
</br>
\n
<hr>
\section SimpleExampleFixedSize Simple example with fixed-size matrices and vectors
\section TutorialCoreSimpleExampleFixedSize Simple example with fixed-size matrices and vectors
By fixed-size, we mean that the number of rows and columns are known at compile-time. In this case, Eigen avoids dynamic memory allocation and unroll loops. This is useful for very small sizes (typically up to 4x4).
@ -40,7 +37,9 @@ output:
\include Tutorial_simple_example_fixed_size.out
</td></tr></table>
<a href="#" class="top">top</a>\section SimpleExampleDynamicSize Simple example with dynamic-size matrices and vectors
<a href="#" class="top">top</a>\section TutorialCoreSimpleExampleDynamicSize Simple example with dynamic-size matrices and vectors
Dynamic-size means that the number of rows and columns are not known at compile-time. In this case, they are stored as runtime variables and the arrays are dynamically allocated.
@ -57,7 +56,7 @@ output:
<a href="#" class="top">top</a>\section MatrixTypes Matrix and vector types
<a href="#" class="top">top</a>\section TutorialCoreMatrixTypes Matrix and vector types
In Eigen, all kinds of dense matrices and vectors are represented by the template class Matrix. In most cases you can simply use one of the \ref matrixtypedefs "several convenient typedefs".
@ -66,7 +65,7 @@ The template class Matrix takes a number of template parameters, but for now it
\code Matrix<Scalar, RowsAtCompileTime, ColsAtCompileTime> \endcode
\li \c Scalar is the scalar type, i.e. the type of the coefficients. That is, if you want a vector of floats, choose \c float here.
\li \c RowsAtCompileTime and \c ColsAtCompileTime are the number of rows and columns of the matrix as known at compile-time.
\li \c RowsAtCompileTime and \c ColsAtCompileTime are the number of rows and columns of the matrix as known at compile-time.
For example, \c Vector3d is a typedef for \code Matrix<double, 3, 1> \endcode
@ -76,8 +75,10 @@ What if the matrix has dynamic-size i.e. the number of rows or cols isn't known
<a href="#" class="top">top</a>\section MatrixInitialization Matrix and vector creation and initialization
<a href="#" class="top">top</a>\section TutorialCoreMatrixInitialization Matrix and vector creation and initialization
\subsection TutorialPredefMat PredefinedMatrix
Eigen offers several methods to create or set matrices with coefficients equals to either a constant value, the identity matrix or even random values:
<table class="tutorial_code">
<tr>
@ -138,6 +139,16 @@ x.setRandom(size);
\endcode
</td>
</tr>
<tr><td colspan="3">Basis vectors \link MatrixBase::Unit [details]\endlink</td></tr>
<tr><td>\code
Vector3f::UnixX() // 1 0 0
Vector3f::UnixY() // 0 1 0
Vector3f::UnixZ() // 0 0 1
\endcode</td><td></td><td>\code
VectorXf::Unit(size,i)
VectorXf::Unit(4,1) == Vector4f(0,1,0,0)
== Vector4f::UnitY()
\endcode
</table>
Here is an usage example:
@ -159,7 +170,25 @@ v = 6 6 6
\endcode
</td></tr></table>
Eigen also offer a comma initializer syntax which allows to set all the coefficients of a matrix to specific values:
\subsection TutorialMap Map
Any memory buffer can be mapped as an Eigen's expression:
<table class="tutorial_code"><tr><td>
\code
std::vector<float> stlarray(10);
Map<VectorXf>(&stlarray[0], stlarray.size()).setOnes();
int data[4] = 1, 2, 3, 4;
Matrix2i mat2x2(data);
MatrixXi mat2x2 = Map<Matrix2i>(data);
MatrixXi mat2x2 = Map<MatrixXi>(data,2,2);
\endcode
</td></tr></table>
\subsection TutorialCommaInit CommaInitializer
Eigen also offer a comma initializer syntax which allows you to set all the coefficients of a matrix to specific values:
<table class="tutorial_code"><tr><td>
\include Tutorial_commainit_01.cpp
</td>
@ -177,16 +206,16 @@ output:
\verbinclude Tutorial_commainit_02.out
</td></tr></table>
<p class="note">\b Side \b note: here .finished() is used to get the actual matrix object once the comma initialization
<span class="note">\b Side \b note: here .finished() is used to get the actual matrix object once the comma initialization
of our temporary submatrix is done. Note that despite the appearant complexity of such an expression
Eigen's comma initializer usually yields to very optimized code without any overhead.</p>
Eigen's comma initializer usually yields to very optimized code without any overhead.</span>
<a href="#" class="top">top</a>\section BasicLinearAlgebra Basic Linear Algebra
<a href="#" class="top">top</a>\section TutorialCoreBasicLinearAlgebra Basic Linear Algebra
In short all mathematically well defined operators can be used right away as in the following example:
\code
@ -221,7 +250,7 @@ outer product</td><td>\code
mat = vec1 * vec2.transpose();\endcode
</td></tr>
<tr><td>
\link MatrixBase::cross() cross product \endcode</td><td>\code
\link MatrixBase::cross() cross product \endlink</td><td>\code
#include <Eigen/Geometry>
vec3 = vec1.cross(vec2);\endcode</td></tr>
</table>
@ -287,18 +316,18 @@ mat3 = mat1.cwise().abs2(mat2);
</table>
</td></tr></table>
<p class="note">\b Side \b note: If you feel the \c .cwise() syntax is too verbose for your taste and don't bother to have non mathematical operator directly available feel free to extend MatrixBase as described \ref ExtendingMatrixBase "here".</p>
<span class="note">\b Side \b note: If you feel the \c .cwise() syntax is too verbose for your taste and don't bother to have non mathematical operator directly available feel free to extend MatrixBase as described \ref ExtendingMatrixBase "here".</span>
<a href="#" class="top">top</a>\section Reductions Reductions
<a href="#" class="top">top</a>\section TutorialCoreReductions Reductions
Eigen provides several several reduction methods such as:
\link Cwise::minCoeff() minCoeff() \endlink, \link Cwise::maxCoeff() maxCoeff() \endlink,
\link Cwise::sum() sum() \endlink, \link Cwise::trace() trace() \endlink,
\link Cwise::norm() norm() \endlink, \link Cwise::norm2() norm2() \endlink,
\link Cwise::all() all() \endlink,and \link Cwise::any() any() \endlink.
\link MatrixBase::minCoeff() minCoeff() \endlink, \link MatrixBase::maxCoeff() maxCoeff() \endlink,
\link MatrixBase::sum() sum() \endlink, \link MatrixBase::trace() trace() \endlink,
\link MatrixBase::norm() norm() \endlink, \link MatrixBase::norm2() norm2() \endlink,
\link MatrixBase::all() all() \endlink,and \link MatrixBase::any() any() \endlink.
All reduction operations can be done matrix-wise,
\link MatrixBase::colwise() column-wise \endlink or
\link MatrixBase::rowwise() row-wise \endlink. Usage example:
@ -316,13 +345,13 @@ mat = 2 7 8
\endcode</td></tr>
</table>
<p class="note">\b Side \b note: The all() and any() functions are especially useful in combinaison with coeff-wise comparison operators (\ref CwiseAll "example").</p>
<span class="note">\b Side \b note: The all() and any() functions are especially useful in combinaison with coeff-wise comparison operators (\ref CwiseAll "example").</span>
<a href="#" class="top">top</a>\section SubMatrix Sub matrices
<a href="#" class="top">top</a>\section TutorialCoreSubMatrix Sub matrices
Read-write access to a \link MatrixBase::col(int) column \endlink
or a \link MatrixBase::row(int) row \endlink of a matrix:
@ -383,7 +412,7 @@ Read-write access to sub-matrices:
<a href="#" class="top">top</a>\section MatrixTransformations Matrix transformations
<a href="#" class="top">top</a>\section TutorialCoreMatrixTransformations Matrix transformations
<table class="tutorial_code">
<tr><td>
@ -411,11 +440,11 @@ mat3x3 = mat4x4.minor(i,j);\endcode
</table>
<a href="#" class="top">top</a>\section TriangularMatrix Dealing with triangular matrices
<a href="#" class="top">top</a>\section TutorialCoreTriangularMatrix Dealing with triangular matrices
todo
<a href="#" class="top">top</a>\section Performance Notes on performances
<a href="#" class="top">top</a>\section TutorialCorePerformance Notes on performances
<table class="tutorial_code">
<tr><td>\code
@ -447,17 +476,70 @@ m4 = m4 * m4.transpose().eval();\endcode</td><td>
forces immediate evaluation of the transpose</td></tr>
</table>
<a href="#" class="top">top</a>\section Geometry Geometry features
*/
maybe a second chapter for that
<a href="#" class="top">top</a>\section AdvancedLinearAlgebra Advanced Linear Algebra
Again, let's do another chapter for that
\subsection LinearSolvers Solving linear problems
\subsection LU LU
\subsection Cholesky Cholesky
\subsection QR QR
\subsection EigenProblems Eigen value problems
/** \page TutorialGeometry Tutorial 2/3 - Geometry
\ingroup Tutorial
<div class="eimainmenu">\ref index "Overview"
| \ref TutorialCore "Core features"
| \b Geometry
| \ref TutorialAdvancedLinearAlgebra "Advanced linear algebra"
</div>
\b Table \b of \b contents
- \ref TutorialGeoRotations
- \ref TutorialGeoTransformation
<a href="#" class="top">top</a>\section TutorialGeoRotations 2D and 3D Rotations
todo
<a href="#" class="top">top</a>\section TutorialGeoTransformation 2D and 3D Transformations
todo
*/
/** \page TutorialAdvancedLinearAlgebra Tutorial 3/3 - Advanced linear algebra
\ingroup Tutorial
<div class="eimainmenu">\ref index "Overview"
| \ref TutorialCore "Core features"
| \ref TutorialGeometry "Geometry"
| \b Advanced \b linear \b algebra
</div>
\b Table \b of \b contents
- \ref TutorialAdvLinearSolvers
- \ref TutorialAdvLU
- \ref TutorialAdvCholesky
- \ref TutorialAdvQR
- \ref TutorialAdvEigenProblems
\section TutorialAdvLinearSolvers Solving linear problems
todo
<a href="#" class="top">top</a>\section TutorialAdvLU LU
todo
<a href="#" class="top">top</a>\section TutorialAdvCholesky Cholesky
todo
<a href="#" class="top">top</a>\section TutorialAdvQR QR
todo
<a href="#" class="top">top</a>\section TutorialAdvEigenProblems Eigen value problems
todo
*/

View File

@ -451,10 +451,17 @@ A.top {
A.top:hover, A.logo:hover {
background-color: transparent;font-weight : bolder;
}
P.note {
SPAN.note {
font-size: 7pt;
}
DIV.navigation {
min-height : 64px;
/* background : url("Eigen_Silly_Professor_64x64.png") no-repeat left; */
padding-left : 80px;
padding-top : 5px;
}
TABLE.noborder {
border-bottom-style : none;
border-left-style : none;
@ -496,5 +503,10 @@ TABLE.tutorial_code TD {
vertical-align: middle;
}
DIV.eimainmenu {
text-align: center;
/* border-top: solid; */
/* border-bottom: solid; */
}

View File

@ -7,4 +7,4 @@
<a name="top"></a>
<a class="logo" href="http://eigen.tuxfamily.org/">
<img src="Eigen_Silly_Professor_64x64.png" width=64 height=64 alt="Eiegn's silly professor"
style="border:none" /></a>
style="position:absolute; border:none" /></a>

View File

@ -2,13 +2,9 @@
DIV.tabs
{
position: absolute;
right: 10pt;
left: 0px;
background : url("tab_b.gif") repeat-x bottom;
margin-bottom : 4px;
margin-left : 100px;
margin-top : -2em;
float: left;
width:100%;
background : url("tab_b.gif") repeat-x bottom;
}
DIV.tabs UL