MADNESS  0.10.1
Public Member Functions | Private Types | Private Attributes | Friends | List of all members
madness::TensorTrain< T > Class Template Reference

#include <tensortrain.h>

Inheritance diagram for madness::TensorTrain< T >:
Inheritance graph
[legend]
Collaboration diagram for madness::TensorTrain< T >:
Collaboration graph
[legend]

Public Member Functions

 TensorTrain ()
 empty constructor More...
 
 TensorTrain (const long &ndims, const long *dims)
 ctor for a TensorTrain, set up only the dimensions, no data More...
 
 TensorTrain (const std::vector< long > &dims)
 ctor for a TensorTrain, set up only the dimensions, no data More...
 
 TensorTrain (const std::vector< Tensor< T > > &c)
 ctor for a TensorTrain, with core tensors explicitly given More...
 
 TensorTrain (const Tensor< T > &t, double eps)
 ctor for a TensorTrain, with the tolerance eps More...
 
 TensorTrain (const Tensor< T > &t, double eps, const std::vector< long > dims)
 ctor for a TensorTrain, with the tolerance eps More...
 
 TensorTrain (const TensorTrain &other)
 copy constructor, shallow More...
 
void decompose (const Tensor< T > &t, double eps, const std::vector< long > &dims)
 decompose the input tensor into a TT representation More...
 
TensorTrain< T > & emul (const TensorTrain< T > &other)
 compute the Hadamard product of two TensorTrains More...
 
void fusedim (const long i)
 merge two dimensions into one More...
 
TensorTrain< T > & gaxpy (const std::array< Slice, TENSOR_MAXDIM > &s1, const TensorTrain< T > &rhs, T beta, const std::array< Slice, TENSOR_MAXDIM > &s2)
 Inplace generalized saxpy with slices and without alpha. More...
 
TensorTrain< T > & gaxpy (T alpha, const TensorTrain< T > &rhs, T beta)
 Inplace generalized saxpy ... this = this*alpha + other*beta. More...
 
Tensor< T > & get_core (const int i)
 reference to the internal core More...
 
const Tensor< T > & get_core (const int i) const
 const reference to the internal core More...
 
bool is_operator () const
 check if this is an operator (r,k',k,r) More...
 
bool is_tensor () const
 check if this is a tensor (r,k,r) More...
 
bool is_zero_rank () const
 if rank is zero More...
 
TensorTrain< T > & make_operator ()
 convert this into an operator representation (r,k',k,r) More...
 
TensorTrain< T > & make_tensor ()
 convert this into a tensor representation (r,k,r) More...
 
float_scalar_type normf () const
 returns the Frobenius norm More...
 
template<class Q >
 operator TensorTrain< Q > () const
 Type conversion makes a deep copy. More...
 
TensorTrain< Toperator* (const T &factor) const
 return this multiplied by a scalar More...
 
TensorTrain< T > & operator+= (const TensorTrain< T > &rhs)
 inplace addition of two Tensortrains; will increase ranks of this More...
 
TensorTrain< T > & operator-= (const TensorTrain< T > &rhs)
 inplace subtraction of two Tensortrains; will increase ranks of this More...
 
TensorTrain< T > & operator= (const T &number)
 assign a number to this tensor More...
 
TensorTrainoperator= (const TensorTrain &other)
 assigment operator More...
 
Tptr (const int ivec=0)
 Returns a pointer to the internal data. More...
 
const Tptr (const int ivec=0) const
 Returns a pointer to the internal data. More...
 
std::vector< long > ranks () const
 return the TT ranks More...
 
long ranks (const int i) const
 return the TT ranks for dimension i (to i+1) More...
 
long real_size () const
 return the size of this instance, including static memory for vectors and such More...
 
Tensor< Treconstruct (const bool flat=false) const
 reconstruct this to a full representation More...
 
void scale (T fac)
 scale this by a number More...
 
template<typename Archive >
void serialize (Archive &ar)
 serialize this More...
 
void set_size_and_dim (std::vector< long > dims)
 
long size () const
 return the number of coefficients in all core tensors More...
 
TensorTrain< Tsplitdim (long idim, long k1, long k2, const double eps) const
 
template<class Q >
std::enable_if<(TensorTypeData< T >::iscomplex or TensorTypeData< Q >::iscomplex), TENSOR_RESULT_TYPE(T, Q)>::type trace (const TensorTrain< Q > &B) const
 
template<class Q >
std::enable_if<!(TensorTypeData< T >::iscomplex or TensorTypeData< Q >::iscomplex), TENSOR_RESULT_TYPE(T, Q)>::type trace (const TensorTrain< Q > &B) const
 Return the trace of two tensors, no complex conjugate involved. More...
 
template<typename R = T>
std::enable_if<!std::is_arithmetic< R >::value, void >::type truncate (double eps)
 recompress and truncate this TT representation More...
 
template<typename R = T>
std::enable_if< std::is_arithmetic< R >::value, void >::type truncate (double eps)
 recompress and truncate this TT representation More...
 
void two_mode_representation (Tensor< T > &U, Tensor< T > &VT, Tensor< typename Tensor< T >::scalar_type > &s) const
 construct a two-mode representation (aka unnormalized SVD) More...
 
bool verify () const
 return the dimensions of this tensor More...
 
void zero_me ()
 turn this into an empty tensor with all cores properly shaped More...
 
void zero_me (const std::vector< long > &dim)
 turn this into an empty tensor with all cores properly shaped More...
 
- Public Member Functions inherited from madness::BaseTensor
 BaseTensor ()
 
virtual ~BaseTensor ()
 
bool conforms (const BaseTensor *t) const
 Returns true if this and *t are the same shape and size. More...
 
long dim (int i) const
 Returns the size of dimension i. More...
 
const long * dims () const
 Returns the array of tensor dimensions. More...
 
long id () const
 Returns the typeid of the tensor (c.f., TensorTypeData<T> ) More...
 
bool iscontiguous () const
 Returns true if the tensor refers to contiguous memory locations. More...
 
long ndim () const
 Returns the number of dimensions in the tensor. More...
 
long size () const
 Returns the number of elements in the tensor. More...
 
long stride (int i) const
 Returns the stride associated with dimension i. More...
 
const long * strides () const
 Returns the array of tensor strides. More...
 

Private Types

typedef TensorTypeData< T >::float_scalar_type float_scalar_type
 C++ typename of the floating point type associated with scalar real type. More...
 
typedef TensorTypeData< T >::scalar_type scalar_type
 C++ typename of the real type associated with a complex type. More...
 
typedef T type
 C++ typename of this tensor. More...
 

Private Attributes

std::vector< Tensor< T > > core
 
bool zero_rank
 true if rank is zero More...
 

Friends

template<typename Q >
class TensorTrain
 
TensorTrain copy (const TensorTrain &other)
 deep copy of the whole tensor More...
 
TensorTrain copy (const TensorTrain &other, const std::array< Slice, TENSOR_MAXDIM > &s)
 deep copy of a slice of the tensor More...
 
template<typename R , typename Q >
TensorTrain< TENSOR_RESULT_TYPE(R, Q)> general_transform (const TensorTrain< R > &t, const Tensor< Q > c[])
 
template<typename R , typename Q >
TensorTrain< TENSOR_RESULT_TYPE(R, Q)> outer (const TensorTrain< R > &t1, const TensorTrain< Q > &t2)
 
template<typename R , typename Q >
TensorTrain< TENSOR_RESULT_TYPE(R, Q)> transform (const TensorTrain< R > &t, const Tensor< Q > &c)
 
template<typename R , typename Q >
TensorTrain< TENSOR_RESULT_TYPE(R, Q)> transform_dir (const TensorTrain< R > &t, const Tensor< Q > &c, const int axis)
 

Additional Inherited Members

- Static Public Member Functions inherited from madness::BaseTensor
static int get_instance_count ()
 Returns the count of all current instances of tensors & slice tensors of all types. More...
 
- Protected Member Functions inherited from madness::BaseTensor
void cycledim_inplace (long shift, long start, long end)
 Cyclic shift of dimensions. More...
 
void flat_inplace ()
 Reshapes the tensor inplace into 1D. More...
 
void fusedim_inplace (long i)
 Fuses dimensions i and i+1. More...
 
void mapdim_inplace (const std::vector< long > &map)
 General permutation of dimensions. More...
 
void reshape_inplace (const std::vector< long > &d)
 Reshapes the tensor inplace. More...
 
void reshape_inplace (int ndimnew, const long *d)
 Reshapes the tensor inplace. More...
 
void set_dims_and_size (long nd, const long d[])
 
void splitdim_inplace (long i, long dimi0, long dimi1)
 Splits dimension i. More...
 
void swapdim_inplace (long i, long j)
 Swaps the dimensions. More...
 
- Protected Attributes inherited from madness::BaseTensor
long _dim [TENSOR_MAXDIM]
 Size of each dimension. More...
 
long _id
 Id from TensorTypeData<T> in type_data.h. More...
 
long _ndim
 Number of dimensions (-1=invalid; 0=no supported; >0=tensor) More...
 
long _size
 Number of elements in the tensor. More...
 
long _stride [TENSOR_MAXDIM]
 Increment between elements in each dimension. More...
 

Detailed Description

template<typename T>
class madness::TensorTrain< T >

A tensor train is a multi-modal representation of a tensor t

t(i,j,k,l) = \sum G^1_{a1,i,a2} G^2_{a2,j,a3} G^3_{a3,k,a4} G^4_{a4,l,a5}
Definition: test_derivative.cc:24
static const Slice _(0,-1, 1)
Function< T, NDIM > sum(World &world, const std::vector< Function< T, NDIM > > &f, bool fence=true)
Returns new function — q = sum_i f[i].
Definition: vmra.h:1421
static const long k
Definition: rk.cc:44
const double a2
Definition: vnucso.cc:86
const double a1
Definition: vnucso.cc:85

The "core" tensors G are connected via a linear index network, where the first index a1 and the last index a5 are boundary indices and are set to 1.

The tensor train representation is suited for any number of dimensions and in general at least as fast as the 2-way decomposition SVD. If the tensor has full rank it will need about twice the storage space of the full tensor

Member Typedef Documentation

◆ float_scalar_type

template<typename T >
typedef TensorTypeData<T>::float_scalar_type madness::TensorTrain< T >::float_scalar_type
private

C++ typename of the floating point type associated with scalar real type.

◆ scalar_type

template<typename T >
typedef TensorTypeData<T>::scalar_type madness::TensorTrain< T >::scalar_type
private

C++ typename of the real type associated with a complex type.

◆ type

template<typename T >
typedef T madness::TensorTrain< T >::type
private

C++ typename of this tensor.

Constructor & Destructor Documentation

◆ TensorTrain() [1/7]

template<typename T >
madness::TensorTrain< T >::TensorTrain ( )
inline

empty constructor

References madness::TensorTrain< T >::set_size_and_dim().

◆ TensorTrain() [2/7]

template<typename T >
madness::TensorTrain< T >::TensorTrain ( const Tensor< T > &  t,
double  eps 
)
inline

ctor for a TensorTrain, with the tolerance eps

The tensor train will represent the input tensor with accuracy || t - this ||_2 < eps

Note that we rely on specific layout of the memory in the tensors, e.g. we pass SliceTensors on to lapack. This will only work if the slices are contiguous.

Parameters
[in]tfull representation of a tensor
[in]epsthe accuracy threshold

References d(), madness::TensorTrain< T >::decompose(), madness::BaseTensor::dim(), madness::BaseTensor::dims(), madness::Tensor< T >::flat(), MADNESS_ASSERT, madness::BaseTensor::ndim(), madness::BaseTensor::set_dims_and_size(), and madness::BaseTensor::size().

◆ TensorTrain() [3/7]

template<typename T >
madness::TensorTrain< T >::TensorTrain ( const Tensor< T > &  t,
double  eps,
const std::vector< long >  dims 
)
inline

ctor for a TensorTrain, with the tolerance eps

The tensor train will represent the input tensor with accuracy || t - this ||_2 < eps

Note that we rely on specific layout of the memory in the tensors, e.g. we pass SliceTensors on to lapack. This will only work if the slices are contiguous.

Parameters
[in]tfull representation of a tensor
[in]epsthe accuracy threshold
[in]dimsthe tt structure

References madness::TensorTrain< T >::decompose(), madness::BaseTensor::dims(), MADNESS_ASSERT, madness::BaseTensor::ndim(), madness::BaseTensor::set_dims_and_size(), and madness::BaseTensor::size().

◆ TensorTrain() [4/7]

template<typename T >
madness::TensorTrain< T >::TensorTrain ( const long &  ndims,
const long *  dims 
)
inline

◆ TensorTrain() [5/7]

template<typename T >
madness::TensorTrain< T >::TensorTrain ( const std::vector< long > &  dims)
inline

◆ TensorTrain() [6/7]

template<typename T >
madness::TensorTrain< T >::TensorTrain ( const std::vector< Tensor< T > > &  c)
inline

ctor for a TensorTrain, with core tensors explicitly given

the core tensors must have the shape (k1,r1) (r1,k2,r2) .. (r2,k3)

Parameters
[in]corevector of core tensors, properly shaped

References c, madness::TensorTrain< T >::core, d(), madness::BaseTensor::dim(), madness::TensorTrain< T >::set_size_and_dim(), madness::TensorTrain< T >::zero_me(), and madness::TensorTrain< T >::zero_rank.

◆ TensorTrain() [7/7]

template<typename T >
madness::TensorTrain< T >::TensorTrain ( const TensorTrain< T > &  other)
inline

Member Function Documentation

◆ decompose()

template<typename T >
void madness::TensorTrain< T >::decompose ( const Tensor< T > &  t,
double  eps,
const std::vector< long > &  dims 
)
inline

decompose the input tensor into a TT representation

Parameters
[in]ttensor in full rank
[in]epsthe precision threshold
[in]dimsthe tt structure

References madness::_(), aa, b, c, madness::copy(), madness::TensorTrain< T >::copy, madness::TensorTrain< T >::core, d(), madness::BaseTensor::dims(), k, max, madness::SRConf< T >::max_sigma(), madness::print(), madness::BaseTensor::size(), madness::svd_result(), T(), u(), and madness::TensorTrain< T >::zero_rank.

Referenced by madness::TensorTrain< T >::TensorTrain().

◆ emul()

template<typename T >
TensorTrain<T>& madness::TensorTrain< T >::emul ( const TensorTrain< T > &  other)
inline

◆ fusedim()

template<typename T >
void madness::TensorTrain< T >::fusedim ( const long  i)
inline

merge two dimensions into one

merge dimension i and i+1 into new dimension i

Parameters
[in]ithe first dimension

(r1=0, k1*k2, r2=0)

References madness::TensorTrain< T >::core, d(), madness::BaseTensor::dim(), madness::BaseTensor::fusedim_inplace(), madness::inner(), and madness::TensorTrain< T >::zero_rank.

Referenced by madness::TensorTrain< T >::emul(), and madness::TensorTrain< T >::make_tensor().

◆ gaxpy() [1/2]

template<typename T >
TensorTrain<T>& madness::TensorTrain< T >::gaxpy ( const std::array< Slice, TENSOR_MAXDIM > &  s1,
const TensorTrain< T > &  rhs,
T  beta,
const std::array< Slice, TENSOR_MAXDIM > &  s2 
)
inline

◆ gaxpy() [2/2]

template<typename T >
TensorTrain<T>& madness::TensorTrain< T >::gaxpy ( T  alpha,
const TensorTrain< T > &  rhs,
T  beta 
)
inline

◆ get_core() [1/2]

template<typename T >
Tensor<T>& madness::TensorTrain< T >::get_core ( const int  i)
inline

reference to the internal core

References madness::TensorTrain< T >::core, MADNESS_ASSERT, and madness::BaseTensor::ndim().

Referenced by madness::apply(), and madness::tt_identity().

◆ get_core() [2/2]

template<typename T >
const Tensor<T>& madness::TensorTrain< T >::get_core ( const int  i) const
inline

const reference to the internal core

References madness::TensorTrain< T >::core, MADNESS_ASSERT, and madness::BaseTensor::ndim().

◆ is_operator()

template<typename T >
bool madness::TensorTrain< T >::is_operator ( ) const
inline

check if this is an operator (r,k',k,r)

References madness::TensorTrain< T >::is_tensor().

Referenced by madness::TensorTrain< T >::make_operator().

◆ is_tensor()

template<typename T >
bool madness::TensorTrain< T >::is_tensor ( ) const
inline

◆ is_zero_rank()

template<typename T >
bool madness::TensorTrain< T >::is_zero_rank ( ) const
inline

if rank is zero

References madness::TensorTrain< T >::zero_rank.

◆ make_operator()

template<typename T >
TensorTrain<T>& madness::TensorTrain< T >::make_operator ( )
inline

◆ make_tensor()

template<typename T >
TensorTrain<T>& madness::TensorTrain< T >::make_tensor ( )
inline

◆ normf()

template<typename T >
float_scalar_type madness::TensorTrain< T >::normf ( ) const
inline

returns the Frobenius norm

References std::abs(), and madness::TensorTrain< T >::trace().

◆ operator TensorTrain< Q >()

template<typename T >
template<class Q >
madness::TensorTrain< T >::operator TensorTrain< Q > ( ) const
inline

◆ operator*()

template<typename T >
TensorTrain<T> madness::TensorTrain< T >::operator* ( const T factor) const
inline

return this multiplied by a scalar

Returns
new tensor

References madness::TensorTrain< T >::copy, and madness::TensorTrain< T >::scale().

◆ operator+=()

template<typename T >
TensorTrain<T>& madness::TensorTrain< T >::operator+= ( const TensorTrain< T > &  rhs)
inline

inplace addition of two Tensortrains; will increase ranks of this

inefficient if many additions are performed, since it requires many calls of new.

Parameters
[in]rhsa TensorTrain to be added

References madness::TensorTrain< T >::gaxpy().

◆ operator-=()

template<typename T >
TensorTrain<T>& madness::TensorTrain< T >::operator-= ( const TensorTrain< T > &  rhs)
inline

inplace subtraction of two Tensortrains; will increase ranks of this

inefficient if many subtractions are performed, since it requires many calls of new.

Parameters
[in]rhsa TensorTrain to be added

References madness::TensorTrain< T >::gaxpy().

◆ operator=() [1/2]

template<typename T >
TensorTrain<T>& madness::TensorTrain< T >::operator= ( const T number)
inline

◆ operator=() [2/2]

template<typename T >
TensorTrain& madness::TensorTrain< T >::operator= ( const TensorTrain< T > &  other)
inline

◆ ptr() [1/2]

template<typename T >
T* madness::TensorTrain< T >::ptr ( const int  ivec = 0)
inline

Returns a pointer to the internal data.

Parameters
[in]ivecindex of core vector to which the return values points

References madness::TensorTrain< T >::core.

◆ ptr() [2/2]

template<typename T >
const T* madness::TensorTrain< T >::ptr ( const int  ivec = 0) const
inline

Returns a pointer to the internal data.

Parameters
[in]ivecindex of core vector to which the return values points

References madness::TensorTrain< T >::core.

◆ ranks() [1/2]

template<typename T >
std::vector<long> madness::TensorTrain< T >::ranks ( ) const
inline

◆ ranks() [2/2]

template<typename T >
long madness::TensorTrain< T >::ranks ( const int  i) const
inline

◆ real_size()

template<typename T >
long madness::TensorTrain< T >::real_size ( ) const
inline

return the size of this instance, including static memory for vectors and such

References madness::TensorTrain< T >::core, madness::TensorTrain< T >::size(), and T().

◆ reconstruct()

template<typename T >
Tensor<T> madness::TensorTrain< T >::reconstruct ( const bool  flat = false) const
inline

reconstruct this to a full representation

Parameters
[in]flatreturn this in flat representation
Returns
this in full rank representation

References madness::BaseTensor::_size, madness::TensorTrain< T >::core, madness::BaseTensor::dim(), madness::BaseTensor::dims(), madness::Tensor< T >::fusedim(), madness::inner(), MADNESS_ASSERT, madness::BaseTensor::ndim(), madness::TensorTrain< T >::size(), and madness::TensorTrain< T >::zero_rank.

Referenced by test_stuff().

◆ scale()

template<typename T >
void madness::TensorTrain< T >::scale ( T  fac)
inline

scale this by a number

Parameters
[in]facthe factor to multiply
Returns
*this * fac

References madness::TensorTrain< T >::core, and madness::TensorTrain< T >::zero_rank.

Referenced by madness::TensorTrain< T >::gaxpy(), and madness::TensorTrain< T >::operator*().

◆ serialize()

template<typename T >
template<typename Archive >
void madness::TensorTrain< T >::serialize ( Archive &  ar)
inline

◆ set_size_and_dim()

template<typename T >
void madness::TensorTrain< T >::set_size_and_dim ( std::vector< long >  dims)
inline

◆ size()

template<typename T >
long madness::TensorTrain< T >::size ( ) const
inline

◆ splitdim()

template<typename T >
TensorTrain<T> madness::TensorTrain< T >::splitdim ( long  idim,
long  k1,
long  k2,
const double  eps 
) const
inline

Returns new view/tensor splitting dimension i as dimi0*dimi1 to produce conforming d+1 dimension tensor

Parameters
[in]idimthe dimension to be split
[in]k1new first dimension of idim
[in]k2new second dimension of idim
[in]epsthreshold for SVD (choose negative to keep all terms)
Returns
new deep copy of this with split dimensions

References madness::TensorTrain< T >::TensorTrain, madness::_(), madness::TensorTrain< T >::copy, madness::TensorTrain< T >::core, madness::BaseTensor::dim(), madness::BaseTensor::dims(), k1, k2, MADNESS_ASSERT, madness::SRConf< T >::max_sigma(), madness::BaseTensor::ndim(), madness::r2(), madness::TensorTrain< T >::ranks(), madness::BaseTensor::splitdim_inplace(), madness::svd(), and madness::TensorTrain< T >::zero_rank.

◆ trace() [1/2]

template<typename T >
template<class Q >
std::enable_if<(TensorTypeData<T>::iscomplex or TensorTypeData<Q>::iscomplex), TENSOR_RESULT_TYPE(T,Q)>::type madness::TensorTrain< T >::trace ( const TensorTrain< Q > &  B) const
inline

◆ trace() [2/2]

template<typename T >
template<class Q >
std::enable_if<!(TensorTypeData<T>::iscomplex or TensorTypeData<Q>::iscomplex), TENSOR_RESULT_TYPE(T,Q)>::type madness::TensorTrain< T >::trace ( const TensorTrain< Q > &  B) const
inline

◆ truncate() [1/2]

template<typename T >
template<typename R = T>
std::enable_if<!std::is_arithmetic<R>::value, void>::type madness::TensorTrain< T >::truncate ( double  eps)
inline

recompress and truncate this TT representation

this in recompressed TT form with optimal rank

Parameters
[in]epsthe truncation threshold

References MADNESS_EXCEPTION.

Referenced by madness::SeparatedConvolution< Q, NDIM >::make_tt_representation().

◆ truncate() [2/2]

template<typename T >
template<typename R = T>
std::enable_if<std::is_arithmetic<R>::value, void>::type madness::TensorTrain< T >::truncate ( double  eps)
inline

◆ two_mode_representation()

template<typename T >
void madness::TensorTrain< T >::two_mode_representation ( Tensor< T > &  U,
Tensor< T > &  VT,
Tensor< typename Tensor< T >::scalar_type > &  s 
) const
inline

construct a two-mode representation (aka unnormalized SVD)

Parameters
[out]UThe left singular vectors, ({i},rank)
[out]VTThe right singular vectors, (rank,{i})
[out]sVector holding 1's, (rank)

References madness::TensorTrain< T >::core, madness::BaseTensor::dim(), madness::inner(), MADNESS_ASSERT, MADNESS_EXCEPTION, madness::BaseTensor::ndim(), and madness::TensorTrain< T >::zero_rank.

◆ verify()

template<typename T >
bool madness::TensorTrain< T >::verify ( ) const
inline

return the dimensions of this tensor

check that the ranks of all core tensors are consistent

References c, madness::TensorTrain< T >::core, d(), madness::BaseTensor::dim(), madness::BaseTensor::ndim(), and madness::TensorTrain< T >::size().

Referenced by madness::TensorTrain< T >::gaxpy(), and madness::TensorTrain< T >::truncate().

◆ zero_me() [1/2]

template<typename T >
void madness::TensorTrain< T >::zero_me ( )
inline

turn this into an empty tensor with all cores properly shaped

References madness::BaseTensor::dims(), and madness::BaseTensor::ndim().

Referenced by madness::TensorTrain< T >::TensorTrain(), and madness::TensorTrain< T >::truncate().

◆ zero_me() [2/2]

template<typename T >
void madness::TensorTrain< T >::zero_me ( const std::vector< long > &  dim)
inline

turn this into an empty tensor with all cores properly shaped

References madness::BaseTensor::dim().

Friends And Related Function Documentation

◆ TensorTrain

template<typename T >
template<typename Q >
friend class TensorTrain
friend

◆ copy [1/2]

template<typename T >
TensorTrain copy ( const TensorTrain< T > &  other)
friend

deep copy of the whole tensor

if argument has zero rank return a zero-rank tensor of the same dimensions

Referenced by madness::TensorTrain< T >::decompose(), madness::TensorTrain< T >::emul(), madness::TensorTrain< T >::operator*(), and madness::TensorTrain< T >::splitdim().

◆ copy [2/2]

template<typename T >
TensorTrain copy ( const TensorTrain< T > &  other,
const std::array< Slice, TENSOR_MAXDIM > &  s 
)
friend

deep copy of a slice of the tensor

this operation does not change the ranks, i.e. the resulting tensor is most likely not in an optimal compression state

Parameters
[in]othertensor to be sliced
[in]svector of slices

◆ general_transform

template<typename T >
template<typename R , typename Q >
TensorTrain<TENSOR_RESULT_TYPE(R,Q)> general_transform ( const TensorTrain< R > &  t,
const Tensor< Q c[] 
)
friend

◆ outer

template<typename T >
template<typename R , typename Q >
TensorTrain<TENSOR_RESULT_TYPE(R,Q)> outer ( const TensorTrain< R > &  t1,
const TensorTrain< Q > &  t2 
)
friend

◆ transform

template<typename T >
template<typename R , typename Q >
TensorTrain<TENSOR_RESULT_TYPE(R,Q)> transform ( const TensorTrain< R > &  t,
const Tensor< Q > &  c 
)
friend

◆ transform_dir

template<typename T >
template<typename R , typename Q >
TensorTrain<TENSOR_RESULT_TYPE(R,Q)> transform_dir ( const TensorTrain< R > &  t,
const Tensor< Q > &  c,
const int  axis 
)
friend

Member Data Documentation

◆ core

template<typename T >
std::vector<Tensor<T> > madness::TensorTrain< T >::core
private

◆ zero_rank

template<typename T >
bool madness::TensorTrain< T >::zero_rank
private

The documentation for this class was generated from the following file: