mpqc::math Namespace Reference

Namespaces

 clustering
 
 detail
 
 groups
 

Classes

struct  compress
 
class  CP3
 Computes 3-way CP decomposition of a order-3 tensor. More...
 
class  DavidsonDiag
 Davidson Algorithm. More...
 
class  DavidsonDiagPred
 
class  DecomposedTensor
 
class  DecompToTaTensor
 
class  DIIS
 Factory for creating DIIS objects. More...
 
class  FiniteDifferenceDerivative
 Computes finite-difference approximation to function derivatives. More...
 
class  Function
 Function maps Parameters to a Value. More...
 
class  FunctionVisitorBase
 
class  Group
 Group is an abstract discrete group. More...
 
struct  mat_to_tile< TA::Tile< DecomposedTensor< T > > >
 
struct  mat_to_tile< Tile, std::enable_if_t< TA::detail::is_contiguous_tensor_v< Tile > > >
 
class  Optimizer
 It is itself a Function, namely Function<Params, std::tuple<>>. More...
 
class  pair_accumulator
 
class  pair_smasher
 
class  PetiteList
 A helper class for symmetric tensor algebra. More...
 
class  QuasiNewtonOptimizer
 It is itself a Function, namely Function<Params, OptimizerParams>. More...
 
class  SingleStateDavidsonDiag
 
class  SymmPetiteList
 
class  TaToDecompTensor
 
class  TaylorExpansionCoefficients
 N-th order Taylor expansion of a function of K variables. More...
 
class  TaylorExpansionFunction
 
class  Tile
 
class  VectorCluster
 

Typedefs

template<typename T >
using Matrix = Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor >
 
template<typename T >
using Vector = Eigen::Matrix< T, Eigen::Dynamic, 1 >
 

Functions

template<typename T1 , typename T2 >
constexpr std::complex< double > Z (T1 re, T2 im)
 
template<typename T >
RowMatrix< T > random_unitary (int size, int seed=42)
 
template<typename T >
RowMatrix< T > complex_random_unitary (int size)
 
template<typename Value >
std::ostream & operator<< (std::ostream &os, const TaylorExpansionCoefficients< Value > &x)
 
std::string to_string (PetiteList::Symmetry symmetry)
 
template<typename Tile , typename Policy >
TA::DistArray< Tile, Policy > cholesky (TA::DistArray< Tile, Policy > const &A)
 
template<typename Tile , typename Policy >
TA::DistArray< Tile, Policy > cholesky_inverse (TA::DistArray< Tile, Policy > const &A)
 
template<typename Tile , typename Policy >
TA::DistArray< Tile, Policy > eigen_inverse (const TA::DistArray< Tile, Policy > &A)
 
template<typename T >
std::tuple< RowMatrix< T >, RowMatrix< T >, size_t, double, double > gensqrtinv (const RowMatrix< T > &S, bool symmetric=false, double max_condition_number=1e8)
 
template<typename Tile , typename Policy >
TA::DistArray< Tile, Policy > conditioned_orthogonalizer (TA::DistArray< Tile, Policy > S_array, double S_condition_number_threshold)
 
template<typename Tile , typename Policy >
TA::DistArray< Tile, Policy > create_diagonal_array_from_eigen (madness::World &world, const TA::TiledRange1 &trange1, const TA::TiledRange1 &trange2, typename Tile::numeric_type val)
 
template<typename Tile , typename Scalar , std::enable_if_t< TA::detail::is_scalar_v< Scalar >> * = nullptr>
void make_diagonal_tile (Tile &tile, Scalar val)
 
template<typename T >
void make_diagonal_tile (TA::Tile< DecomposedTensor< T >> &tile, T val)
 
template<typename Tile , typename Policy >
TiledArray::DistArray< Tile, typename std::enable_if< std::is_same< Policy, TA::SparsePolicy >::value, TA::SparsePolicy >::type > create_diagonal_matrix (TiledArray::DistArray< Tile, Policy > const &model, double val)
 
template<typename Tile , typename Policy >
TiledArray::DistArray< Tile, typename std::enable_if< std::is_same< Policy, TA::DensePolicy >::value, TA::DensePolicy >::type > create_diagonal_matrix (TiledArray::DistArray< Tile, Policy > const &model, double val)
 
template<typename Tile >
TiledArray::DistArray< Tile, TiledArray::SparsePolicy > diagonal_matrix (TiledArray::TiledRange const &trange, double val, madness::World &world)
 takes a TiledArray::TiledRange and a value and returns a diagonal matrix. More...
 
template<typename T , unsigned int N, typename Tile >
TiledArray::Array< T, N, Tile, TiledArray::DensePolicy > create_diagonal_matrix (TiledArray::Array< T, N, Tile, TiledArray::DensePolicy > const &model, double val)
 
template<typename It >
double min_eval_guess (It first, It second)
 
template<typename Array >
std::pair< double, double > symmetric_min_max_evals (Array const &S)
 
template<typename Array >
Array init_X (Array const &S)
 
template<typename Array >
Array invert (Array const &S)
 
template<typename Array >
double min_eval_est (Array const &H, Array const &S)
 
template<typename D >
void gram_schmidt (std::vector< D > &V, double threshold, std::size_t start=0)
 
template<typename D >
void gram_schmidt (const std::vector< D > &V1, std::vector< D > &V2, double threshold)
 
template<typename Tile , typename Policy >
TA::DistArray< Tile, Policy > pseudoinverse (TA::DistArray< Tile, Policy > const &A, bool HPSD=false)
 
template<typename Tile , std::enable_if_t< TiledArray::detail::is_contiguous_tensor_v< Tile >> * = nullptr>
void eigen_estimator (std::array< Tile, 2 > &result, Tile const &tile)
 
template<typename Array >
double max_eval_est (Array const &S)
 
template<typename Scalar , typename Tensor , std::enable_if_t< TA::detail::is_contiguous_tensor_v< Tensor >> * = nullptr>
void add_to_diag_tile (Scalar val, Tensor &tile)
 
void add_to_diag_tile (double val, TA::Tile< DecomposedTensor< double >> &tile)
 
template<typename Array >
void add_to_diag (Array &A, double val)
 
template<typename Array >
std::array< typename Array::value_type::numeric_type, 2 > eval_guess (Array const &A)
 computes {min,max} eigenvalue bounds for a square matrix A using Gerschgorin circle theorem More...
 
template<typename Array >
void third_order_update (Array const &S, Array &Z)
 
template<typename Array >
Array inverse_sqrt (Array const &S)
 
void gauss_legendre (int N, Eigen::VectorXd &w, Eigen::VectorXd &x)
 this function uses orthogonal polynomial approach for obtaining quadrature roots and weights. Taken from paper: Gaussian Quadrature and the Eigenvalue Problem by John A. Gubner. http://gubner.ece.wisc.edu/gaussquad.pdf The code is outlined at Example 15: Legendre polynomials More...
 
template<typename T >
Matrix< T > tile_to_eigen (TA::Tensor< T > const &t)
 
template<typename T >
Matrix< T > tile_to_eigen (TA::Tile< DecomposedTensor< T >> const &t)
 
template<typename Tile , typename Policy >
Matrix< typename Tile::value_typearray_to_eigen (TA::DistArray< Tile, Policy > const &A)
 converts a TiledArray::Array to a row-major Eigen Matrix More...
 
template<typename Tile , typename Policy >
std::vector< Matrix< typename Tile::value_type > > array_to_eigen (TA::DistArrayVector< Tile, Policy > const &A)
 converts a TiledArray::ArrayVector to a row-major Eigen Matrix More...
 
template<typename Tile , typename Policy >
TA::DistArray< Tile, typename std::enable_if< std::is_same< Policy, TA::SparsePolicy >::value, TA::SparsePolicy >::type > eigen_to_array (madness::World &world, Matrix< typename Tile::numeric_type > const &M, TA::TiledRange1 tr0, TA::TiledRange1 tr1, double cut=1e-7)
 
template<typename Tile , typename Policy >
TA::DistArray< Tile, typename std::enable_if< std::is_same< Policy, TA::DensePolicy >::value, TA::DensePolicy >::type > eigen_to_array (madness::World &world, Matrix< typename Tile::numeric_type > const &M, TA::TiledRange1 tr0, TA::TiledRange1 tr1, double cut=1e-7)
 
template<typename T >
std::ostream & operator<< (std::ostream &os, DecomposedTensor< T > const &t)
 
template<typename T >
unsigned long tile_clr_storage (TiledArray::Tile< math::DecomposedTensor< T >> const &tile)
 
template<typename T >
DecomposedTensor< T > add_order2 (DecomposedTensor< T > const &l, DecomposedTensor< T > const &r, T factor)
 
template<typename T >
DecomposedTensor< T > add (DecomposedTensor< T > const &l, DecomposedTensor< T > const &r)
 
template<typename T , typename Perm , typename = std::enable_if_t<TiledArray::detail::is_permutation_v<Perm>>>
DecomposedTensor< T > add (DecomposedTensor< T > const &l, DecomposedTensor< T > const &r, Perm const &p)
 
template<typename T >
DecomposedTensor< T > add (DecomposedTensor< T > const &l, DecomposedTensor< T > const &r, const T factor)
 
template<typename T , typename Perm , typename = std::enable_if_t<TiledArray::detail::is_permutation_v<Perm>>>
DecomposedTensor< T > add (DecomposedTensor< T > const &l, DecomposedTensor< T > const &r, const T factor, Perm const &p)
 
template<typename T >
DecomposedTensor< T > add (DecomposedTensor< T > const &l, const T factor)
 
template<typename T , typename Perm , typename = std::enable_if_t<TiledArray::detail::is_permutation_v<Perm>>>
DecomposedTensor< T > add (DecomposedTensor< T > const &l, const T factor, Perm const &p)
 
template<typename T >
DecomposedTensor< T > & add_to_order2 (DecomposedTensor< T > &l, DecomposedTensor< T > const &r)
 
template<typename T >
DecomposedTensor< T > & add_to (DecomposedTensor< T > &l, DecomposedTensor< T > const &r)
 
template<typename T >
DecomposedTensor< T > & add_to (DecomposedTensor< T > &l, DecomposedTensor< T > const &r, const T factor)
 
template<typename T >
DecomposedTensor< T > & add_to (DecomposedTensor< T > &l, const T factor)
 
void recompress (DecomposedTensor< double > &t)
 Currently modifies input data regardless could cause some loss of accuracy. More...
 
DecomposedTensor< double > two_way_decomposition (DecomposedTensor< double > const &t)
 Returns an empty DecomposedTensor if the compression rank was to large. More...
 
TA::Tensor< double > combine (DecomposedTensor< double > const &t)
 
void piv_cholesky (Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > &a)
 performs the pivoted Cholesky decomposition More...
 
template<typename T , typename Op >
DecomposedTensor< T > binary (DecomposedTensor< T > const &l, DecomposedTensor< T > const &r, Op &&op)
 
template<typename T , typename Op , typename Perm , typename = std::enable_if_t<TiledArray::detail::is_permutation_v<Perm>>>
DecomposedTensor< T > binary (DecomposedTensor< T > const &l, DecomposedTensor< T > const &r, Op &&op, Perm const &p)
 
template<typename T , typename Op >
DecomposedTensor< T > & inplace_binary (DecomposedTensor< T > &l, DecomposedTensor< T > const &r, Op &&op)
 
template<typename T >
TA::TensorD gemm (DecomposedTensor< T > const &a, TA::TensorD const &b, const T factor, TA::math::GemmHelper const &gh)
 
template<typename T >
TA::TensorD gemm (TA::TensorD &c, DecomposedTensor< T > const &a, TA::TensorD const &b, const T factor, TA::math::GemmHelper const &gh)
 
template<typename T >
DecomposedTensor< T > gemm (DecomposedTensor< T > const &a, DecomposedTensor< T > const &b, const T factor, TA::math::GemmHelper const &gh)
 
template<typename T >
DecomposedTensor< T > & gemm (DecomposedTensor< T > &c, DecomposedTensor< T > const &a, DecomposedTensor< T > const &b, const T factor, TA::math::GemmHelper const &gh)
 
template<typename T >
DecomposedTensor< T > mult (DecomposedTensor< T > const &l, DecomposedTensor< T > const &r)
 
template<typename T , typename Perm , typename = std::enable_if_t<TiledArray::detail::is_permutation_v<Perm>>>
DecomposedTensor< T > mult (DecomposedTensor< T > const &l, DecomposedTensor< T > const &r, Perm const &p)
 
template<typename T , typename Perm , typename = std::enable_if_t<TiledArray::detail::is_permutation_v<Perm>>>
DecomposedTensor< T > mult (DecomposedTensor< T > const &l, DecomposedTensor< T > const &r, Perm const &p, const T factor)
 
template<typename T >
DecomposedTensor< T > mult (DecomposedTensor< T > const &l, DecomposedTensor< T > const &r, const T factor)
 
template<typename T >
DecomposedTensor< T > & mult_to (DecomposedTensor< T > &l, DecomposedTensor< T > const &r)
 
template<typename T >
DecomposedTensor< T > & mult_to (DecomposedTensor< T > &l, DecomposedTensor< T > const &r, const T factor)
 
template<typename T >
DecomposedTensor< T > subt (DecomposedTensor< T > const &l, DecomposedTensor< T > const &r)
 
template<typename T , typename Perm , typename = std::enable_if_t<TiledArray::detail::is_permutation_v<Perm>>>
DecomposedTensor< T > subt (DecomposedTensor< T > const &l, DecomposedTensor< T > const &r, Perm const &p)
 
template<typename T >
DecomposedTensor< T > subt (DecomposedTensor< T > const &l, DecomposedTensor< T > const &r, const T factor)
 
template<typename T , typename Perm , typename = std::enable_if_t<TiledArray::detail::is_permutation_v<Perm>>>
DecomposedTensor< T > subt (DecomposedTensor< T > const &l, DecomposedTensor< T > const &r, const T factor, Perm const &p)
 
template<typename T >
DecomposedTensor< T > subt (DecomposedTensor< T > const &l, const T factor)
 
template<typename T , typename Perm , typename = std::enable_if_t<TiledArray::detail::is_permutation_v<Perm>>>
DecomposedTensor< T > subt (DecomposedTensor< T > const &l, const T factor, Perm const &p)
 
template<typename T >
DecomposedTensor< T > & subt_to (DecomposedTensor< T > &l, DecomposedTensor< T > const &r)
 
template<typename T , typename F >
DecomposedTensor< T > & subt_to (DecomposedTensor< T > &l, DecomposedTensor< T > const &r, const F factor)
 
template<typename T >
DecomposedTensor< T > & subt_to (DecomposedTensor< T > &l, const T factor)
 
template<typename T >
trace (DecomposedTensor< T > const &t)
 Returns the trace of a decompose tensor Currently not recommended due to implementation details. More...
 
template<typename T >
sum (DecomposedTensor< T > const &t)
 Returns the sum of a decomposed tensor tile Currently not recommended due to implementation details. More...
 
template<typename T >
min (DecomposedTensor< T > const &t)
 Returns the min of a decomposed tensor tile Currently not recommended due to implementation details. More...
 
template<typename T >
max (DecomposedTensor< T > const &t)
 Returns the max of a decomposed tensor tile Currently not recommended due to implementation details. More...
 
template<typename T >
abs_min (DecomposedTensor< T > const &t)
 Returns the abs_min of a decomposed tensor tile Currently not recommended due to implementation details. More...
 
template<typename T >
abs_max (DecomposedTensor< T > const &t)
 Returns the abs_max of a decomposed tensor tile Currently not recommended due to implementation details. More...
 
template<typename T >
product (DecomposedTensor< T > const &t)
 Returns the product of a decomposed tensor tile Currently not recommended due to implementation details. More...
 
template<typename T >
norm (DecomposedTensor< T > const &t)
 
template<typename T , typename R >
void norm (DecomposedTensor< T > const &t, R &result)
 
template<typename T >
squared_norm (DecomposedTensor< T > const &t)
 
template<typename T , typename Perm , typename = std::enable_if_t<TiledArray::detail::is_permutation_v<Perm>>>
DecomposedTensor< T > permute (DecomposedTensor< T > const &t, Perm const &p)
 
template<typename T >
DecomposedTensor< T > clone (DecomposedTensor< T > const &t)
 
template<typename T , typename F >
DecomposedTensor< T > scale (DecomposedTensor< T > const &t, F factor)
 
template<typename T , typename F , typename Perm , typename = std::enable_if_t<TiledArray::detail::is_permutation_v<Perm>>>
DecomposedTensor< T > scale (DecomposedTensor< T > const &t, F factor, Perm const &p)
 
template<typename T >
DecomposedTensor< T > neg (DecomposedTensor< T > const &t)
 
template<typename T , typename Perm , typename = std::enable_if_t<TiledArray::detail::is_permutation_v<Perm>>>
DecomposedTensor< T > neg (DecomposedTensor< T > const &t, Perm const &p)
 
template<typename T , typename Perm , typename = std::enable_if_t<TiledArray::detail::is_permutation_v<Perm>>>
DecomposedTensor< T > & neg_to (DecomposedTensor< T > &t, Perm const &p)
 
template<typename T >
DecomposedTensor< T > & neg_to (DecomposedTensor< T > &t)
 
template<typename T , typename F >
DecomposedTensor< T > & scale_to (DecomposedTensor< T > &t, F factor)
 
template<typename T , typename F , typename Perm , typename = std::enable_if_t<TiledArray::detail::is_permutation_v<Perm>>>
DecomposedTensor< T > & scale_to (DecomposedTensor< T > &t, F factor, Perm const &p)
 
template<typename T >
DecomposedTensor< T > compress (DecomposedTensor< T > const &t, double cut)
 
template<typename T >
bool empty (DecomposedTensor< T > const &t)
 
template<typename T , typename Op >
DecomposedTensor< T > unary (DecomposedTensor< T > const &l, Op &&op)
 
template<typename T , typename Op , typename Perm , typename = std::enable_if_t<TiledArray::detail::is_permutation_v<Perm>>>
DecomposedTensor< T > binary (DecomposedTensor< T > const &l, Op &&op, Perm const &p)
 
template<typename T , typename Op >
DecomposedTensor< T > & inplace_unary (DecomposedTensor< T > &l, Op &&op)
 
template<typename T >
std::ostream & operator<< (std::ostream &os, Tile< T > const &tile)
 
template<typename T >
unsigned long tile_clr_storage (Tile< DecomposedTensor< T >> const &tile)
 
template<typename T >
void attach_to_closest (Matrix< T > const &D, std::vector< VectorCluster< T >> &clusters, bool init=false)
 
template<typename T >
std::vector< VectorCluster< T > > init_rows (Matrix< T > const &D, unsigned long num_clusters)
 
template<typename T >
Vector< unsigned long > get_pivots (std::vector< VectorCluster< T >> const &clusters)
 
template<typename T >
auto all_have_elems (std::vector< T > &vs)
 
template<typename T >
TA::TiledRange1 localize_vectors_with_kmeans (Matrix< T > const &xyz, Matrix< T > &D, unsigned long num_clusters)
 
template<typename Array , class conv_class >
Array cp_als_rank (const Array &reference, std::vector< Array > &factor_matrices, std::vector< size_t > &symmetries, conv_class &conv, int block_size=1, int rank=0, double max_ALS=5e3, bool fast_pI=true, bool direct=true, bool recompose=false)
 
template<typename Array , class conv_class >
Array cp_rals_rank (const Array &reference, std::vector< Array > &factor_matrices, std::vector< size_t > symmetries, conv_class &conv, int block_size=1, int rank=0, double max_ALS=5e3, bool fast_pI=true, bool direct=true, bool recompose=false)
 
template<typename Array , class conv_class >
void cp_df_als_rank (const Array &referenceL, const Array &referenceR, const Array &ref_full, std::vector< Array > &factor_matrices, conv_class &conv, std::vector< size_t > &symm, int block_size=1, int rank=0, double max_ALS=5e3, bool fast_pI=true)
 
template<typename Array >
void coupled_cp_als (const Array &referenceL, const Array &referenceR, std::vector< Array > &factor_matrices, int block_size=1, int rank=0, double max_ALS=5e3, double tcutALS=5e-3, int SVD_initial_guess=false, int SVD_rank=0, int step=1)
 
template<typename Array >
void cp_exchange_formation (const Array &referenceL, const Array &referenceR, std::vector< Array > &factor_matrices, int block_size=1, int rank=0, double max_ALS=5e3, double tcutALS=1e-2, bool fast_pI=true)
 
template<typename Array >
void cp_decomp_for_LT_T (const Array &Xab, const Array &Xai, std::vector< Array > &factor_matrices, int block_size=1, int rank=0, double max_ALS=5e3, double tcutALS=1e-2, bool fast_pI=true)
 
template<typename Tile , typename Policy >
std::tuple< std::vector< TiledArray::DistArray< Tile, Policy > >, std::vector< TiledArray::DistArray< Tile, Policy > >, std::vector< TiledArray::DistArray< Tile, Policy > > > compute_cp3_df_int_decomp (const TA::DistArray< Tile, Policy > &Xab, double rank_block_factor, int unocc_block_size, size_t cp_rank, double cp_precision, const TA::DistArray< Tile, Policy > &Had_reblock_matrix=TA::DistArray< Tile, Policy >(), bool cp_ps=false, bool cp_robust=false, bool verbose=false, bool ta_cp3=false, const TA::DistArray< Tile, Policy > &Xai=TA::DistArray< Tile, Policy >())
 
template<typename integral , typename Tile , typename Policy >
void generate_cc_cp3_decomp (integral &ints, int had_block_size_, double cp3_rank_, bool reduced_abcd_memory_, double rank_block_factor_, int vir_block_size, double cp3_precision_, bool ta_als_)
 
template<typename integral , typename TArray , typename Tile , typename Policy >
void generate_cc_cp4_decomp (integral &ints, TArray Xma, TArray Cma, TArray Cmi, int had_block_size_, double cp3_rank_, double cp4_rank_, double rank_block_factor_, int vir_block_size, double cp3_precision_, double cp4_precision_)
 

Typedef Documentation

◆ Matrix

template<typename T >
using mpqc::math::Matrix = typedef Eigen::Matrix<T, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>

◆ Vector

template<typename T >
using mpqc::math::Vector = typedef Eigen::Matrix<T, Eigen::Dynamic, 1>

Function Documentation

◆ abs_max()

template<typename T >
T mpqc::math::abs_max ( DecomposedTensor< T > const &  t)

Returns the abs_max of a decomposed tensor tile Currently not recommended due to implementation details.

◆ abs_min()

template<typename T >
T mpqc::math::abs_min ( DecomposedTensor< T > const &  t)

Returns the abs_min of a decomposed tensor tile Currently not recommended due to implementation details.

◆ add() [1/6]

template<typename T >
DecomposedTensor<T> mpqc::math::add ( DecomposedTensor< T > const &  l,
const T  factor 
)

◆ add() [2/6]

template<typename T , typename Perm , typename = std::enable_if_t<TiledArray::detail::is_permutation_v<Perm>>>
DecomposedTensor<T> mpqc::math::add ( DecomposedTensor< T > const &  l,
const T  factor,
Perm const &  p 
)

◆ add() [3/6]

template<typename T >
DecomposedTensor<T> mpqc::math::add ( DecomposedTensor< T > const &  l,
DecomposedTensor< T > const &  r 
)

◆ add() [4/6]

template<typename T >
DecomposedTensor<T> mpqc::math::add ( DecomposedTensor< T > const &  l,
DecomposedTensor< T > const &  r,
const T  factor 
)

◆ add() [5/6]

template<typename T , typename Perm , typename = std::enable_if_t<TiledArray::detail::is_permutation_v<Perm>>>
DecomposedTensor<T> mpqc::math::add ( DecomposedTensor< T > const &  l,
DecomposedTensor< T > const &  r,
const T  factor,
Perm const &  p 
)

◆ add() [6/6]

template<typename T , typename Perm , typename = std::enable_if_t<TiledArray::detail::is_permutation_v<Perm>>>
DecomposedTensor<T> mpqc::math::add ( DecomposedTensor< T > const &  l,
DecomposedTensor< T > const &  r,
Perm const &  p 
)

◆ add_order2()

template<typename T >
DecomposedTensor<T> mpqc::math::add_order2 ( DecomposedTensor< T > const &  l,
DecomposedTensor< T > const &  r,
factor 
)

◆ add_to() [1/3]

template<typename T >
DecomposedTensor<T>& mpqc::math::add_to ( DecomposedTensor< T > &  l,
const T  factor 
)

◆ add_to() [2/3]

template<typename T >
DecomposedTensor<T>& mpqc::math::add_to ( DecomposedTensor< T > &  l,
DecomposedTensor< T > const &  r 
)

◆ add_to() [3/3]

template<typename T >
DecomposedTensor<T>& mpqc::math::add_to ( DecomposedTensor< T > &  l,
DecomposedTensor< T > const &  r,
const T  factor 
)

◆ add_to_diag()

template<typename Array >
void mpqc::math::add_to_diag ( Array A,
double  val 
)

◆ add_to_diag_tile() [1/2]

void mpqc::math::add_to_diag_tile ( double  val,
TA::Tile< DecomposedTensor< double >> &  tile 
)
inline

◆ add_to_diag_tile() [2/2]

template<typename Scalar , typename Tensor , std::enable_if_t< TA::detail::is_contiguous_tensor_v< Tensor >> * = nullptr>
void mpqc::math::add_to_diag_tile ( Scalar  val,
Tensor &  tile 
)
inline

◆ add_to_order2()

template<typename T >
DecomposedTensor<T>& mpqc::math::add_to_order2 ( DecomposedTensor< T > &  l,
DecomposedTensor< T > const &  r 
)

◆ all_have_elems()

template<typename T >
auto mpqc::math::all_have_elems ( std::vector< T > &  vs)

◆ array_to_eigen() [1/2]

template<typename Tile , typename Policy >
Matrix<typename Tile::value_type> mpqc::math::array_to_eigen ( TA::DistArray< Tile, Policy > const &  A)

converts a TiledArray::Array to a row-major Eigen Matrix

Parameters
[in]Aan order-2 tensor
Warning
this is a collective operation over the World object returned by A.world()

◆ array_to_eigen() [2/2]

template<typename Tile , typename Policy >
std::vector<Matrix<typename Tile::value_type> > mpqc::math::array_to_eigen ( TA::DistArrayVector< Tile, Policy > const &  A)

converts a TiledArray::ArrayVector to a row-major Eigen Matrix

Parameters
[in]Aan order-2 tensor
Warning
this is a collective operation over the World object returned by A.world()

◆ attach_to_closest()

template<typename T >
void mpqc::math::attach_to_closest ( Matrix< T > const &  D,
std::vector< VectorCluster< T >> &  clusters,
bool  init = false 
)

◆ binary() [1/3]

template<typename T , typename Op >
DecomposedTensor<T> mpqc::math::binary ( DecomposedTensor< T > const &  l,
DecomposedTensor< T > const &  r,
Op &&  op 
)

◆ binary() [2/3]

template<typename T , typename Op , typename Perm , typename = std::enable_if_t<TiledArray::detail::is_permutation_v<Perm>>>
DecomposedTensor<T> mpqc::math::binary ( DecomposedTensor< T > const &  l,
DecomposedTensor< T > const &  r,
Op &&  op,
Perm const &  p 
)

◆ binary() [3/3]

template<typename T , typename Op , typename Perm , typename = std::enable_if_t<TiledArray::detail::is_permutation_v<Perm>>>
DecomposedTensor<T> mpqc::math::binary ( DecomposedTensor< T > const &  l,
Op &&  op,
Perm const &  p 
)

◆ cholesky()

template<typename Tile , typename Policy >
TA::DistArray<Tile, Policy> mpqc::math::cholesky ( TA::DistArray< Tile, Policy > const &  A)

◆ cholesky_inverse()

template<typename Tile , typename Policy >
TA::DistArray<Tile, Policy> mpqc::math::cholesky_inverse ( TA::DistArray< Tile, Policy > const &  A)

◆ clone()

template<typename T >
DecomposedTensor<T> mpqc::math::clone ( DecomposedTensor< T > const &  t)

◆ combine()

TiledArray::Tensor< double > mpqc::math::combine ( DecomposedTensor< double > const &  t)

◆ complex_random_unitary()

template<typename T >
RowMatrix<T> mpqc::math::complex_random_unitary ( int  size)

◆ compress()

template<typename T >
DecomposedTensor<T> mpqc::math::compress ( DecomposedTensor< T > const &  t,
double  cut 
)

◆ compute_cp3_df_int_decomp()

template<typename Tile , typename Policy >
std::tuple<std::vector<TiledArray::DistArray<Tile, Policy> >, std::vector<TiledArray::DistArray<Tile, Policy> >, std::vector<TiledArray::DistArray<Tile, Policy> > > mpqc::math::compute_cp3_df_int_decomp ( const TA::DistArray< Tile, Policy > &  Xab,
double  rank_block_factor,
int  unocc_block_size,
size_t  cp_rank,
double  cp_precision,
const TA::DistArray< Tile, Policy > &  Had_reblock_matrix = TA::DistArray<Tile, Policy>(),
bool  cp_ps = false,
bool  cp_robust = false,
bool  verbose = false,
bool  ta_cp3 = false,
const TA::DistArray< Tile, Policy > &  Xai = TA::DistArray<Tile, Policy>() 
)

Function to compute the CP3 decomposition of DF integrals in any method can use PS ( $ Xab_{cp} Xab $), DF ( $ Xab_{cp} Xab_{cp}$), or rCP (robust) ( $ Xab_{cp} (2.0 Xab - Xab_{cp} $)

Returns
a tuple first the CP factor matrices for a Xab tensor, then tensor (X * T) where X is the factor matrix of the DF mode and T is the correct Xab depending on if using PS, DF or rCP and, finally, a contraction of (X * Xai) if Xai is provided.
Parameters
[in]XabThe reference DF like tensor being decomposed.
[in]rank_block_factorblocking info for the rank mode block size is determined by unnoc_block_size * rank_block_factor
[in]cp_rankrank of the cp decomposition
[in]cp_precisionALS precision for the CP decomposition
[in]Had_reblock_matrixDistArray<Tile, Policy>() A matrix to reblock the unnoc dimensions for the mixed contraction in CC methods.
[in]cp_psfalse Using the CP PS approach?
[in]cp_robustfalse Using the rCP-DF approach?
[in]verbosefalse Print CP information
[in]ta_cp3false Currently unavailable using TA to compute CP decomposition
[in]XaiDistArray<Tile, Policy>() A second DF tensor, used in a full CP3 + CP4 implementation of CC methods.

◆ conditioned_orthogonalizer()

template<typename Tile , typename Policy >
TA::DistArray<Tile, Policy> mpqc::math::conditioned_orthogonalizer ( TA::DistArray< Tile, Policy >  S_array,
double  S_condition_number_threshold 
)

◆ coupled_cp_als()

template<typename Array >
void mpqc::math::coupled_cp_als ( const Array referenceL,
const Array referenceR,
std::vector< Array > &  factor_matrices,
int  block_size = 1,
int  rank = 0,
double  max_ALS = 5e3,
double  tcutALS = 5e-3,
int  SVD_initial_guess = false,
int  SVD_rank = 0,
int  step = 1 
)

REQUIRES BTAS Calculates the canonical product (CP) decomposition of two reference tensors where the first factor of both tensors are eqiuvalent and minimize both loss functions at a specific rank using alternating least squares (ALS) method in BTAS

Parameters
[in]referenceLIn: First reference tensor to be decomposed.
[in]referenceRIn: Second reference tensor to be decomposed.
[in,out]factor_matricesIn: empty vector. Out: Set of $ N_{left} + *N_{right} - 1 $ CP-factor matrices.
[in]convconvergence test used to determine if ALS optmization is finished
[in]symmA list which tells solver which modes are equivalent should be ordered such that later modes equal earlier modes. If empty all modes will be optimized independently
[in]block_size1 Block size of the rank dimension for the factor matrices
[in]rank0 what is the CP decomposition rank
[in]step1 if SVD_rank < rank how many columns to build rank per optimization.

◆ cp_als_rank()

template<typename Array , class conv_class >
Array mpqc::math::cp_als_rank ( const Array reference,
std::vector< Array > &  factor_matrices,
std::vector< size_t > &  symmetries,
conv_class &  conv,
int  block_size = 1,
int  rank = 0,
double  max_ALS = 5e3,
bool  fast_pI = true,
bool  direct = true,
bool  recompose = false 
)

REQUIRES BTAS Calculates the canonical product (CP) decomposition of /c reference to a specific rank using alternating least squares (ALS) method in BTAS

Parameters
[in]referenceIn: tensor to be decomposed.
[in,out]factor_matricesIn: empty vector. Out: Set of N CP-factor matrices.
[in]symmetriesA list which tells solver which modes are equivalent should be ordered such that later modes equal earlier modes. If empty all modes will be optimized independently
[in]convconvergence test used to determine if ALS optmization is finished
[in]block_size1 Block size of the rank dimension for the factor matrices
[in]rank0 what is the CP decomposition rank
[in]max_als1e3 max number of iterations to optimize rank r factors using ALS
[in]fast_pItrue Should ALS use a fast pseudo inverse scheme
[in]directtrue Should the method compute the ALS without the Khatri-Rao product
[in]recomposefalse Return the CP approximation of tensor reference
Returns
CP reference approximate tensor.

◆ cp_decomp_for_LT_T()

template<typename Array >
void mpqc::math::cp_decomp_for_LT_T ( const Array Xab,
const Array Xai,
std::vector< Array > &  factor_matrices,
int  block_size = 1,
int  rank = 0,
double  max_ALS = 5e3,
double  tcutALS = 1e-2,
bool  fast_pI = true 
)

REQUIRES BTAS Specialized version of cp_df_als_rank for LT (T)

Parameters
[in]XabIn: Left side of reference tensor product to be decomposed.
[in]XaiIn: Right side of reference tensor product to be decomposed.
[in,out]factor_matricesIn: empty vector. Out: Set of N CP-factor matrices.
[in]block_size1 Block size of the rank dimension for the factor matrices
[in]rank0 what is the CP decomposition rank
[in]max_als5e3 max number of iterations to optimize rank r factors using ALS
[in]tcutALS1e-2 ALS stopping parameter using BTAS FitCheck
[in]fast_pItrue Should ALS use a fast pseudo inverse scheme

◆ cp_df_als_rank()

template<typename Array , class conv_class >
void mpqc::math::cp_df_als_rank ( const Array referenceL,
const Array referenceR,
const Array ref_full,
std::vector< Array > &  factor_matrices,
conv_class &  conv,
std::vector< size_t > &  symm,
int  block_size = 1,
int  rank = 0,
double  max_ALS = 5e3,
bool  fast_pI = true 
)

REQUIRES BTAS Calculates the canonical product (CP) decomposition of a reference ( created via a tensor contraction of referenceL and referenceLR over the 2 references first mode) ata specific rank using alternating least squares (ALS) method in BTAS

Parameters
[in]referenceLIn: Left side of reference tensor product to be decomposed.
[in]referenceRIn: Right side of reference tensor product to be decomposed.
[in,out]factor_matricesIn: empty vector. Out: Set of N CP-factor matrices.
[in]convconvergence test used to determine if ALS optmization is finished
[in]symmA list which tells solver which modes are equivalent should be ordered such that later modes equal earlier modes. If empty all modes will be optimized independently
[in]block_size1 Block size of the rank dimension for the factor matrices
[in]rank0 what is the CP decomposition rank
[in]max_als1e3 max number of iterations to optimize rank r factors using ALS
[in]fast_pItrue Should ALS use a fast pseudo inverse scheme

◆ cp_exchange_formation()

template<typename Array >
void mpqc::math::cp_exchange_formation ( const Array referenceL,
const Array referenceR,
std::vector< Array > &  factor_matrices,
int  block_size = 1,
int  rank = 0,
double  max_ALS = 5e3,
double  tcutALS = 1e-2,
bool  fast_pI = true 
)

REQUIRES BTAS Specialized version of cp_df_als_rank Assumes one is using 2 order-3 reference tensors after decomposition forms order-3 factors in the order (referenceL $ \in \mathbf{R}^{X *\times a \times i} $ referenceR $ \in \mathbf{R}^{X \times b \times j} *$) $ f_1 \in \mathbf{R}^{a\times b \times r}$ $ f_2 \in *\mathbf{R}^{a\times j \times r}$ $ f_3 \in \mathbf{R}^{i\times b \times *r}$ $ f_4 \in \mathbf{R}^{i\times j \times r}$

Parameters
[in]referenceLIn: Left side of reference tensor product to be decomposed.
[in]referenceRIn: Right side of reference tensor product to be decomposed.
[in,out]factor_matricesIn: empty vector. Out: Set of N CP-factor matrices.
[in]block_size1 Block size of the rank dimension for the factor matrices
[in]rank0 what is the CP decomposition rank
[in]max_als5e3 max number of iterations to optimize rank r factors using ALS
[in]tcutALS1e-2 ALS stopping parameter using BTAS FitCheck
[in]fast_pItrue Should ALS use a fast pseudo inverse scheme

◆ cp_rals_rank()

template<typename Array , class conv_class >
Array mpqc::math::cp_rals_rank ( const Array reference,
std::vector< Array > &  factor_matrices,
std::vector< size_t >  symmetries,
conv_class &  conv,
int  block_size = 1,
int  rank = 0,
double  max_ALS = 5e3,
bool  fast_pI = true,
bool  direct = true,
bool  recompose = false 
)

REQUIRES BTAS Calculates the canonical product (CP) decomposition of /c reference to a specific rank using a regularized alternating least squares (RALS) method in BTAS

Parameters
[in]referenceIn: tensor to be decomposed.
[in,out]factor_matricesIn: empty vector. Out: Set of N CP-factor matrices.
[in]symmetriesA list which tells solver which modes are equivalent should be ordered such that later modes equal earlier modes. If empty all modes will be optimized independently
[in]convconvergence test used to determine if ALS optmization is finished
[in]block_size1 Block size of the rank dimension for the factor matrices
[in]rank0 what is the CP decomposition rank
[in]max_als1e3 max number of iterations to optimize rank r factors using ALS
[in]fast_pItrue Should ALS use a fast pseudo inverse scheme
[in]directtrue Should the method compute the ALS without the Khatri-Rao product
[in]recomposefalse Return the CP approximation of tensor reference
Returns
CP reference approximate tensor.

◆ create_diagonal_array_from_eigen()

template<typename Tile , typename Policy >
TA::DistArray<Tile, Policy> mpqc::math::create_diagonal_array_from_eigen ( madness::World &  world,
const TA::TiledRange1 &  trange1,
const TA::TiledRange1 &  trange2,
typename Tile::numeric_type  val 
)

◆ create_diagonal_matrix() [1/3]

template<typename T , unsigned int N, typename Tile >
TiledArray::Array<T, N, Tile, TiledArray::DensePolicy> mpqc::math::create_diagonal_matrix ( TiledArray::Array< T, N, Tile, TiledArray::DensePolicy > const &  model,
double  val 
)

◆ create_diagonal_matrix() [2/3]

template<typename Tile , typename Policy >
TiledArray::DistArray< Tile, typename std::enable_if<std::is_same<Policy, TA::SparsePolicy>::value, TA::SparsePolicy>::type> mpqc::math::create_diagonal_matrix ( TiledArray::DistArray< Tile, Policy > const &  model,
double  val 
)

create sparse diagonal TA::DistArray matrix

Template Parameters
Tile
PolicyTA::SparsePolicy
Parameters
modelmodel used to construct result array
valvalue
Returns

◆ create_diagonal_matrix() [3/3]

template<typename Tile , typename Policy >
TiledArray::DistArray< Tile, typename std::enable_if<std::is_same<Policy, TA::DensePolicy>::value, TA::DensePolicy>::type> mpqc::math::create_diagonal_matrix ( TiledArray::DistArray< Tile, Policy > const &  model,
double  val 
)

create sparse diagonal TA::DistArray matrix

Template Parameters
Tile
PolicyTA::SparsePolicy
Parameters
modelmodel used to construct result array
valvalue
Returns

◆ diagonal_matrix()

template<typename Tile >
TiledArray::DistArray<Tile, TiledArray::SparsePolicy> mpqc::math::diagonal_matrix ( TiledArray::TiledRange const &  trange,
double  val,
madness::World &  world 
)

takes a TiledArray::TiledRange and a value and returns a diagonal matrix.

This function creates a diagonal matrix given a TiledArray::TiledRange and a value. The template parameters must be a numeric type followed by a tile type. Finally there must be a create_diagonal_tile overload for the tile type that gets passed in.

By default this function will create the array using the default madness::World, but you can pass in a world as an optional third parameter.

Todo:
Finish diagonal matrix this is a little trickier than previous identity functions because it needs to gracefully handle non symmetric tiling.

◆ eigen_estimator()

template<typename Tile , std::enable_if_t< TiledArray::detail::is_contiguous_tensor_v< Tile >> * = nullptr>
void mpqc::math::eigen_estimator ( std::array< Tile, 2 > &  result,
Tile const &  tile 
)

◆ eigen_inverse()

template<typename Tile , typename Policy >
TA::DistArray<Tile, Policy> mpqc::math::eigen_inverse ( const TA::DistArray< Tile, Policy > &  A)

◆ eigen_to_array() [1/2]

template<typename Tile , typename Policy >
TA::DistArray< Tile, typename std::enable_if<std::is_same<Policy, TA::SparsePolicy>::value, TA::SparsePolicy>::type> mpqc::math::eigen_to_array ( madness::World &  world,
Matrix< typename Tile::numeric_type > const &  M,
TA::TiledRange1  tr0,
TA::TiledRange1  tr1,
double  cut = 1e-7 
)

convert eigen(replicated) to sparse TA::DistArray

Template Parameters
Tile
Policy
Parameters
world
MEigen matrix, must be replicated on all nodes
tr0the TiledArray::TiledRange1 object for the row dimension of the result
tr1the TiledArray::TiledRange1 object for the column dimension of the result
cut
Returns

◆ eigen_to_array() [2/2]

template<typename Tile , typename Policy >
TA::DistArray< Tile, typename std::enable_if<std::is_same<Policy, TA::DensePolicy>::value, TA::DensePolicy>::type> mpqc::math::eigen_to_array ( madness::World &  world,
Matrix< typename Tile::numeric_type > const &  M,
TA::TiledRange1  tr0,
TA::TiledRange1  tr1,
double  cut = 1e-7 
)

convert eigen(replicated) to dense TA::DistArray

Template Parameters
Tile
PolicyTA::DensePolicy
Parameters
world
MEigen matrix, must be replicated on all nodes
tr0the TiledArray::TiledRange1 object for the row dimension of the result
tr1the TiledArray::TiledRange1 object for the column dimension of the result
cut
Returns

◆ empty()

template<typename T >
bool mpqc::math::empty ( DecomposedTensor< T > const &  t)

◆ eval_guess()

template<typename Array >
std::array<typename Array::value_type::numeric_type, 2> mpqc::math::eval_guess ( Array const &  A)

computes {min,max} eigenvalue bounds for a square matrix A using Gerschgorin circle theorem

Note
https://en.wikipedia.org/wiki/Gershgorin_circle_theorem
Template Parameters
Arraya DistArray type
Parameters
[in]Aa square matrix
Returns
{min,max} eigenavlue bounds

◆ gauss_legendre()

void mpqc::math::gauss_legendre ( int  N,
Eigen::VectorXd &  w,
Eigen::VectorXd &  x 
)

this function uses orthogonal polynomial approach for obtaining quadrature roots and weights. Taken from paper: Gaussian Quadrature and the Eigenvalue Problem by John A. Gubner. http://gubner.ece.wisc.edu/gaussquad.pdf The code is outlined at Example 15: Legendre polynomials

Parameters
N[in]is the number of quadrature points
w[out]will return the weights
x[out]will return the roots

◆ gemm() [1/4]

template<typename T >
DecomposedTensor<T>& mpqc::math::gemm ( DecomposedTensor< T > &  c,
DecomposedTensor< T > const &  a,
DecomposedTensor< T > const &  b,
const T  factor,
TA::math::GemmHelper const &  gh 
)

◆ gemm() [2/4]

template<typename T >
DecomposedTensor<T> mpqc::math::gemm ( DecomposedTensor< T > const &  a,
DecomposedTensor< T > const &  b,
const T  factor,
TA::math::GemmHelper const &  gh 
)

◆ gemm() [3/4]

template<typename T >
TA::TensorD mpqc::math::gemm ( DecomposedTensor< T > const &  a,
TA::TensorD const &  b,
const T  factor,
TA::math::GemmHelper const &  gh 
)

◆ gemm() [4/4]

template<typename T >
TA::TensorD mpqc::math::gemm ( TA::TensorD &  c,
DecomposedTensor< T > const &  a,
TA::TensorD const &  b,
const T  factor,
TA::math::GemmHelper const &  gh 
)

◆ generate_cc_cp3_decomp()

template<typename integral , typename Tile , typename Policy >
void mpqc::math::generate_cc_cp3_decomp ( integral &  ints,
int  had_block_size_,
double  cp3_rank_,
bool  reduced_abcd_memory_,
double  rank_block_factor_,
int  vir_block_size,
double  cp3_precision_,
bool  ta_als_ 
)

Generates all necessary integrals for cc cp3 currently used in CCSD and CCSDT

Parameters
[in,out]intsIn: CC integrals will access Xab and possibly Xai Out: Will compute X * T (as described above) and will compute cp factor matrices
[in]had_block_size_size of the unocc blocking for mixed contraction
[in]cp3_rank_rank of the CP decomposition of Xab
[in]reduced_abcd_memory_If reduced_abcd_memory_ = false will also need Xai contraction
[in]rank_block_factor_Used to compute the block size for rank mode rank block = rank_block_factor_ * vir_block_size
[in]cp3_precision_ALS optimization stopping condition parameter.
[in]ta_als_using TA to compute CP decomposition ** currently disabled **

◆ generate_cc_cp4_decomp()

template<typename integral , typename TArray , typename Tile , typename Policy >
void mpqc::math::generate_cc_cp4_decomp ( integral &  ints,
TArray  Xma,
TArray  Cma,
TArray  Cmi,
int  had_block_size_,
double  cp3_rank_,
double  cp4_rank_,
double  rank_block_factor_,
int  vir_block_size,
double  cp3_precision_,
double  cp4_precision_ 
)

Generates all necessary integrals for cc cp3 + cp4 implementation currently used in CCSD and CCSDT approximation computes the CP3 decomposition of Xab and the CP4 decomposition of $ G^{a \mu}_{ib} $

Parameters
[in,out]intsIn: CC integrals will access Xab Out: Will compute X * T (as described above) and will compute cp factor matrices
[in]had_block_size_size of the unocc blocking for mixed contraction
[in]cp3_rank_rank of the CP decomposition of Xab
[in]reduced_abcd_memory_If reduced_abcd_memory_ = false will also need Xai contraction
[in]rank_block_factor_Used to compute the block size for rank mode rank block = rank_block_factor_ * vir_block_size
[in]cp3_precision_ALS optimization stopping condition parameter.
[in]ta_als_using TA to compute CP decomposition ** currently disabled **

◆ gensqrtinv()

template<typename T >
std::tuple<RowMatrix<T>, RowMatrix<T>, size_t, double, double> mpqc::math::gensqrtinv ( const RowMatrix< T > &  S,
bool  symmetric = false,
double  max_condition_number = 1e8 
)

◆ get_pivots()

template<typename T >
Vector<unsigned long> mpqc::math::get_pivots ( std::vector< VectorCluster< T >> const &  clusters)

◆ gram_schmidt() [1/2]

template<typename D >
void mpqc::math::gram_schmidt ( const std::vector< D > &  V1,
std::vector< D > &  V2,
double  threshold 
)

vector V1 is already orthonormalized orthonormalize V2 with respect to V1

◆ gram_schmidt() [2/2]

template<typename D >
void mpqc::math::gram_schmidt ( std::vector< D > &  V,
double  threshold,
std::size_t  start = 0 
)

◆ init_rows()

template<typename T >
std::vector<VectorCluster<T> > mpqc::math::init_rows ( Matrix< T > const &  D,
unsigned long  num_clusters 
)

◆ init_X()

template<typename Array >
Array mpqc::math::init_X ( Array const &  S)

◆ inplace_binary()

template<typename T , typename Op >
DecomposedTensor<T>& mpqc::math::inplace_binary ( DecomposedTensor< T > &  l,
DecomposedTensor< T > const &  r,
Op &&  op 
)

◆ inplace_unary()

template<typename T , typename Op >
DecomposedTensor<T>& mpqc::math::inplace_unary ( DecomposedTensor< T > &  l,
Op &&  op 
)

◆ inverse_sqrt()

template<typename Array >
Array mpqc::math::inverse_sqrt ( Array const &  S)

◆ invert()

template<typename Array >
Array mpqc::math::invert ( Array const &  S)

◆ localize_vectors_with_kmeans()

template<typename T >
TA::TiledRange1 mpqc::math::localize_vectors_with_kmeans ( Matrix< T > const &  xyz,
Matrix< T > &  D,
unsigned long  num_clusters 
)

◆ make_diagonal_tile() [1/2]

template<typename T >
void mpqc::math::make_diagonal_tile ( TA::Tile< DecomposedTensor< T >> &  tile,
val 
)

◆ make_diagonal_tile() [2/2]

template<typename Tile , typename Scalar , std::enable_if_t< TA::detail::is_scalar_v< Scalar >> * = nullptr>
void mpqc::math::make_diagonal_tile ( Tile tile,
Scalar  val 
)

◆ max()

template<typename T >
T mpqc::math::max ( DecomposedTensor< T > const &  t)

Returns the max of a decomposed tensor tile Currently not recommended due to implementation details.

◆ max_eval_est()

template<typename Array >
double mpqc::math::max_eval_est ( Array const &  S)

◆ min()

template<typename T >
T mpqc::math::min ( DecomposedTensor< T > const &  t)

Returns the min of a decomposed tensor tile Currently not recommended due to implementation details.

◆ min_eval_est()

template<typename Array >
double mpqc::math::min_eval_est ( Array const &  H,
Array const &  S 
)

◆ min_eval_guess()

template<typename It >
double mpqc::math::min_eval_guess ( It  first,
It  second 
)

◆ mult() [1/4]

template<typename T >
DecomposedTensor<T> mpqc::math::mult ( DecomposedTensor< T > const &  l,
DecomposedTensor< T > const &  r 
)

◆ mult() [2/4]

template<typename T >
DecomposedTensor<T> mpqc::math::mult ( DecomposedTensor< T > const &  l,
DecomposedTensor< T > const &  r,
const T  factor 
)

◆ mult() [3/4]

template<typename T , typename Perm , typename = std::enable_if_t<TiledArray::detail::is_permutation_v<Perm>>>
DecomposedTensor<T> mpqc::math::mult ( DecomposedTensor< T > const &  l,
DecomposedTensor< T > const &  r,
Perm const &  p 
)

◆ mult() [4/4]

template<typename T , typename Perm , typename = std::enable_if_t<TiledArray::detail::is_permutation_v<Perm>>>
DecomposedTensor<T> mpqc::math::mult ( DecomposedTensor< T > const &  l,
DecomposedTensor< T > const &  r,
Perm const &  p,
const T  factor 
)

◆ mult_to() [1/2]

template<typename T >
DecomposedTensor<T>& mpqc::math::mult_to ( DecomposedTensor< T > &  l,
DecomposedTensor< T > const &  r 
)

◆ mult_to() [2/2]

template<typename T >
DecomposedTensor<T>& mpqc::math::mult_to ( DecomposedTensor< T > &  l,
DecomposedTensor< T > const &  r,
const T  factor 
)

◆ neg() [1/2]

template<typename T >
DecomposedTensor<T> mpqc::math::neg ( DecomposedTensor< T > const &  t)

◆ neg() [2/2]

template<typename T , typename Perm , typename = std::enable_if_t<TiledArray::detail::is_permutation_v<Perm>>>
DecomposedTensor<T> mpqc::math::neg ( DecomposedTensor< T > const &  t,
Perm const &  p 
)

◆ neg_to() [1/2]

template<typename T >
DecomposedTensor<T>& mpqc::math::neg_to ( DecomposedTensor< T > &  t)

◆ neg_to() [2/2]

template<typename T , typename Perm , typename = std::enable_if_t<TiledArray::detail::is_permutation_v<Perm>>>
DecomposedTensor<T>& mpqc::math::neg_to ( DecomposedTensor< T > &  t,
Perm const &  p 
)

◆ norm() [1/2]

template<typename T >
T mpqc::math::norm ( DecomposedTensor< T > const &  t)

◆ norm() [2/2]

template<typename T , typename R >
void mpqc::math::norm ( DecomposedTensor< T > const &  t,
R &  result 
)

◆ operator<<() [1/3]

template<typename Value >
std::ostream& mpqc::math::operator<< ( std::ostream &  os,
const TaylorExpansionCoefficients< Value > &  x 
)

◆ operator<<() [2/3]

template<typename T >
std::ostream& mpqc::math::operator<< ( std::ostream &  os,
DecomposedTensor< T > const &  t 
)

◆ operator<<() [3/3]

template<typename T >
std::ostream& mpqc::math::operator<< ( std::ostream &  os,
Tile< T > const &  tile 
)
inline

◆ permute()

template<typename T , typename Perm , typename = std::enable_if_t<TiledArray::detail::is_permutation_v<Perm>>>
DecomposedTensor<T> mpqc::math::permute ( DecomposedTensor< T > const &  t,
Perm const &  p 
)

◆ piv_cholesky()

void mpqc::math::piv_cholesky ( Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > &  a)

performs the pivoted Cholesky decomposition

This function will take a matrix that is symmetric semi-positive definate and over write it with the permuted and truncated Cholesky vectors corresponding to the matrix.

◆ product()

template<typename T >
T mpqc::math::product ( DecomposedTensor< T > const &  t)

Returns the product of a decomposed tensor tile Currently not recommended due to implementation details.

◆ pseudoinverse()

template<typename Tile , typename Policy >
TA::DistArray<Tile, Policy> mpqc::math::pseudoinverse ( TA::DistArray< Tile, Policy > const &  A,
bool  HPSD = false 
)

◆ random_unitary()

template<typename T >
RowMatrix<T> mpqc::math::random_unitary ( int  size,
int  seed = 42 
)

◆ recompress()

void mpqc::math::recompress ( DecomposedTensor< double > &  t)

Currently modifies input data regardless could cause some loss of accuracy.

◆ scale() [1/2]

template<typename T , typename F >
DecomposedTensor<T> mpqc::math::scale ( DecomposedTensor< T > const &  t,
factor 
)

◆ scale() [2/2]

template<typename T , typename F , typename Perm , typename = std::enable_if_t<TiledArray::detail::is_permutation_v<Perm>>>
DecomposedTensor<T> mpqc::math::scale ( DecomposedTensor< T > const &  t,
factor,
Perm const &  p 
)

◆ scale_to() [1/2]

template<typename T , typename F >
DecomposedTensor<T>& mpqc::math::scale_to ( DecomposedTensor< T > &  t,
factor 
)

◆ scale_to() [2/2]

template<typename T , typename F , typename Perm , typename = std::enable_if_t<TiledArray::detail::is_permutation_v<Perm>>>
DecomposedTensor<T>& mpqc::math::scale_to ( DecomposedTensor< T > &  t,
factor,
Perm const &  p 
)

◆ squared_norm()

template<typename T >
T mpqc::math::squared_norm ( DecomposedTensor< T > const &  t)

◆ subt() [1/6]

template<typename T >
DecomposedTensor<T> mpqc::math::subt ( DecomposedTensor< T > const &  l,
const T  factor 
)

◆ subt() [2/6]

template<typename T , typename Perm , typename = std::enable_if_t<TiledArray::detail::is_permutation_v<Perm>>>
DecomposedTensor<T> mpqc::math::subt ( DecomposedTensor< T > const &  l,
const T  factor,
Perm const &  p 
)

◆ subt() [3/6]

template<typename T >
DecomposedTensor<T> mpqc::math::subt ( DecomposedTensor< T > const &  l,
DecomposedTensor< T > const &  r 
)

◆ subt() [4/6]

template<typename T >
DecomposedTensor<T> mpqc::math::subt ( DecomposedTensor< T > const &  l,
DecomposedTensor< T > const &  r,
const T  factor 
)

◆ subt() [5/6]

template<typename T , typename Perm , typename = std::enable_if_t<TiledArray::detail::is_permutation_v<Perm>>>
DecomposedTensor<T> mpqc::math::subt ( DecomposedTensor< T > const &  l,
DecomposedTensor< T > const &  r,
const T  factor,
Perm const &  p 
)

◆ subt() [6/6]

template<typename T , typename Perm , typename = std::enable_if_t<TiledArray::detail::is_permutation_v<Perm>>>
DecomposedTensor<T> mpqc::math::subt ( DecomposedTensor< T > const &  l,
DecomposedTensor< T > const &  r,
Perm const &  p 
)

◆ subt_to() [1/3]

template<typename T >
DecomposedTensor<T>& mpqc::math::subt_to ( DecomposedTensor< T > &  l,
const T  factor 
)

◆ subt_to() [2/3]

template<typename T >
DecomposedTensor<T>& mpqc::math::subt_to ( DecomposedTensor< T > &  l,
DecomposedTensor< T > const &  r 
)

◆ subt_to() [3/3]

template<typename T , typename F >
DecomposedTensor<T>& mpqc::math::subt_to ( DecomposedTensor< T > &  l,
DecomposedTensor< T > const &  r,
const F  factor 
)

◆ sum()

template<typename T >
T mpqc::math::sum ( DecomposedTensor< T > const &  t)

Returns the sum of a decomposed tensor tile Currently not recommended due to implementation details.

◆ symmetric_min_max_evals()

template<typename Array >
std::pair<double, double> mpqc::math::symmetric_min_max_evals ( Array const &  S)

◆ third_order_update()

template<typename Array >
void mpqc::math::third_order_update ( Array const &  S,
Array Z 
)

◆ tile_clr_storage() [1/2]

template<typename T >
unsigned long mpqc::math::tile_clr_storage ( Tile< DecomposedTensor< T >> const &  tile)

◆ tile_clr_storage() [2/2]

template<typename T >
unsigned long mpqc::math::tile_clr_storage ( TiledArray::Tile< math::DecomposedTensor< T >> const &  tile)

◆ tile_to_eigen() [1/2]

template<typename T >
Matrix<T> mpqc::math::tile_to_eigen ( TA::Tensor< T > const &  t)

◆ tile_to_eigen() [2/2]

template<typename T >
Matrix<T> mpqc::math::tile_to_eigen ( TA::Tile< DecomposedTensor< T >> const &  t)

◆ to_string()

std::string mpqc::math::to_string ( PetiteList::Symmetry  symmetry)

converts PetiteList::Symmetry to std::string

Parameters
symmetrya PetiteList::Symmetry object
Returns
string representation of symmetry

◆ trace()

template<typename T >
T mpqc::math::trace ( DecomposedTensor< T > const &  t)

Returns the trace of a decompose tensor Currently not recommended due to implementation details.

◆ two_way_decomposition()

DecomposedTensor< double > mpqc::math::two_way_decomposition ( DecomposedTensor< double > const &  t)

Returns an empty DecomposedTensor if the compression rank was to large.

◆ unary()

template<typename T , typename Op >
DecomposedTensor<T> mpqc::math::unary ( DecomposedTensor< T > const &  l,
Op &&  op 
)

◆ Z()

template<typename T1 , typename T2 >
constexpr std::complex<double> mpqc::math::Z ( T1  re,
T2  im 
)
inlineconstexpr