MADNESS
0.10.1
|
#include <tensortrain.h>
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< T > | operator* (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... | |
TensorTrain & | operator= (const TensorTrain &other) |
assigment operator More... | |
T * | ptr (const int ivec=0) |
Returns a pointer to the internal data. More... | |
const T * | ptr (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< T > | reconstruct (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< T > | splitdim (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... | |
A tensor train is a multi-modal representation of a tensor t
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
|
private |
C++ typename of the floating point type associated with scalar real type.
|
private |
C++ typename of the real type associated with a complex type.
|
private |
C++ typename of this tensor.
|
inline |
empty constructor
References madness::TensorTrain< T >::set_size_and_dim().
|
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.
[in] | t | full representation of a tensor |
[in] | eps | the 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().
|
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.
[in] | t | full representation of a tensor |
[in] | eps | the accuracy threshold |
[in] | dims | the 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().
|
inline |
ctor for a TensorTrain, set up only the dimensions, no data
References madness::BaseTensor::dims(), madness::BaseTensor::set_dims_and_size(), and madness::TensorTrain< T >::zero_rank.
|
inline |
ctor for a TensorTrain, set up only the dimensions, no data
References madness::BaseTensor::dims(), madness::TensorTrain< T >::set_size_and_dim(), and madness::TensorTrain< T >::zero_rank.
|
inline |
ctor for a TensorTrain, with core tensors explicitly given
the core tensors must have the shape (k1,r1) (r1,k2,r2) .. (r2,k3)
[in] | core | vector 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.
|
inline |
copy constructor, shallow
References madness::BaseTensor::dims(), madness::BaseTensor::ndim(), and madness::BaseTensor::set_dims_and_size().
|
inline |
decompose the input tensor into a TT representation
[in] | t | tensor in full rank |
[in] | eps | the precision threshold |
[in] | dims | the 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().
|
inline |
compute the Hadamard product of two TensorTrains
see Sec. 4.2 of the TT paper
References madness::_(), a1, madness::TensorTrain< T >::copy, madness::TensorTrain< T >::core, madness::BaseTensor::dim(), madness::BaseTensor::dims(), madness::TensorTrain< T >::fusedim(), k, MADNESS_ASSERT, MADNESS_EXCEPTION, madness::BaseTensor::ndim(), madness::TensorTrain< T >::outer, madness::TensorTrain< T >::ranks(), madness::Tensor< T >::reshape(), madness::TensorTrain< T >::size(), and madness::TensorTrain< T >::zero_rank.
|
inline |
merge two dimensions into one
merge dimension i and i+1 into new dimension i
[in] | i | the 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().
|
inline |
Inplace generalized saxpy with slices and without alpha.
return this = this(s1) + other(s2) * beta
References madness::_(), beta, madness::TensorTrain< T >::core, d(), madness::BaseTensor::dim(), k, MADNESS_ASSERT, MADNESS_EXCEPTION, madness::BaseTensor::ndim(), madness::TensorTrain< T >::verify(), and madness::TensorTrain< T >::zero_rank.
|
inline |
Inplace generalized saxpy ... this = this*alpha + other*beta.
References madness::_(), alpha, beta, madness::TensorTrain< T >::core, d(), k, MADNESS_ASSERT, MADNESS_EXCEPTION, madness::BaseTensor::ndim(), madness::TensorTrain< T >::scale(), madness::TensorTrain< T >::verify(), and madness::TensorTrain< T >::zero_rank.
Referenced by madness::TensorTrain< T >::operator+=(), and madness::TensorTrain< T >::operator-=().
|
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().
|
inline |
const reference to the internal core
References madness::TensorTrain< T >::core, MADNESS_ASSERT, and madness::BaseTensor::ndim().
|
inline |
check if this is an operator (r,k',k,r)
References madness::TensorTrain< T >::is_tensor().
Referenced by madness::TensorTrain< T >::make_operator().
|
inline |
check if this is a tensor (r,k,r)
References madness::TensorTrain< T >::core, and madness::BaseTensor::ndim().
Referenced by madness::TensorTrain< T >::is_operator(), and madness::TensorTrain< T >::make_tensor().
|
inline |
if rank is zero
References madness::TensorTrain< T >::zero_rank.
|
inline |
convert this into an operator representation (r,k',k,r)
References madness::TensorTrain< T >::core, madness::TensorTrain< T >::is_operator(), k0, k2, MADNESS_ASSERT, and madness::BaseTensor::ndim().
Referenced by madness::SeparatedConvolution< Q, NDIM >::make_tt_representation().
|
inline |
convert this into a tensor representation (r,k,r)
References madness::TensorTrain< T >::core, madness::TensorTrain< T >::fusedim(), madness::TensorTrain< T >::is_tensor(), and madness::BaseTensor::ndim().
Referenced by madness::SeparatedConvolution< Q, NDIM >::make_tt_representation().
|
inline |
returns the Frobenius norm
References std::abs(), and madness::TensorTrain< T >::trace().
|
inline |
Type conversion makes a deep copy.
References c, madness::TensorTrain< T >::core, madness::BaseTensor::dims(), and madness::TensorTrain< T >::zero_rank.
|
inline |
return this multiplied by a scalar
References madness::TensorTrain< T >::copy, and madness::TensorTrain< T >::scale().
|
inline |
inplace addition of two Tensortrains; will increase ranks of this
inefficient if many additions are performed, since it requires many calls of new.
[in] | rhs | a TensorTrain to be added |
References madness::TensorTrain< T >::gaxpy().
|
inline |
inplace subtraction of two Tensortrains; will increase ranks of this
inefficient if many subtractions are performed, since it requires many calls of new.
[in] | rhs | a TensorTrain to be added |
References madness::TensorTrain< T >::gaxpy().
|
inline |
assign a number to this tensor
References madness::TensorTrain< T >::core, madness::BaseTensor::dim(), madness::BaseTensor::ndim(), and madness::TensorTrain< T >::zero_rank.
|
inline |
assigment operator
References madness::TensorTrain< T >::core, madness::BaseTensor::dims(), madness::BaseTensor::ndim(), madness::BaseTensor::set_dims_and_size(), and madness::TensorTrain< T >::zero_rank.
|
inline |
Returns a pointer to the internal data.
[in] | ivec | index of core vector to which the return values points |
References madness::TensorTrain< T >::core.
|
inline |
Returns a pointer to the internal data.
[in] | ivec | index of core vector to which the return values points |
References madness::TensorTrain< T >::core.
|
inline |
return the TT ranks
References madness::TensorTrain< T >::core, madness::BaseTensor::dim(), madness::BaseTensor::ndim(), and madness::TensorTrain< T >::zero_rank.
Referenced by madness::apply(), madness::TensorTrain< T >::emul(), madness::TensorTrain< T >::splitdim(), test_stuff(), and madness::TensorTrain< T >::truncate().
|
inline |
return the TT ranks for dimension i (to i+1)
References madness::TensorTrain< T >::core, MADNESS_EXCEPTION, madness::BaseTensor::ndim(), madness::print(), and madness::TensorTrain< T >::zero_rank.
|
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().
|
inline |
reconstruct this to a full representation
[in] | flat | return this in flat 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().
|
inline |
scale this by a number
[in] | fac | the factor to multiply |
References madness::TensorTrain< T >::core, and madness::TensorTrain< T >::zero_rank.
Referenced by madness::TensorTrain< T >::gaxpy(), and madness::TensorTrain< T >::operator*().
|
inline |
|
inline |
References madness::BaseTensor::dims(), and madness::BaseTensor::set_dims_and_size().
Referenced by madness::TensorTrain< T >::TensorTrain().
|
inline |
return the number of coefficients in all core tensors
References madness::TensorTrain< T >::core, and madness::TensorTrain< T >::zero_rank.
Referenced by madness::TensorTrain< T >::emul(), madness::TensorTrain< T >::real_size(), madness::TensorTrain< T >::reconstruct(), madness::TensorTrain< T >::truncate(), and madness::TensorTrain< T >::verify().
|
inline |
Returns new view/tensor splitting dimension i
as dimi0*dimi1
to produce conforming d+1 dimension tensor
[in] | idim | the dimension to be split |
[in] | k1 | new first dimension of idim |
[in] | k2 | new second dimension of idim |
[in] | eps | threshold for SVD (choose negative to keep all terms) |
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.
|
inline |
References MADNESS_EXCEPTION.
Referenced by madness::TensorTrain< T >::normf().
|
inline |
Return the trace of two tensors, no complex conjugate involved.
References d(), madness::BaseTensor::dim(), madness::Tensor< T >::fusedim(), madness::inner_result(), MADNESS_ASSERT, MADNESS_EXCEPTION, max, madness::BaseTensor::ndim(), Q(), madness::Tensor< T >::reshape(), T(), madness::TENSOR_RESULT_TYPE(), and madness::Tensor< T >::trace().
|
inline |
recompress and truncate this TT representation
this in recompressed TT form with optimal rank
[in] | eps | the truncation threshold |
References MADNESS_EXCEPTION.
Referenced by madness::SeparatedConvolution< Q, NDIM >::make_tt_representation().
|
inline |
recompress and truncate this TT representation
this in recompressed TT form with optimal rank
[in] | eps | the truncation threshold |
References madness::_(), madness::copy(), madness::TensorTrain< T >::core, d(), madness::BaseTensor::dim(), madness::Tensor< T >::flat(), madness::inner(), k, madness::kmax, L, madness::lq_result(), m, MADNESS_EXCEPTION, max, madness::SRConf< T >::max_sigma(), madness::BaseTensor::ndim(), madness::TensorTrain< T >::ranks(), madness::Tensor< T >::reshape(), madness::TensorTrain< T >::size(), madness::svd_result(), TENSOR_MAXDIM, madness::TensorTrain< T >::verify(), madness::TensorTrain< T >::zero_me(), and madness::TensorTrain< T >::zero_rank.
|
inline |
construct a two-mode representation (aka unnormalized SVD)
[out] | U | The left singular vectors, ({i},rank) |
[out] | VT | The right singular vectors, (rank,{i}) |
[out] | s | Vector 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.
|
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().
|
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().
|
inline |
turn this into an empty tensor with all cores properly shaped
References madness::BaseTensor::dim().
Referenced by madness::TensorTrain< T >::splitdim().
|
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().
|
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
[in] | other | tensor to be sliced |
[in] | s | vector of slices |
|
friend |
|
friend |
Referenced by madness::TensorTrain< T >::emul().
|
friend |
|
friend |
|
private |
holding the core tensors of a tensor train the tensors have the shape (k,r0) (r0,k,r1) (r1,k,r2) .. (rn-1,k)
Referenced by madness::TensorTrain< T >::TensorTrain(), madness::TensorTrain< T >::decompose(), madness::TensorTrain< T >::emul(), madness::TensorTrain< T >::fusedim(), madness::TensorTrain< T >::gaxpy(), madness::general_transform(), madness::TensorTrain< T >::get_core(), madness::TensorTrain< T >::is_tensor(), madness::TensorTrain< T >::make_operator(), madness::TensorTrain< T >::make_tensor(), madness::TensorTrain< T >::operator TensorTrain< Q >(), madness::TensorTrain< T >::operator=(), madness::outer(), madness::TensorTrain< T >::ptr(), madness::TensorTrain< T >::ranks(), madness::TensorTrain< T >::real_size(), madness::TensorTrain< T >::reconstruct(), madness::TensorTrain< T >::scale(), madness::TensorTrain< T >::serialize(), madness::TensorTrain< T >::size(), madness::TensorTrain< T >::splitdim(), madness::transform(), madness::transform_dir(), madness::TensorTrain< T >::truncate(), madness::TensorTrain< T >::two_mode_representation(), and madness::TensorTrain< T >::verify().
|
private |
true if rank is zero
Referenced by madness::TensorTrain< T >::TensorTrain(), madness::TensorTrain< T >::decompose(), madness::TensorTrain< T >::emul(), madness::TensorTrain< T >::fusedim(), madness::TensorTrain< T >::gaxpy(), madness::general_transform(), madness::TensorTrain< T >::is_zero_rank(), madness::TensorTrain< T >::operator TensorTrain< Q >(), madness::TensorTrain< T >::operator=(), madness::outer(), madness::TensorTrain< T >::ranks(), madness::TensorTrain< T >::reconstruct(), madness::TensorTrain< T >::scale(), madness::TensorTrain< T >::serialize(), madness::TensorTrain< T >::size(), madness::TensorTrain< T >::splitdim(), madness::transform(), madness::transform_dir(), madness::TensorTrain< T >::truncate(), and madness::TensorTrain< T >::two_mode_representation().