mirror of
https://gitlab.com/libeigen/eigen.git
synced 2025-05-31 10:34:02 +08:00

* use them (big simplification in Assign.h) * axe (Inner|Outer)StrideAtCompileTime that were just introduced * ei_int_if_dynamic now asserts that the size is the expected one: adapt to that in Block.h * add rowStride() / colStride() in DenseBase * implement innerStride() / outerStride() everywhere needed
119 lines
3.8 KiB
C++
119 lines
3.8 KiB
C++
// This file is part of Eigen, a lightweight C++ template library
|
|
// for linear algebra.
|
|
//
|
|
// Copyright (C) 2008 Gael Guennebaud <g.gael@free.fr>
|
|
// Copyright (C) 2006-2008 Benoit Jacob <jacob.benoit.1@gmail.com>
|
|
//
|
|
// 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_NESTBYVALUE_H
|
|
#define EIGEN_NESTBYVALUE_H
|
|
|
|
/** \class NestByValue
|
|
*
|
|
* \brief Expression which must be nested by value
|
|
*
|
|
* \param ExpressionType the type of the object of which we are requiring nesting-by-value
|
|
*
|
|
* This class is the return type of MatrixBase::nestByValue()
|
|
* and most of the time this is the only way it is used.
|
|
*
|
|
* \sa MatrixBase::nestByValue()
|
|
*/
|
|
template<typename ExpressionType>
|
|
struct ei_traits<NestByValue<ExpressionType> > : public ei_traits<ExpressionType>
|
|
{};
|
|
|
|
template<typename ExpressionType> class NestByValue
|
|
: public ExpressionType::template MakeBase< NestByValue<ExpressionType> >::Type
|
|
{
|
|
public:
|
|
|
|
typedef typename ExpressionType::template MakeBase<NestByValue<ExpressionType> >::Type Base;
|
|
EIGEN_DENSE_PUBLIC_INTERFACE(NestByValue)
|
|
|
|
inline NestByValue(const ExpressionType& matrix) : m_expression(matrix) {}
|
|
|
|
inline int rows() const { return m_expression.rows(); }
|
|
inline int cols() const { return m_expression.cols(); }
|
|
inline int outerStride() const { return m_expression.outerStride(); }
|
|
inline int innerStride() const { return m_expression.innerStride(); }
|
|
|
|
inline const CoeffReturnType coeff(int row, int col) const
|
|
{
|
|
return m_expression.coeff(row, col);
|
|
}
|
|
|
|
inline Scalar& coeffRef(int row, int col)
|
|
{
|
|
return m_expression.const_cast_derived().coeffRef(row, col);
|
|
}
|
|
|
|
inline const CoeffReturnType coeff(int index) const
|
|
{
|
|
return m_expression.coeff(index);
|
|
}
|
|
|
|
inline Scalar& coeffRef(int index)
|
|
{
|
|
return m_expression.const_cast_derived().coeffRef(index);
|
|
}
|
|
|
|
template<int LoadMode>
|
|
inline const PacketScalar packet(int row, int col) const
|
|
{
|
|
return m_expression.template packet<LoadMode>(row, col);
|
|
}
|
|
|
|
template<int LoadMode>
|
|
inline void writePacket(int row, int col, const PacketScalar& x)
|
|
{
|
|
m_expression.const_cast_derived().template writePacket<LoadMode>(row, col, x);
|
|
}
|
|
|
|
template<int LoadMode>
|
|
inline const PacketScalar packet(int index) const
|
|
{
|
|
return m_expression.template packet<LoadMode>(index);
|
|
}
|
|
|
|
template<int LoadMode>
|
|
inline void writePacket(int index, const PacketScalar& x)
|
|
{
|
|
m_expression.const_cast_derived().template writePacket<LoadMode>(index, x);
|
|
}
|
|
|
|
operator const ExpressionType&() const { return m_expression; }
|
|
|
|
protected:
|
|
const ExpressionType m_expression;
|
|
};
|
|
|
|
/** \returns an expression of the temporary version of *this.
|
|
*/
|
|
template<typename Derived>
|
|
inline const NestByValue<Derived>
|
|
DenseBase<Derived>::nestByValue() const
|
|
{
|
|
return NestByValue<Derived>(derived());
|
|
}
|
|
|
|
#endif // EIGEN_NESTBYVALUE_H
|