xerus
a general purpose tensor library
xerus Namespace Reference

The main namespace of xerus. More...

Namespaces

 blasWrapper
 In this namespace the minimal wrappers for the BLAS and LAPACK functions are collected.
 
 internal
 Namespace for function and classes designated only for internal use.
 
 misc
 Collection of classes and functions that provide elementary functionality that is not special to xerus as a tensor library.
 

Classes

class  ADFVariant
 Wrapper class for all ADF variants. More...
 
class  ALSVariant
 Wrapper class for all ALS variants (dmrg etc.) More...
 
class  CQ
 Helper class to allow an intuitive syntax for an rank revealing orthogonal factorisation. More...
 
class  GeometricCGVariant
 Wrapper class for all geometric (ie. Riemannian) CG variants. More...
 
struct  HOSVDRetraction
 retraction that performs a HOSVD to project back onto the Manifold More...
 
class  Index
 Class used to represent indices that can be used to write tensor calculations in index notation. More...
 
class  InternalSolver
 
class  MeasurmentComparator
 
class  PerformanceData
 Storage class for the performance data collected during an algorithm (typically iteration count, time and residual) More...
 
class  QC
 Helper class to allow an intuitive syntax for an rank revealing orthogonal factorisation. More...
 
class  QR
 Helper class to allow an intuitive syntax for QR factorisations. More...
 
class  RankOneMeasurementSet
 
class  RQ
 Helper class to allow an intuitive syntax for RQ factorisations. More...
 
class  SinglePointMeasurementSet
 Class used to represent a single point measurments. More...
 
class  SteepestDescentVariant
 Wrapper class for all steepest descent variants. More...
 
class  SVD
 Helper class to allow an intuitive syntax for SVD factorisations. More...
 
class  Tensor
 Class that handles simple (non-decomposed) tensors in a dense or sparse representation. More...
 
class  TensorFactorisation
 Abstract super class for all tensor factorisations. More...
 
class  TensorNetwork
 Very general class used to represent arbitary tensor networks. More...
 
class  TTNetwork
 Specialized TensorNetwork class used to represent TTTensor and TToperators. More...
 
class  TTTangentVector
 class to compactly represent tangent vectors of the manifold of constant TT-rank More...
 
class  UQMeasurementSet
 

Typedefs

typedef uint8_t byte
 unsigned int type of exactly 8 bit More...
 
typedef int16_t int16
 
typedef int32_t int32
 
typedef int64_t int64
 
typedef int8_t int8
 
typedef TTNetwork< true > TTOperator
 
using TTRetractionI = std::function< void(TTTensor &, const TTTangentVector &)>
 
using TTRetractionII = std::function< void(TTTensor &, const TTTensor &)>
 
typedef TTNetwork< false > TTTensor
 
using TTVectorTransport = std::function< void(const TTTensor &, TTTangentVector &)>
 
typedef unsigned int uint
 
typedef uint16_t uint16
 
typedef uint32_t uint32
 
typedef uint64_t uint64
 
typedef uint8_t uint8
 
typedef unsigned long ulong
 
typedef unsigned short ushort
 
typedef double value_t
 The type of values to be used by xerus. More...
 

Functions

const ADFVariant ADF (0, 1e-8, 0.999)
 
const ALSVariant ALS (1, 0, ALSVariant::lapack_solver, false)
 
const ALSVariant ALS_SPD (1, 0, ALSVariant::lapack_solver, true)
 
void ALSRetractionI (TTTensor &_U, const TTTangentVector &_change)
 retraction that performs an ALS half-sweep to project back onto the Manifold. Automatically retains the ranks of _U More...
 
void ALSRetractionII (TTTensor &_U, const TTTensor &_change)
 retraction that performs an ALS half-sweep to project back onto the Manifold. Automatically retains the ranks of _U More...
 
bool approx_entrywise_equal (const xerus::Tensor &_a, const xerus::Tensor &_b, const xerus::value_t _eps=EPSILON)
 Checks whether two Tensors are approximately entrywise equal. More...
 
bool approx_entrywise_equal (const xerus::Tensor &_tensor, const std::vector< value_t > &_values, const xerus::value_t _eps=EPSILON)
 Compares the given Tensor entriewise to the given data. More...
 
bool approx_equal (const TensorNetwork &_a, const TensorNetwork &_b, const value_t _eps=EPSILON)
 Checks whether two TensorNetworks are approximately equal. More...
 
bool approx_equal (const TensorNetwork &_a, const Tensor &_b, const value_t _eps=EPSILON)
 Convinience wrapper, casts the the given TensorNetwork _a to Tensor and calls the Tensor function. More...
 
bool approx_equal (const Tensor &_a, const TensorNetwork &_b, const value_t _eps=EPSILON)
 Convinience wrapper, casts the the given TensorNetwork _b to Tensor and calls the Tensor function. More...
 
bool approx_equal (const xerus::Tensor &_a, const xerus::Tensor &_b, const xerus::value_t _eps=EPSILON)
 Checks whether two tensors are approximately equal. More...
 
const ALSVariant ASD (1, 0, ALSVariant::ASD_solver, false)
 
const ALSVariant ASD_SPD (1, 0, ALSVariant::ASD_solver, true)
 
template<class MeasurmentSet >
void calc_measurment_norm (double *_norms, const MeasurmentSet &_measurments)
 
template<>
void calc_measurment_norm< RankOneMeasurementSet > (double *_norms, const RankOneMeasurementSet &_measurments)
 
template<>
void calc_measurment_norm< SinglePointMeasurementSet > (double *_norms, const SinglePointMeasurementSet &_measurments)
 
void calculate_cq (Tensor &_C, Tensor &_Q, Tensor _input, const size_t _splitPos)
 Low-Level CQ calculation of a given Tensor _input = _C _Q. More...
 
XERUS_force_inline std::tuple< size_t, size_t, size_t > calculate_factorization_sizes (const Tensor &_input, const size_t _splitPos)
 
void calculate_qc (Tensor &_Q, Tensor &_C, Tensor _input, const size_t _splitPos)
 Low-Level QC calculation of a given Tensor _input = _Q _C. More...
 
void calculate_qr (Tensor &_Q, Tensor &_R, Tensor _input, const size_t _splitPos)
 Low-Level QR calculation of a given Tensor _input = _Q _R. More...
 
void calculate_rq (Tensor &_R, Tensor &_Q, Tensor _input, const size_t _splitPos)
 Low-Level RQ calculation of a given Tensor _input = _R _Q. More...
 
void calculate_svd (Tensor &_U, Tensor &_S, Tensor &_Vt, Tensor _input, const size_t _splitPos, const size_t _maxRank, const value_t _eps)
 Low-Level SVD calculation of a given Tensor _input = _U _S _Vt. More...
 
void contract (Tensor &_result, const Tensor &_lhs, const bool _lhsTrans, const Tensor &_rhs, const bool _rhsTrans, const size_t _numModes)
 Low-level contraction between Tensors. More...
 
Tensor contract (const Tensor &_lhs, const bool _lhsTrans, const Tensor &_rhs, const bool _rhsTrans, const size_t _numModes)
 
XERUS_force_inline void contract (Tensor &_result, const Tensor &_lhs, const Tensor &_rhs, const size_t _numModes)
 Low-level contraction between Tensors. More...
 
XERUS_force_inline Tensor contract (const Tensor &_lhs, const Tensor &_rhs, const size_t _numModes)
 
void decomposition_als (TTTensor &_x, const Tensor &_b, const double _eps=EPSILON, const size_t _maxIterations=1000)
 
const ALSVariant DMRG (2, 0, ALSVariant::lapack_solver, false)
 
const ALSVariant DMRG_SPD (2, 0, ALSVariant::lapack_solver, true)
 
template<bool isOperator>
TTNetwork< isOperator > dyadic_product (const TTNetwork< isOperator > &_lhs, const TTNetwork< isOperator > &_rhs)
 Computes the dyadic product of _lhs and _rhs. More...
 
template<bool isOperator>
TTNetwork< isOperator > dyadic_product (const std::vector< TTNetwork< isOperator >> &_tensors)
 Computes the dyadic product of all given TTNetworks. More...
 
template TTNetwork< true > dyadic_product (const TTNetwork< true > &_lhs, const TTNetwork< true > &_rhs)
 
template TTNetwork< false > dyadic_product (const TTNetwork< false > &_lhs, const TTNetwork< false > &_rhs)
 
template TTNetwork< true > dyadic_product (const std::vector< TTNetwork< true >> &_tensors)
 
template TTNetwork< false > dyadic_product (const std::vector< TTNetwork< false >> &_tensors)
 
template<bool isOperator>
TTNetwork< isOperator > entrywise_product (const TTNetwork< isOperator > &_A, const TTNetwork< isOperator > &_B)
 Calculates the componentwise product of two tensors given in the TT format. More...
 
Tensor entrywise_product (const Tensor &_A, const Tensor &_B)
 calculates the entrywise product of two Tensors More...
 
template TTNetwork< false > entrywise_product (const TTNetwork< false > &_A, const TTNetwork< false > &_B)
 
template TTNetwork< true > entrywise_product (const TTNetwork< true > &_A, const TTNetwork< true > &_B)
 
template<bool isOperator>
size_t find_largest_entry (const TTNetwork< isOperator > &_T, double _accuracy, value_t _lowerBound=0.0)
 Finds the position of the approximately largest entry. More...
 
template size_t find_largest_entry (const TTNetwork< true > &, double, value_t)
 
template size_t find_largest_entry (const TTNetwork< false > &, double, value_t)
 
static XERUS_force_inline value_t frob_norm (const TensorNetwork &_network)
 Calculates the frobenious norm of the given TensorNetwork. More...
 
static XERUS_force_inline value_t frob_norm (const Tensor &_tensor)
 Calculates the frobenius norm of the given tensor. More...
 
const GeometricCGVariant GeometricCG (0, 1e-8, false, SubmanifoldRetractionI, ProjectiveVectorTransport)
 
template<bool isOperator>
std::vector< std::vector< std::tuple< size_t, size_t, value_t > > > get_grouped_entries (const Tensor &_component)
 
std::vector< size_t > get_step_sizes (const Tensor::DimensionTuple &_dimensions)
 
void HOSVDRetractionI (TTTensor &_U, const TTTangentVector &_change)
 retraction that performs a HOSVD to project back onto the Manifold More...
 
void HOSVDRetractionII (TTTensor &_U, const TTTensor &_change)
 
double IHT (TTTensor &_x, const SinglePointMeasurementSet &_measurments, PerformanceData &_perfData=NoPerfData)
 
void increase_indices (const size_t _i, value_t *&_outPosition, const std::vector< size_t > &_stepSizes, const std::vector< size_t > &_baseDimensions)
 
void line_search (TTTensor &_x, value_t &_alpha, const TTTangentVector &_direction, value_t _derivative, value_t &_residual, std::function< void(TTTensor &, const TTTangentVector &)> _retraction, std::function< value_t()> _calculate_residual, value_t _changeInAlpha=0.5)
 
void matrix_matrix_product (double *const _C, const size_t _leftDim, const size_t _rightDim, const double _alpha, const std::map< size_t, double > &_A, const bool _transposeA, const size_t _midDim, const double *const _B, const bool _transposeB)
 
void matrix_matrix_product (double *const _C, const size_t _leftDim, const size_t _rightDim, const double _alpha, const double *const _A, const bool _transposeA, const size_t _midDim, const std::map< size_t, double > &_B, const bool _transposeB)
 
void matrix_matrix_product (std::map< size_t, double > &_C, const size_t _leftDim, const size_t _rightDim, const double _alpha, const std::map< size_t, double > &_A, const bool _transposeA, const size_t _midDim, const double *const _B, const bool _transposeB)
 
void matrix_matrix_product (std::map< size_t, double > &_C, const size_t _leftDim, const size_t _rightDim, const double _alpha, const double *const _A, const bool _transposeA, const size_t _midDim, const std::map< size_t, double > &_B, const bool _transposeB)
 
void matrix_matrix_product (double *const _C, const size_t _leftDim, const size_t _rightDim, const double _alpha, const std::map< size_t, double > &_A, const bool _transposeA, const size_t _midDim, const double *const _B)
 
void matrix_matrix_product (std::map< size_t, double > &_C, const size_t, const size_t _rightDim, const double _alpha, const std::map< size_t, double > &_A, const size_t _midDim, const double *const _B)
 
PerformanceData NoPerfData (false, false)
 
static XERUS_force_inline value_t one_norm (const Tensor &_tensor)
 Calculates the 1-norm of the given tensor. More...
 
bool operator!= (const Index &_a, const Index &_b)
 Two Indices are equal if their valueId coincides. Fixed indices are never equal. More...
 
TTTangentVector operator* (value_t _alpha, const TTTangentVector &_rhs)
 
TensorNetwork operator* (TensorNetwork &_lhs, value_t _factor)
 
TensorNetwork operator* (value_t _factor, TensorNetwork &_rhs)
 
template<bool isOperator>
TTNetwork< isOperator > operator* (TTNetwork< isOperator > _network, const value_t _factor)
 Calculates the entrywise multiplication of the given TTNetwork _network with a constant _factor. More...
 
template<bool isOperator>
TTNetwork< isOperator > operator* (const value_t _factor, TTNetwork< isOperator > _network)
 Calculates the entrywise multiplication of the given TTNetwork _network with a constant _factor. More...
 
Tensor operator* (const value_t _factor, Tensor _tensor)
 Calculates the entrywise multiplication of the Tensor _tensor with the constant _factor. More...
 
Tensor operator* (Tensor _tensor, const value_t _factor)
 Calculates the entrywise multiplication of the Tensor _tensor with the constant _factor. More...
 
template TTNetwork< false > operator* (TTNetwork< false > _network, const value_t _factor)
 
template TTNetwork< true > operator* (TTNetwork< true > _network, const value_t _factor)
 
template TTNetwork< false > operator* (const value_t _factor, TTNetwork< false > _network)
 
template TTNetwork< true > operator* (const value_t _factor, TTNetwork< true > _network)
 
template<bool isOperator>
TTNetwork< isOperator > operator+ (TTNetwork< isOperator > _lhs, const TTNetwork< isOperator > &_rhs)
 Calculates the entrywise sum of the given TTNetworks _lhs and _rhs. More...
 
Tensor operator+ (Tensor _lhs, const Tensor &_rhs)
 Calculates the entrywise sum of _lhs and _rhs. More...
 
template TTNetwork< false > operator+ (TTNetwork< false > _lhs, const TTNetwork< false > &_rhs)
 
template TTNetwork< true > operator+ (TTNetwork< true > _lhs, const TTNetwork< true > &_rhs)
 
template<bool isOperator>
TTNetwork< isOperator > operator- (TTNetwork< isOperator > _lhs, const TTNetwork< isOperator > &_rhs)
 Calculates the entrywise difference of the given TTNetworks _lhs and _rhs. More...
 
Tensor operator- (Tensor _lhs, const Tensor &_rhs)
 Calculates the entrywise difference between _lhs and _rhs. More...
 
template TTNetwork< false > operator- (TTNetwork< false > _lhs, const TTNetwork< false > &_rhs)
 
template TTNetwork< true > operator- (TTNetwork< true > _lhs, const TTNetwork< true > &_rhs)
 
internal::IndexedTensorMoveable< Tensoroperator/ (internal::IndexedTensorReadOnly< Tensor > &&_b, internal::IndexedTensorReadOnly< Tensor > &&_A)
 
TensorNetwork operator/ (TensorNetwork &_lhs, value_t _factor)
 
template<bool isOperator>
TTNetwork< isOperator > operator/ (TTNetwork< isOperator > _network, const value_t _divisor)
 Calculates the entrywise divison of this TTNetwork by a constant _divisor. More...
 
Tensor operator/ (Tensor _tensor, const value_t _divisor)
 Calculates the entrywise divison of the Tensor _tensor with the constant _divisor. More...
 
template TTNetwork< false > operator/ (TTNetwork< false > _network, const value_t _divisor)
 
template TTNetwork< true > operator/ (TTNetwork< true > _network, const value_t _divisor)
 
bool operator< (const Index &_a, const Index &_b)
 The Comparision operator is needed for indices to be orderable in std::set, the valueId is used. More...
 
std::ostream & operator<< (std::ostream &_out, const xerus::Index &_idx)
 Allows to pretty print indices, giving the valueId and span. More...
 
std::ostream & operator<< (std::ostream &_out, const TensorNetwork::Link &_rhs)
 
std::ostream & operator<< (std::ostream &_out, const Tensor &_tensor)
 Prints the Tensor to the given outStream. More...
 
bool operator== (const Index &_a, const Index &_b)
 Two Indices are equal if their valueId coincides. Fixed indices are never equal. More...
 
template<bool isOperator>
void perform_component_product (Tensor &_newComponent, const Tensor &_componentA, const Tensor &_componentB)
 
template<class MeasurmentSet >
size_t position_or_zero (const MeasurmentSet &_measurments, const size_t _meas, const size_t _corePosition)
 
template<>
size_t position_or_zero< RankOneMeasurementSet > (const RankOneMeasurementSet &, const size_t, const size_t)
 
template<>
size_t position_or_zero< SinglePointMeasurementSet > (const SinglePointMeasurementSet &_measurments, const size_t _meas, const size_t _corePosition)
 
XERUS_force_inline void prepare_factorization_output (Tensor &_lhs, Tensor &_rhs, const Tensor &_input, const size_t _splitPos, const size_t _rank, Tensor::Representation _rep)
 
std::unique_ptr< Tensorprepare_split (size_t &_lhsSize, size_t &_rhsSize, size_t &_rank, size_t &_splitPos, std::vector< Index > &_lhsPreliminaryIndices, std::vector< Index > &_rhsPreliminaryIndices, internal::IndexedTensorReadOnly< Tensor > &&_base, internal::IndexedTensor< Tensor > &&_lhs, internal::IndexedTensor< Tensor > &&_rhs)
 
void ProjectiveVectorTransport (const TTTensor &_newBase, TTTangentVector &_tangentVector)
 simple vector transport by projecting onto the new tangent plane More...
 
void pseudo_inverse (Tensor &_inverse, const Tensor &_input, const size_t _splitPos)
 Low-Level calculation of the pseudo inverse of a given Tensor. More...
 
Tensor pseudo_inverse (const Tensor &_input, const size_t _splitPos)
 
Tensor randVar_to_position (const double _v, const size_t _polyDegree)
 
void reshuffle (Tensor &_out, const Tensor &_base, const std::vector< size_t > &_shuffle)
 : Performs a simple reshuffle. Much less powerfull then a full evaluate, but more efficient. More...
 
Tensor reshuffle (const Tensor &_base, const std::vector< size_t > &_shuffle)
 
XERUS_force_inline void set_factorization_output (Tensor &_lhs, std::unique_ptr< value_t[]> &&_lhsData, Tensor &_rhs, std::unique_ptr< value_t[]> &&_rhsData, const Tensor &_input, const size_t _splitPos, const size_t _rank)
 
XERUS_force_inline void set_factorization_output (Tensor &_lhs, std::map< size_t, value_t > &&_lhsData, Tensor &_rhs, std::map< size_t, value_t > &&_rhsData, const Tensor &_input, const size_t _splitPos, const size_t _rank)
 
void solve (internal::IndexedTensorWritable< Tensor > &&_x, internal::IndexedTensorReadOnly< Tensor > &&_A, internal::IndexedTensorReadOnly< Tensor > &&_b)
 
void solve (Tensor &_X, const Tensor &_A, const Tensor &_B, size_t _extraDegree=0)
 Solves the equation Ax = b for x. If the solution is not unique, the output need not be the minimal norm solution. More...
 
void solve_least_squares (Tensor &_X, const Tensor &_A, const Tensor &_B, const size_t _extraDegree=0)
 Solves the least squares problem ||_A _X - _B||_F. More...
 
const SteepestDescentVariant SteepestDescent (0, 1e-8, false, SubmanifoldRetractionII)
 
void SubmanifoldRetractionI (TTTensor &_U, const TTTangentVector &_change)
 retraction that performs componentwise addition of \( U_i \) and \( W_i \) where \( W_i \) is the i-th component of the riemannian tangential vector representation More...
 
void SubmanifoldRetractionII (TTTensor &_U, const TTTensor &_change)
 retraction that performs componentwise addition of \( U_i \) and \( W_i \) where \( W_i \) is the i-th component of the riemannian tangential vector representation More...
 
XERUS_force_inline void transpose (double *const __restrict _out, const double *const __restrict _in, const size_t _leftDim, const size_t _rightDim)
 
XERUS_force_inline std::unique_ptr< double[]> transpose (const double *const _A, const size_t _leftDim, const size_t _rightDim)
 
XERUS_force_inline void transpose (std::map< size_t, double > &__restrict _out, const std::map< size_t, double > &__restrict _in, const size_t _leftDim, const size_t _rightDim)
 
XERUS_force_inline std::map< size_t, double > transpose (const std::map< size_t, double > &_A, const size_t _leftDim, const size_t _rightDim)
 
template<bool isOperator>
void transpose_if_operator (TTNetwork< isOperator > &_ttNetwork)
 
template<>
void transpose_if_operator< false > (TTNetwork< false > &)
 
template<>
void transpose_if_operator< true > (TTNetwork< true > &_ttNetwork)
 
void uq_adf (TTTensor &_x, const std::vector< std::vector< double >> &_randomVariables, const std::vector< Tensor > &_solutions)
 
TTTensor uq_adf (const UQMeasurementSet &_measurments, const TTTensor &_guess)
 
Tensor uq_avg (const TTTensor &_x, const size_t _N, const size_t _numSpecial)
 
std::pair< std::vector< std::vector< double > >, std::vector< Tensor > > uq_mc (const TTTensor &_x, const size_t _N, const size_t _numSpecial)
 

Variables

const ADFVariant ADF
 Default variant of the ADF algorithm. More...
 
const ALSVariant ALS
 default variant of the single-site ALS algorithm for non-symmetric operators using the lapack solver More...
 
const ALSVariant ALS_SPD
 default variant of the single-site ALS algorithm for symmetric positive-definite operators using the lapack solver More...
 
const ALSVariant ASD
 default variant of the alternating steepest descent for non-symmetric operators More...
 
const ALSVariant ASD_SPD
 default variant of the alternating steepest descent for symmetric positive-definite operators More...
 
const ALSVariant DMRG
 default variant of the two-site DMRG algorithm for non-symmetric operators using the lapack solver More...
 
const ALSVariant DMRG_SPD
 default variant of the two-site DMRG algorithm for symmetric positive-definite operators using the lapack solver More...
 
constexpr const value_t EPSILON = 8*std::numeric_limits<value_t>::epsilon()
 The default epsilon value in xerus. More...
 
const GeometricCGVariant GeometricCG
 default variant of the steepest descent algorithm using the lapack solver More...
 
PerformanceData NoPerfData
 
const SteepestDescentVariant SteepestDescent
 default variant of the steepest descent algorithm using the lapack solver More...
 
const int VERSION_COMMIT = XERUS_VERSION_COMMIT
 
const int VERSION_MAJOR = XERUS_VERSION_MAJOR
 
const int VERSION_MINOR = XERUS_VERSION_MINOR
 
const int VERSION_REVISION = XERUS_VERSION_REVISION
 

Detailed Description

The main namespace of xerus.

With very few exceptions all classes, functions and variables declared by the xerus are placed in this namespace.

Typedef Documentation

◆ byte

typedef uint8_t xerus::byte

unsigned int type of exactly 8 bit

Definition at line 49 of file standard.h.

◆ int16

typedef int16_t xerus::int16

Definition at line 56 of file standard.h.

◆ int32

typedef int32_t xerus::int32

Definition at line 57 of file standard.h.

◆ int64

typedef int64_t xerus::int64

Definition at line 58 of file standard.h.

◆ int8

typedef int8_t xerus::int8

Definition at line 55 of file standard.h.

◆ TTOperator

typedef TTNetwork< true > xerus::TTOperator

Definition at line 40 of file performanceData.h.

◆ TTRetractionI

using xerus::TTRetractionI = typedef std::function<void(TTTensor &, const TTTangentVector &)>

Definition at line 70 of file retractions.h.

◆ TTRetractionII

using xerus::TTRetractionII = typedef std::function<void(TTTensor &, const TTTensor &)>

Definition at line 71 of file retractions.h.

◆ TTTensor

typedef TTNetwork< false > xerus::TTTensor

Definition at line 37 of file performanceData.h.

◆ TTVectorTransport

using xerus::TTVectorTransport = typedef std::function<void(const TTTensor &, TTTangentVector &)>

Definition at line 72 of file retractions.h.

◆ uint

typedef unsigned int xerus::uint

Definition at line 51 of file standard.h.

◆ uint16

typedef uint16_t xerus::uint16

Definition at line 61 of file standard.h.

◆ uint32

typedef uint32_t xerus::uint32

Definition at line 62 of file standard.h.

◆ uint64

typedef uint64_t xerus::uint64

Definition at line 63 of file standard.h.

◆ uint8

typedef uint8_t xerus::uint8

Definition at line 60 of file standard.h.

◆ ulong

typedef unsigned long xerus::ulong

Definition at line 52 of file standard.h.

◆ ushort

typedef unsigned short xerus::ushort

Definition at line 50 of file standard.h.

◆ value_t

typedef double xerus::value_t

The type of values to be used by xerus.

In future versions this should be allowed to be float, double, or complex. In the current version however this is fixed to double.

Definition at line 43 of file basic.h.

Function Documentation

◆ ADF()

const ADFVariant xerus::ADF ( ,
1e-  8,
0.  999 
)

◆ ALS()

const ALSVariant xerus::ALS ( ,
,
ALSVariant::lapack_solver  ,
false   
)

◆ ALS_SPD()

const ALSVariant xerus::ALS_SPD ( ,
,
ALSVariant::lapack_solver  ,
true   
)

◆ ALSRetractionI()

void xerus::ALSRetractionI ( TTTensor _U,
const TTTangentVector _change 
)

retraction that performs an ALS half-sweep to project back onto the Manifold. Automatically retains the ranks of _U

Definition at line 69 of file retractions.cpp.

◆ ALSRetractionII()

void xerus::ALSRetractionII ( TTTensor _U,
const TTTensor _change 
)

retraction that performs an ALS half-sweep to project back onto the Manifold. Automatically retains the ranks of _U

Definition at line 62 of file retractions.cpp.

◆ approx_entrywise_equal() [1/2]

bool xerus::approx_entrywise_equal ( const xerus::Tensor _a,
const xerus::Tensor _b,
const xerus::value_t  _eps = EPSILON 
)

Checks whether two Tensors are approximately entrywise equal.

Check whether |_a[i] - _b[i] |/(|_a[i] | + |_b[i] | < _eps is fullfilled for all i.

Parameters
_athe first test candidate.
_bthe second test candidate
_epsthe maximal relative difference between the entries of _a and _b.
Returns
TRUE if _a and _b are determined to be approximately equal, FALSE otherwise.

Definition at line 1745 of file tensor.cpp.

◆ approx_entrywise_equal() [2/2]

bool xerus::approx_entrywise_equal ( const xerus::Tensor _tensor,
const std::vector< value_t > &  _values,
const xerus::value_t  _eps = EPSILON 
)

Compares the given Tensor entriewise to the given data.

Parameters
_tensorthe tensor.
_valuesthe data to compare to, must be of the same size as the Tensor has entries (i.e. size).
_eps(optional) epsilon used to determine whether an entry is equal to a data point.
Returns
TRUE if approx_equal() is true for all entries/data points, FALSE otherwise.

Definition at line 1765 of file tensor.cpp.

◆ approx_equal() [1/4]

bool xerus::approx_equal ( const TensorNetwork _a,
const TensorNetwork _b,
const value_t  _eps = EPSILON 
)

Checks whether two TensorNetworks are approximately equal.

Check whether ||_a - _b ||/(||a ||/2 + ||_b ||/2) < _eps, i.e. whether the relative difference in the frobenius norm is sufficently small.

Parameters
_athe first test candidate.
_bthe second test candidate
_epsthe maximal relative difference between _a and _b.
Returns
TRUE if _a and _b are determined to be approximately equal, FALSE otherwise.

Definition at line 1411 of file tensorNetwork.cpp.

◆ approx_equal() [2/4]

bool xerus::approx_equal ( const TensorNetwork _a,
const Tensor _b,
const value_t  _eps = EPSILON 
)

Convinience wrapper, casts the the given TensorNetwork _a to Tensor and calls the Tensor function.

Definition at line 1418 of file tensorNetwork.cpp.

◆ approx_equal() [3/4]

bool xerus::approx_equal ( const Tensor _a,
const TensorNetwork _b,
const value_t  _eps = EPSILON 
)

Convinience wrapper, casts the the given TensorNetwork _b to Tensor and calls the Tensor function.

Definition at line 1423 of file tensorNetwork.cpp.

◆ approx_equal() [4/4]

bool xerus::approx_equal ( const xerus::Tensor _a,
const xerus::Tensor _b,
const xerus::value_t  _eps = EPSILON 
)

Checks whether two tensors are approximately equal.

Check whether ||_a - _b ||/(||a ||/2 + ||_b ||/2) < _eps, i.e. whether the relative difference in the frobenius norm is sufficently small.

Parameters
_athe first test candidate.
_bthe second test candidate
_epsthe maximal relative difference between _a and _b.
Returns
TRUE if _a and _b are determined to be approximately equal, FALSE otherwise.

Definition at line 1738 of file tensor.cpp.

◆ ASD()

const ALSVariant xerus::ASD ( ,
,
ALSVariant::ASD_solver  ,
false   
)

◆ ASD_SPD()

const ALSVariant xerus::ASD_SPD ( ,
,
ALSVariant::ASD_solver  ,
true   
)

◆ calc_measurment_norm()

template<class MeasurmentSet >
void xerus::calc_measurment_norm ( double *  _norms,
const MeasurmentSet &  _measurments 
)
inline

◆ calc_measurment_norm< RankOneMeasurementSet >()

template<>
void xerus::calc_measurment_norm< RankOneMeasurementSet > ( double *  _norms,
const RankOneMeasurementSet _measurments 
)
inline

Definition at line 556 of file adf.cpp.

◆ calc_measurment_norm< SinglePointMeasurementSet >()

template<>
void xerus::calc_measurment_norm< SinglePointMeasurementSet > ( double *  _norms,
const SinglePointMeasurementSet _measurments 
)
inline

Definition at line 549 of file adf.cpp.

◆ calculate_cq()

void xerus::calculate_cq ( Tensor _C,
Tensor _Q,
Tensor  _input,
const size_t  _splitPos 
)

Low-Level CQ calculation of a given Tensor _input = _C _Q.

This is a rank revealing RQ decomposition with coloum pivoting. In contrast to an RQ the C is not nessecarily upper triangular.

Parameters
_COutput Tensor for the resulting R.
_QOutput Tensor for the resulting Q.
_inputinput Tensor of which the CQ shall be calculated.
_splitPosindex position at defining the matrification for which the CQ is calculated.

Definition at line 1548 of file tensor.cpp.

◆ calculate_factorization_sizes()

XERUS_force_inline std::tuple<size_t, size_t, size_t> xerus::calculate_factorization_sizes ( const Tensor _input,
const size_t  _splitPos 
)

Definition at line 1361 of file tensor.cpp.

◆ calculate_qc()

void xerus::calculate_qc ( Tensor _Q,
Tensor _C,
Tensor  _input,
const size_t  _splitPos 
)

Low-Level QC calculation of a given Tensor _input = _Q _C.

This is a rank revealing QR decomposition with coloum pivoting. In contrast to an QR the C is not nessecarily upper triangular.

Parameters
_QOutput Tensor for the resulting Q.
_COutput Tensor for the resulting R.
_inputinput Tensor of which the QC shall be calculated.
_splitPosindex position at defining the matrification for which the QC is calculated.

Definition at line 1528 of file tensor.cpp.

◆ calculate_qr()

void xerus::calculate_qr ( Tensor _Q,
Tensor _R,
Tensor  _input,
const size_t  _splitPos 
)

Low-Level QR calculation of a given Tensor _input = _Q _R.

Parameters
_QOutput Tensor for the resulting Q.
_ROutput Tensor for the resulting R.
_inputinput Tensor of which the QR shall be calculated.
_splitPosindex position at defining the matrification for which the QR is calculated.

Definition at line 1492 of file tensor.cpp.

◆ calculate_rq()

void xerus::calculate_rq ( Tensor _R,
Tensor _Q,
Tensor  _input,
const size_t  _splitPos 
)

Low-Level RQ calculation of a given Tensor _input = _R _Q.

Parameters
_ROutput Tensor for the resulting R.
_QOutput Tensor for the resulting Q.
_inputinput Tensor of which the RQ shall be calculated.
_splitPosindex position at defining the matrification for which the RQ is calculated.

Definition at line 1511 of file tensor.cpp.

◆ calculate_svd()

void xerus::calculate_svd ( Tensor _U,
Tensor _S,
Tensor _Vt,
Tensor  _input,
const size_t  _splitPos,
const size_t  _maxRank,
const value_t  _eps 
)

Low-Level SVD calculation of a given Tensor _input = _U _S _Vt.

Parameters
_UOutput Tensor for the resulting U.
_SOutput Tensor for the resulting S.
_VtOutput Tensor for the resulting Vt.
_inputinput Tensor of which the SVD shall be calculated.
_splitPosindex position at defining the matrification for which the SVD is calculated.

Definition at line 1424 of file tensor.cpp.

◆ contract() [1/4]

void xerus::contract ( Tensor _result,
const Tensor _lhs,
const bool  _lhsTrans,
const Tensor _rhs,
const bool  _rhsTrans,
const size_t  _numModes 
)

Low-level contraction between Tensors.

Parameters
_resultOutput for the result of the contraction.
_lhsleft hand side of the contraction.
_lhsTransFlags whether the LHS should be transposed (in the matrifications sense).
_rhsright hand side of the contraction.
_rhsTransFlags whether the RHS should be transposed (in the matrifications sense).
_numModesnumber of modes that shall be contracted.

Definition at line 1252 of file tensor.cpp.

◆ contract() [2/4]

Tensor xerus::contract ( const Tensor _lhs,
const bool  _lhsTrans,
const Tensor _rhs,
const bool  _rhsTrans,
const size_t  _numModes 
)

Definition at line 1354 of file tensor.cpp.

◆ contract() [3/4]

XERUS_force_inline void xerus::contract ( Tensor _result,
const Tensor _lhs,
const Tensor _rhs,
const size_t  _numModes 
)

Low-level contraction between Tensors.

Parameters
_resultOutput for the result of the contraction.
_lhsleft hand side of the contraction.
_rhsright hand side of the contraction.
_numModesnumber of indices that shall be contracted.

Definition at line 869 of file tensor.h.

◆ contract() [4/4]

XERUS_force_inline Tensor xerus::contract ( const Tensor _lhs,
const Tensor _rhs,
const size_t  _numModes 
)

Definition at line 873 of file tensor.h.

◆ decomposition_als()

void xerus::decomposition_als ( TTTensor _x,
const Tensor _b,
const double  _eps = EPSILON,
const size_t  _maxIterations = 1000 
)

Definition at line 36 of file decompositionAls.cpp.

◆ DMRG()

const ALSVariant xerus::DMRG ( ,
,
ALSVariant::lapack_solver  ,
false   
)

◆ DMRG_SPD()

const ALSVariant xerus::DMRG_SPD ( ,
,
ALSVariant::lapack_solver  ,
true   
)

◆ dyadic_product() [1/6]

template<bool isOperator>
TTNetwork< isOperator > xerus::dyadic_product ( const TTNetwork< isOperator > &  _lhs,
const TTNetwork< isOperator > &  _rhs 
)

Computes the dyadic product of _lhs and _rhs.

This function is currently needed to keep the resulting network in the TTNetwork class. Apart from that the result is the same as result(i^d1, j^d2) = _lhs(i^d1)*_rhs(j^d2).

Returns
the dyadic product as a TTNetwork.

Definition at line 1319 of file ttNetwork.cpp.

◆ dyadic_product() [2/6]

template<bool isOperator>
TTNetwork< isOperator > xerus::dyadic_product ( const std::vector< TTNetwork< isOperator >> &  _tensors)

Computes the dyadic product of all given TTNetworks.

This is nothing but the repeated application of dyadic_product() for the given TTNetworks.

Returns
the dyadic product as a TTNetwork.

Definition at line 1435 of file ttNetwork.cpp.

◆ dyadic_product() [3/6]

template TTNetwork<true> xerus::dyadic_product ( const TTNetwork< true > &  _lhs,
const TTNetwork< true > &  _rhs 
)

◆ dyadic_product() [4/6]

template TTNetwork<false> xerus::dyadic_product ( const TTNetwork< false > &  _lhs,
const TTNetwork< false > &  _rhs 
)

◆ dyadic_product() [5/6]

template TTNetwork<true> xerus::dyadic_product ( const std::vector< TTNetwork< true >> &  _tensors)

◆ dyadic_product() [6/6]

template TTNetwork<false> xerus::dyadic_product ( const std::vector< TTNetwork< false >> &  _tensors)

◆ entrywise_product() [1/4]

template<bool isOperator>
TTNetwork< isOperator > xerus::entrywise_product ( const TTNetwork< isOperator > &  _A,
const TTNetwork< isOperator > &  _B 
)

Calculates the componentwise product of two tensors given in the TT format.

In general the resulting rank = rank(A)*rank(B). Retains the core position of _A

Definition at line 1275 of file ttNetwork.cpp.

◆ entrywise_product() [2/4]

Tensor xerus::entrywise_product ( const Tensor _A,
const Tensor _B 
)

calculates the entrywise product of two Tensors

Definition at line 1708 of file tensor.cpp.

◆ entrywise_product() [3/4]

template TTNetwork<false> xerus::entrywise_product ( const TTNetwork< false > &  _A,
const TTNetwork< false > &  _B 
)

◆ entrywise_product() [4/4]

template TTNetwork<true> xerus::entrywise_product ( const TTNetwork< true > &  _A,
const TTNetwork< true > &  _B 
)

◆ find_largest_entry() [1/3]

template<bool isOperator>
size_t xerus::find_largest_entry ( const TTNetwork< isOperator > &  _T,
double  _accuracy,
value_t  _lowerBound = 0.0 
)

Finds the position of the approximately largest entry.

Finds an entry that is at least of size _accuracy * X_max in absolute value, where X_max is the largest entry of the tensor. The smaller _accuracy, the faster the algorithm will work.

Parameters
_accuracyfactor that determains the maximal deviation of the returned entry from the true largest entry.
_lowerBounda lower bound for the largest entry, i.e. there must be an entry in the tensor which is at least of this size (in absolute value). The algorithm may fail completely if this is not fullfilled, but will work using its own approximation if no value (i.e. 0.0) is given.
Returns
the position of the entry found.

Definition at line 6 of file largestEntry.cpp.

◆ find_largest_entry() [2/3]

template size_t xerus::find_largest_entry ( const TTNetwork< true > &  ,
double  ,
value_t   
)

◆ find_largest_entry() [3/3]

template size_t xerus::find_largest_entry ( const TTNetwork< false > &  ,
double  ,
value_t   
)

◆ frob_norm() [1/2]

static XERUS_force_inline value_t xerus::frob_norm ( const TensorNetwork _network)
static

Calculates the frobenious norm of the given TensorNetwork.

Parameters
_networkthe TensorNetwork of which the frobenious norm shall be calculated.
Returns
the frobenious norm.

Definition at line 534 of file tensorNetwork.h.

◆ frob_norm() [2/2]

static XERUS_force_inline value_t xerus::frob_norm ( const Tensor _tensor)
static

Calculates the frobenius norm of the given tensor.

Parameters
_tensorthe Tensor of which the frobenious norm shall be calculated.
Returns
the frobenius norm .

Definition at line 931 of file tensor.h.

◆ GeometricCG()

const GeometricCGVariant xerus::GeometricCG ( ,
1e-  8,
false  ,
SubmanifoldRetractionI  ,
ProjectiveVectorTransport   
)

◆ get_grouped_entries()

template<bool isOperator>
std::vector<std::vector<std::tuple<size_t, size_t, value_t> > > xerus::get_grouped_entries ( const Tensor _component)

Definition at line 496 of file ttNetwork.cpp.

◆ get_step_sizes()

std::vector<size_t> xerus::get_step_sizes ( const Tensor::DimensionTuple _dimensions)
inline

Definition at line 958 of file tensor.cpp.

◆ HOSVDRetractionI()

void xerus::HOSVDRetractionI ( TTTensor _U,
const TTTangentVector _change 
)

retraction that performs a HOSVD to project back onto the Manifold

Definition at line 56 of file retractions.cpp.

◆ HOSVDRetractionII()

void xerus::HOSVDRetractionII ( TTTensor _U,
const TTTensor _change 
)

Definition at line 50 of file retractions.cpp.

◆ IHT()

double xerus::IHT ( TTTensor _x,
const SinglePointMeasurementSet _measurments,
PerformanceData _perfData = NoPerfData 
)

Definition at line 29 of file iht.cpp.

◆ increase_indices()

void xerus::increase_indices ( const size_t  _i,
value_t *&  _outPosition,
const std::vector< size_t > &  _stepSizes,
const std::vector< size_t > &  _baseDimensions 
)

Definition at line 39 of file indexedTensor_tensor_evaluate.cpp.

◆ line_search()

void xerus::line_search ( TTTensor _x,
value_t _alpha,
const TTTangentVector _direction,
value_t  _derivative,
value_t _residual,
std::function< void(TTTensor &, const TTTangentVector &)>  _retraction,
std::function< value_t()>  _calculate_residual,
value_t  _changeInAlpha = 0.5 
)

Definition at line 37 of file steepestDescent.cpp.

◆ matrix_matrix_product() [1/6]

void xerus::matrix_matrix_product ( double *const  _C,
const size_t  _leftDim,
const size_t  _rightDim,
const double  _alpha,
const std::map< size_t, double > &  _A,
const bool  _transposeA,
const size_t  _midDim,
const double *const  _B,
const bool  _transposeB 
)

Definition at line 97 of file sparseTimesFullContraction.cpp.

◆ matrix_matrix_product() [2/6]

void xerus::matrix_matrix_product ( double *const  _C,
const size_t  _leftDim,
const size_t  _rightDim,
const double  _alpha,
const double *const  _A,
const bool  _transposeA,
const size_t  _midDim,
const std::map< size_t, double > &  _B,
const bool  _transposeB 
)

Definition at line 114 of file sparseTimesFullContraction.cpp.

◆ matrix_matrix_product() [3/6]

void xerus::matrix_matrix_product ( std::map< size_t, double > &  _C,
const size_t  _leftDim,
const size_t  _rightDim,
const double  _alpha,
const std::map< size_t, double > &  _A,
const bool  _transposeA,
const size_t  _midDim,
const double *const  _B,
const bool  _transposeB 
)

Definition at line 183 of file sparseTimesFullContraction.cpp.

◆ matrix_matrix_product() [4/6]

void xerus::matrix_matrix_product ( std::map< size_t, double > &  _C,
const size_t  _leftDim,
const size_t  _rightDim,
const double  _alpha,
const double *const  _A,
const bool  _transposeA,
const size_t  _midDim,
const std::map< size_t, double > &  _B,
const bool  _transposeB 
)

Definition at line 210 of file sparseTimesFullContraction.cpp.

◆ matrix_matrix_product() [5/6]

void xerus::matrix_matrix_product ( double *const  _C,
const size_t  _leftDim,
const size_t  _rightDim,
const double  _alpha,
const std::map< size_t, double > &  _A,
const bool  _transposeA,
const size_t  _midDim,
const double *const  _B 
)

Definition at line 66 of file sparseTimesFullContraction.cpp.

◆ matrix_matrix_product() [6/6]

void xerus::matrix_matrix_product ( std::map< size_t, double > &  _C,
const size_t  ,
const size_t  _rightDim,
const double  _alpha,
const std::map< size_t, double > &  _A,
const size_t  _midDim,
const double *const  _B 
)

Definition at line 132 of file sparseTimesFullContraction.cpp.

◆ NoPerfData()

PerformanceData xerus::NoPerfData ( false  ,
false   
)

◆ one_norm()

static XERUS_force_inline value_t xerus::one_norm ( const Tensor _tensor)
static

Calculates the 1-norm of the given tensor.

Parameters
_tensorthe Tensor of which the norm shall be calculated.
Returns
the 1-norm

Definition at line 938 of file tensor.h.

◆ operator!=()

bool xerus::operator!= ( const Index _a,
const Index _b 
)

Two Indices are equal if their valueId coincides. Fixed indices are never equal.

Definition at line 167 of file index.cpp.

◆ operator*() [1/11]

TTTangentVector xerus::operator* ( value_t  _alpha,
const TTTangentVector _rhs 
)

Definition at line 161 of file retractions.cpp.

◆ operator*() [2/11]

TensorNetwork xerus::operator* ( TensorNetwork _lhs,
value_t  _factor 
)

Definition at line 1392 of file tensorNetwork.cpp.

◆ operator*() [3/11]

TensorNetwork xerus::operator* ( value_t  _factor,
TensorNetwork _rhs 
)

Definition at line 1398 of file tensorNetwork.cpp.

◆ operator*() [4/11]

template<bool isOperator>
TTNetwork< isOperator > xerus::operator* ( TTNetwork< isOperator >  _network,
const value_t  _factor 
)

Calculates the entrywise multiplication of the given TTNetwork _network with a constant _factor.

Internally this only results in a change in the global factor.

Parameters
_networkthe TTNetwork,
_factorthe factor,
Returns
the resulting scaled TTNetwork.

Definition at line 1175 of file ttNetwork.cpp.

◆ operator*() [5/11]

template<bool isOperator>
TTNetwork< isOperator > xerus::operator* ( const value_t  _factor,
TTNetwork< isOperator >  _network 
)

Calculates the entrywise multiplication of the given TTNetwork _network with a constant _factor.

Internally this only results in a change in the global factor.

Parameters
_factorthe factor,
_networkthe TTNetwork,
Returns
the resulting scaled TTNetwork.

Definition at line 1182 of file ttNetwork.cpp.

◆ operator*() [6/11]

Tensor xerus::operator* ( const value_t  _factor,
Tensor  _tensor 
)

Calculates the entrywise multiplication of the Tensor _tensor with the constant _factor.

Internally this only results in a change in the global factor.

Parameters
_factorthe factor to be used.
_tensorthe Tensor that shall be scaled.
Returns
the resulting scaled Tensor.

Definition at line 1233 of file tensor.cpp.

◆ operator*() [7/11]

Tensor xerus::operator* ( Tensor  _tensor,
const value_t  _factor 
)

Calculates the entrywise multiplication of the Tensor _tensor with the constant _factor.

Internally this only results in a change in the global factor.

Parameters
_tensorthe Tensor that shall be scaled.
_factorthe factor to be used.
Returns
the resulting scaled Tensor.

Definition at line 1239 of file tensor.cpp.

◆ operator*() [8/11]

template TTNetwork<false> xerus::operator* ( TTNetwork< false >  _network,
const value_t  _factor 
)

◆ operator*() [9/11]

template TTNetwork<true> xerus::operator* ( TTNetwork< true >  _network,
const value_t  _factor 
)

◆ operator*() [10/11]

template TTNetwork<false> xerus::operator* ( const value_t  _factor,
TTNetwork< false >  _network 
)

◆ operator*() [11/11]

template TTNetwork<true> xerus::operator* ( const value_t  _factor,
TTNetwork< true >  _network 
)

◆ operator+() [1/4]

template<bool isOperator>
TTNetwork< isOperator > xerus::operator+ ( TTNetwork< isOperator >  _lhs,
const TTNetwork< isOperator > &  _rhs 
)

Calculates the entrywise sum of the given TTNetworks _lhs and _rhs.

To be well-defined it is required that the dimensions of _lhs and _rhs coincide. The rank of the result are in general the entrywise sum of the ranks of _lhs and _rhs.

Parameters
_lhsthe first summand.
_rhsthe second summand.
Returns
the sum.

Definition at line 1161 of file ttNetwork.cpp.

◆ operator+() [2/4]

Tensor xerus::operator+ ( Tensor  _lhs,
const Tensor _rhs 
)

Calculates the entrywise sum of _lhs and _rhs.

To be well-defined it is required that the dimensions of _lhs and _rhs coincide.

Parameters
_lhsthe first summand.
_rhsthe second summand.
Returns
the sum.

Definition at line 1221 of file tensor.cpp.

◆ operator+() [3/4]

template TTNetwork<false> xerus::operator+ ( TTNetwork< false >  _lhs,
const TTNetwork< false > &  _rhs 
)

◆ operator+() [4/4]

template TTNetwork<true> xerus::operator+ ( TTNetwork< true >  _lhs,
const TTNetwork< true > &  _rhs 
)

◆ operator-() [1/4]

template<bool isOperator>
TTNetwork< isOperator > xerus::operator- ( TTNetwork< isOperator >  _lhs,
const TTNetwork< isOperator > &  _rhs 
)

Calculates the entrywise difference of the given TTNetworks _lhs and _rhs.

To be well-defined it is required that the dimensions of _lhs and _rhs coincide. The rank of the result are in general the entrywise sum of the ranks of _lhs and _rhs.

Parameters
_lhsthe minuend.
_rhsthe subtrahend.
Returns
the difference.

Definition at line 1168 of file ttNetwork.cpp.

◆ operator-() [2/4]

Tensor xerus::operator- ( Tensor  _lhs,
const Tensor _rhs 
)

Calculates the entrywise difference between _lhs and _rhs.

To be well-defined it is required that the dimensions of _lhs and _rhs coincide.

Parameters
_lhsthe minuend.
_rhsthe subtrahend.
Returns
the difference.

Definition at line 1227 of file tensor.cpp.

◆ operator-() [3/4]

template TTNetwork<false> xerus::operator- ( TTNetwork< false >  _lhs,
const TTNetwork< false > &  _rhs 
)

◆ operator-() [4/4]

template TTNetwork<true> xerus::operator- ( TTNetwork< true >  _lhs,
const TTNetwork< true > &  _rhs 
)

◆ operator/() [1/6]

◆ operator/() [2/6]

TensorNetwork xerus::operator/ ( TensorNetwork _lhs,
value_t  _factor 
)

Definition at line 1404 of file tensorNetwork.cpp.

◆ operator/() [3/6]

template<bool isOperator>
TTNetwork< isOperator > xerus::operator/ ( TTNetwork< isOperator >  _network,
const value_t  _divisor 
)

Calculates the entrywise divison of this TTNetwork by a constant _divisor.

Internally this only results in a change in the global factor.

Parameters
_divisorthe divisor,
Returns
the resulting scaled TTNetwork.

Definition at line 1189 of file ttNetwork.cpp.

◆ operator/() [4/6]

Tensor xerus::operator/ ( Tensor  _tensor,
const value_t  _divisor 
)

Calculates the entrywise divison of the Tensor _tensor with the constant _divisor.

Internally this only results in a change in the global factor.

Parameters
_tensorthe Tensor that shall be scaled.
_divisorthe factor to be used.
Returns
the resulting scaled Tensor.

Definition at line 1245 of file tensor.cpp.

◆ operator/() [5/6]

template TTNetwork<false> xerus::operator/ ( TTNetwork< false >  _network,
const value_t  _divisor 
)

◆ operator/() [6/6]

template TTNetwork<true> xerus::operator/ ( TTNetwork< true >  _network,
const value_t  _divisor 
)

◆ operator<()

bool xerus::operator< ( const Index _a,
const Index _b 
)

The Comparision operator is needed for indices to be orderable in std::set, the valueId is used.

Definition at line 171 of file index.cpp.

◆ operator<<() [1/3]

std::ostream & xerus::operator<< ( std::ostream &  _out,
const xerus::Index _idx 
)

Allows to pretty print indices, giving the valueId and span.

Definition at line 175 of file index.cpp.

◆ operator<<() [2/3]

std::ostream & xerus::operator<< ( std::ostream &  _out,
const TensorNetwork::Link _rhs 
)

Definition at line 35 of file link.cpp.

◆ operator<<() [3/3]

std::ostream & xerus::operator<< ( std::ostream &  _out,
const Tensor _tensor 
)

Prints the Tensor to the given outStream.

Parameters
_outthe outstream to be printed to.
_tensorthe tensor.
Returns
_out.

Definition at line 1773 of file tensor.cpp.

◆ operator==()

bool xerus::operator== ( const Index _a,
const Index _b 
)

Two Indices are equal if their valueId coincides. Fixed indices are never equal.

Definition at line 163 of file index.cpp.

◆ perform_component_product()

template<bool isOperator>
void xerus::perform_component_product ( Tensor _newComponent,
const Tensor _componentA,
const Tensor _componentB 
)

Definition at line 1210 of file ttNetwork.cpp.

◆ position_or_zero()

template<class MeasurmentSet >
size_t xerus::position_or_zero ( const MeasurmentSet &  _measurments,
const size_t  _meas,
const size_t  _corePosition 
)
inline

◆ position_or_zero< RankOneMeasurementSet >()

template<>
size_t xerus::position_or_zero< RankOneMeasurementSet > ( const RankOneMeasurementSet ,
const size_t  ,
const size_t   
)
inline

Definition at line 407 of file adf.cpp.

◆ position_or_zero< SinglePointMeasurementSet >()

template<>
size_t xerus::position_or_zero< SinglePointMeasurementSet > ( const SinglePointMeasurementSet _measurments,
const size_t  _meas,
const size_t  _corePosition 
)
inline

Definition at line 402 of file adf.cpp.

◆ prepare_factorization_output()

XERUS_force_inline void xerus::prepare_factorization_output ( Tensor _lhs,
Tensor _rhs,
const Tensor _input,
const size_t  _splitPos,
const size_t  _rank,
Tensor::Representation  _rep 
)

Definition at line 1371 of file tensor.cpp.

◆ prepare_split()

std::unique_ptr<Tensor> xerus::prepare_split ( size_t &  _lhsSize,
size_t &  _rhsSize,
size_t &  _rank,
size_t &  _splitPos,
std::vector< Index > &  _lhsPreliminaryIndices,
std::vector< Index > &  _rhsPreliminaryIndices,
internal::IndexedTensorReadOnly< Tensor > &&  _base,
internal::IndexedTensor< Tensor > &&  _lhs,
internal::IndexedTensor< Tensor > &&  _rhs 
)

Definition at line 36 of file indexedTensor_tensor_factorisations.cpp.

◆ ProjectiveVectorTransport()

void xerus::ProjectiveVectorTransport ( const TTTensor _newBase,
TTTangentVector _tangentVector 
)

simple vector transport by projecting onto the new tangent plane

Definition at line 283 of file retractions.cpp.

◆ pseudo_inverse() [1/2]

void xerus::pseudo_inverse ( Tensor _inverse,
const Tensor _input,
const size_t  _splitPos 
)

Low-Level calculation of the pseudo inverse of a given Tensor.

Currently doen by calculation of the SVD.

Parameters
_inverseOutput the pseudo inverse.
_inputinput Tensor of which the CQ shall be calculated.
_splitPosindex position at defining the matrification for which the pseudo inverse is calculated.

Definition at line 1568 of file tensor.cpp.

◆ pseudo_inverse() [2/2]

Tensor xerus::pseudo_inverse ( const Tensor _input,
const size_t  _splitPos 
)

Definition at line 1576 of file tensor.cpp.

◆ randVar_to_position()

Tensor xerus::randVar_to_position ( const double  _v,
const size_t  _polyDegree 
)

Definition at line 40 of file uqAdf.cpp.

◆ reshuffle() [1/2]

void xerus::reshuffle ( Tensor _out,
const Tensor _base,
const std::vector< size_t > &  _shuffle 
)

: Performs a simple reshuffle. Much less powerfull then a full evaluate, but more efficient.

_shuffle shall be a vector that gives for every old index, its new position.

Definition at line 55 of file indexedTensor_tensor_evaluate.cpp.

◆ reshuffle() [2/2]

Tensor xerus::reshuffle ( const Tensor _base,
const std::vector< size_t > &  _shuffle 
)

Definition at line 139 of file indexedTensor_tensor_evaluate.cpp.

◆ set_factorization_output() [1/2]

XERUS_force_inline void xerus::set_factorization_output ( Tensor _lhs,
std::unique_ptr< value_t[]> &&  _lhsData,
Tensor _rhs,
std::unique_ptr< value_t[]> &&  _rhsData,
const Tensor _input,
const size_t  _splitPos,
const size_t  _rank 
)

Definition at line 1398 of file tensor.cpp.

◆ set_factorization_output() [2/2]

XERUS_force_inline void xerus::set_factorization_output ( Tensor _lhs,
std::map< size_t, value_t > &&  _lhsData,
Tensor _rhs,
std::map< size_t, value_t > &&  _rhsData,
const Tensor _input,
const size_t  _splitPos,
const size_t  _rank 
)

Definition at line 1411 of file tensor.cpp.

◆ solve() [1/2]

◆ solve() [2/2]

void xerus::solve ( Tensor _X,
const Tensor _A,
const Tensor _B,
size_t  _extraDegree = 0 
)

Solves the equation Ax = b for x. If the solution is not unique, the output need not be the minimal norm solution.

Parameters
_XOutput Tensor for the result
_Ainput Operator A
_Binput right-hand-side b
_extraDegreenumber of modes that _x and _B sharefor which the solution should be computed independently.

Definition at line 1654 of file tensor.cpp.

◆ solve_least_squares()

void xerus::solve_least_squares ( Tensor _X,
const Tensor _A,
const Tensor _B,
const size_t  _extraDegree = 0 
)

Solves the least squares problem ||_A _X - _B||_F.

Parameters
_XOutput Tensor for the resulting X.
_Ainput Tensor A.
_Binput Tensor b.
_extraDegreenumber of modes that _X and _B share and for which the least squares problem is independently solved.

Definition at line 1583 of file tensor.cpp.

◆ SteepestDescent()

const SteepestDescentVariant xerus::SteepestDescent ( ,
1e-  8,
false  ,
SubmanifoldRetractionII   
)

◆ SubmanifoldRetractionI()

void xerus::SubmanifoldRetractionI ( TTTensor _U,
const TTTangentVector _change 
)

retraction that performs componentwise addition of \( U_i \) and \( W_i \) where \( W_i \) is the i-th component of the riemannian tangential vector representation

Definition at line 269 of file retractions.cpp.

◆ SubmanifoldRetractionII()

void xerus::SubmanifoldRetractionII ( TTTensor _U,
const TTTensor _change 
)

retraction that performs componentwise addition of \( U_i \) and \( W_i \) where \( W_i \) is the i-th component of the riemannian tangential vector representation

Definition at line 264 of file retractions.cpp.

◆ transpose() [1/4]

XERUS_force_inline void xerus::transpose ( double *const __restrict  _out,
const double *const __restrict  _in,
const size_t  _leftDim,
const size_t  _rightDim 
)

Definition at line 35 of file sparseTimesFullContraction.cpp.

◆ transpose() [2/4]

XERUS_force_inline std::unique_ptr<double[]> xerus::transpose ( const double *const  _A,
const size_t  _leftDim,
const size_t  _rightDim 
)

Definition at line 43 of file sparseTimesFullContraction.cpp.

◆ transpose() [3/4]

XERUS_force_inline void xerus::transpose ( std::map< size_t, double > &__restrict  _out,
const std::map< size_t, double > &__restrict  _in,
const size_t  _leftDim,
const size_t  _rightDim 
)

Definition at line 50 of file sparseTimesFullContraction.cpp.

◆ transpose() [4/4]

XERUS_force_inline std::map<size_t, double> xerus::transpose ( const std::map< size_t, double > &  _A,
const size_t  _leftDim,
const size_t  _rightDim 
)

Definition at line 58 of file sparseTimesFullContraction.cpp.

◆ transpose_if_operator()

template<bool isOperator>
void xerus::transpose_if_operator ( TTNetwork< isOperator > &  _ttNetwork)

◆ transpose_if_operator< false >()

template<>
void xerus::transpose_if_operator< false > ( TTNetwork< false > &  )

Definition at line 974 of file ttNetwork.cpp.

◆ transpose_if_operator< true >()

template<>
void xerus::transpose_if_operator< true > ( TTNetwork< true > &  _ttNetwork)

Definition at line 977 of file ttNetwork.cpp.

◆ uq_adf() [1/2]

void xerus::uq_adf ( TTTensor _x,
const std::vector< std::vector< double >> &  _randomVariables,
const std::vector< Tensor > &  _solutions 
)

Definition at line 327 of file uqAdf.cpp.

◆ uq_adf() [2/2]

TTTensor xerus::uq_adf ( const UQMeasurementSet _measurments,
const TTTensor _guess 
)

Definition at line 334 of file uqAdf.cpp.

◆ uq_avg()

Tensor xerus::uq_avg ( const TTTensor _x,
const size_t  _N,
const size_t  _numSpecial 
)

Definition at line 451 of file uqAdf.cpp.

◆ uq_mc()

std::pair< std::vector< std::vector< double > >, std::vector< Tensor > > xerus::uq_mc ( const TTTensor _x,
const size_t  _N,
const size_t  _numSpecial 
)

Definition at line 425 of file uqAdf.cpp.

Variable Documentation

◆ ADF

const ADFVariant xerus::ADF(0, 1e-8, 0.999)

Default variant of the ADF algorithm.

◆ ALS

const ALSVariant xerus::ALS(1, 0, ALSVariant::lapack_solver, false)

default variant of the single-site ALS algorithm for non-symmetric operators using the lapack solver

◆ ALS_SPD

const ALSVariant xerus::ALS_SPD(1, 0, ALSVariant::lapack_solver, true)

default variant of the single-site ALS algorithm for symmetric positive-definite operators using the lapack solver

◆ ASD

const ALSVariant xerus::ASD(1, 0, ALSVariant::ASD_solver, false)

default variant of the alternating steepest descent for non-symmetric operators

◆ ASD_SPD

const ALSVariant xerus::ASD_SPD(1, 0, ALSVariant::ASD_solver, true)

default variant of the alternating steepest descent for symmetric positive-definite operators

◆ DMRG

const ALSVariant xerus::DMRG(2, 0, ALSVariant::lapack_solver, false)

default variant of the two-site DMRG algorithm for non-symmetric operators using the lapack solver

◆ DMRG_SPD

const ALSVariant xerus::DMRG_SPD(2, 0, ALSVariant::lapack_solver, true)

default variant of the two-site DMRG algorithm for symmetric positive-definite operators using the lapack solver

◆ EPSILON

constexpr const value_t xerus::EPSILON = 8*std::numeric_limits<value_t>::epsilon()

The default epsilon value in xerus.

This value is used as a default value by all xerus routines that work up to a certain accuracy, e.g. the singluar value decompostion or the ALS algorithm.

Definition at line 50 of file basic.h.

◆ GeometricCG

const GeometricCGVariant xerus::GeometricCG(0, 1e-8, false, SubmanifoldRetractionI, ProjectiveVectorTransport)

default variant of the steepest descent algorithm using the lapack solver

◆ NoPerfData

PerformanceData xerus::NoPerfData(false, false)

◆ SteepestDescent

const SteepestDescentVariant xerus::SteepestDescent(0, 1e-8, false, SubmanifoldRetractionII)

default variant of the steepest descent algorithm using the lapack solver

◆ VERSION_COMMIT

const int xerus::VERSION_COMMIT = XERUS_VERSION_COMMIT

Definition at line 31 of file standard.cpp.

◆ VERSION_MAJOR

const int xerus::VERSION_MAJOR = XERUS_VERSION_MAJOR

The version of the compiled xerus library

Definition at line 28 of file standard.cpp.

◆ VERSION_MINOR

const int xerus::VERSION_MINOR = XERUS_VERSION_MINOR

Definition at line 29 of file standard.cpp.

◆ VERSION_REVISION

const int xerus::VERSION_REVISION = XERUS_VERSION_REVISION

Definition at line 30 of file standard.cpp.