10 #ifndef EIGEN_ASSIGNMENT_FUNCTORS_H    11 #define EIGEN_ASSIGNMENT_FUNCTORS_H    21 template<
typename DstScalar,
typename SrcScalar> 
struct assign_op {
    26   template<
int Alignment, 
typename Packet>
    28   { internal::pstoret<DstScalar,Packet,Alignment>(a,b); }
    32 template<
typename DstScalar> 
struct assign_op<DstScalar,void> {};
    34 template<
typename DstScalar,
typename SrcScalar>
    51   template<
int Alignment, 
typename Packet>
    53   { internal::pstoret<DstScalar,Packet,Alignment>(a,
internal::padd(internal::ploadt<Packet,Alignment>(a),b)); }
    55 template<
typename DstScalar,
typename SrcScalar>
    72   template<
int Alignment, 
typename Packet>
    74   { internal::pstoret<DstScalar,Packet,Alignment>(a,
internal::psub(internal::ploadt<Packet,Alignment>(a),b)); }
    76 template<
typename DstScalar,
typename SrcScalar>
    88 template<
typename DstScalar, 
typename SrcScalar=DstScalar>
    94   template<
int Alignment, 
typename Packet>
    96   { internal::pstoret<DstScalar,Packet,Alignment>(a,
internal::pmul(internal::ploadt<Packet,Alignment>(a),b)); }
    98 template<
typename DstScalar, 
typename SrcScalar>
   110 template<
typename DstScalar, 
typename SrcScalar=DstScalar> 
struct div_assign_op {
   115   template<
int Alignment, 
typename Packet>
   117   { internal::pstoret<DstScalar,Packet,Alignment>(a,
internal::pdiv(internal::ploadt<Packet,Alignment>(a),b)); }
   119 template<
typename DstScalar, 
typename SrcScalar>
   149     Scalar t=b; 
const_cast<Scalar&
>(b)=a; a=t;
   152     swap(a,const_cast<Scalar&>(b));
   156 template<
typename Scalar>
   168 #endif // EIGEN_ASSIGNMENT_FUNCTORS_H 
EIGEN_STRONG_INLINE void assignPacket(DstScalar *a, const Packet &b) const
#define EIGEN_STRONG_INLINE
EIGEN_STRONG_INLINE void assignPacket(DstScalar *a, const Packet &b) const
#define EIGEN_EMPTY_STRUCT_CTOR(X)
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
EIGEN_DEVICE_FUNC Packet padd(const Packet &a, const Packet &b)
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void assignCoeff(DstScalar &a, const SrcScalar &b) const
EIGEN_STRONG_INLINE void assignPacket(DstScalar *a, const Packet &b) const
EIGEN_DEVICE_FUNC Packet pdiv(const Packet &a, const Packet &b)
EIGEN_STRONG_INLINE void assignPacket(DstScalar *a, const Packet &b) const
EIGEN_DEVICE_FUNC Packet psub(const Packet &a, const Packet &b)
EIGEN_STRONG_INLINE void assignPacket(DstScalar *a, const Packet &b) const
EIGEN_DEVICE_FUNC Packet pmul(const Packet &a, const Packet &b)
void swap(scoped_array< T > &a, scoped_array< T > &b)