mpqc Namespace Reference

The top-level namespace for all Massively Parallel Quantum Chemistry package. More...

Namespaces

 cc
 
 constants
 
 cuda
 
 detail
 
 lcao
 
 libintx
 
 math
 
 meta
 
 molecule
 
 pbc
 
 pure
 
 python
 
 string
 
 TT
 
 util
 
 utility
 

Classes

class  AlgorithmException
 
class  AssertionFailed
 
class  Atom
 A class which holds the basic information for an atom. More...
 
class  AtomBasedCluster
 is the unit that holds a collection of atom based clusterables that go together. More...
 
class  AtomBasedClusterable
 The AtomBasedClusterable is a class that holds any clusterable type that is built up from atoms. More...
 
class  AtomBasedClusterConcept
 
class  AtomBasedClusterModel
 
class  AtomicData
 This singleton provides versioned access to the atomic data. More...
 
class  CartMolecularCoordinates
 
class  Cluster
 is the unit that holds a collection of clusterables that go together. More...
 
class  Clusterable
 The Clusterable is a class that holds any clusterable type. More...
 
class  ClusterConcept
 
class  ClusterModel
 
struct  CreIndicesUnoccupied
 checks if created indices are unoccupied in a string (or a pair of strings) More...
 
class  Debugger
 
class  DescribedClass
 
class  DirectArray
 
class  DirectTile
 A direct tile for integral construction. More...
 
class  enable_shared_from_this
 just like std::enable_shared_from_this but can appear multiple times in the subhierarchy of a class More...
 
class  Energy
 Taylor expansion of the molecular energy computed by a Wavefunction. More...
 
class  Exception
 
class  ExcitationEnergy
 
class  ExEnv
 Describes the execution environment of the program. More...
 
class  FeatureDisabled
 ‍** This is thrown when a verbose assertion fails. It is used by the More...
 
class  FeatureNotImplemented
 
class  FermionOccupationBlockString
 
class  FermionOccupationDBitString
 
class  FermionOccupationNBitString
 
class  FermionStringCompleteSet
 complete set of strings representable in a fixed number of states More...
 
class  FermionStringDenseSet
 contiguous container of strings More...
 
class  FermionStringSparseSet
 hashmap-based container of strings More...
 
class  FileOperationFailed
 
class  ForceLink
 
class  FormIO
 This utility class is used to print only on node 0 and to provide attractive indentation of output. More...
 
class  Formula
 Formula parses a string representation of quantum mechanical matrix elements and related expressions. More...
 
class  FormulaRegistry
 map Formula to Value object More...
 
class  FormulaRegistryBase
 
struct  FullCreAnnGenerator
 
struct  FullCreAnnPairGenerator
 
struct  FullDeterminantReplacementGenerator
 
class  FullFermionStringSetBuild
 
struct  FullStringReplacementGenerator
 
struct  FundamentalConstants
 
class  GetLongOpt
 Parse command line options. More...
 
class  GFRealPole
 
struct  HBCompositeGenerator
 Combines a generator producing tuples with HBTupleGenerator. More...
 
struct  HBDeterminantReplacementGenerator
 
class  HBSparseTensor
 Heat-bath sparse representation of a tensor of rank RI+RO More...
 
struct  HBStringReplacementGenerator
 
class  HBTupleGenerator
 
class  InputError
 
struct  InRangeFullDeterminantReplacementGenerator
 
class  KeyVal
 KeyVal specifies C++ primitive data (booleans, integers, reals, string) and user-defined objects obtained from JSON/XML/INFO input or by programmatic construction. More...
 
class  LimitExceeded
 
struct  Log
 
class  MaxIterExceeded
 
class  MemAllocFailed
 
class  MolecularCoordinates
 
class  MolecularFiniteDifferenceDerivative
 Evaluates FiniteDifferenceDerivative with respect to molecular coordinates. More...
 
class  MolecularFormula
 
class  Molecule
 Molecule is a class which contains a vector of AtomBasedClusterables. More...
 
class  MPQCInit
 This helper singleton class simplifies initialization of MPQC. More...
 
class  mpqcprintf
 
class  MPQCTask
 An MPQC computation. More...
 
class  NIST_v41_AtomicData
 
class  Operator
 Class to represent "logical" operators used in Formula. More...
 
struct  OpersInStringRange
 
struct  OrbitalTupleGenerator
 
class  PartitionedSparseOrbitalRange
 
class  PrimitiveOperator
 
struct  PrimitiveOperatorGenerator
 
struct  PrimitiveOperatorPairGenerator
 
class  ProgrammingError
 
class  Property
 
class  Provides
 Base for classes that provide Properties . More...
 
class  Provides< Property >
 
class  Registry
 a syntactically-sweet wrapper of std::map, not thread_safe for writing More...
 
struct  SlaterDeterminant
 
struct  SlaterDeterminant< String, std::enable_if_t< meta::is_string_v< String > > >
 
struct  SlaterDeterminant< StringRange, std::enable_if_t<!meta::is_string_v< StringRange > > >
 
class  SlaterDeterminantSet
 
struct  SlaterDeterminantSparseMatrix
 
class  SparseOrbitalRange
 a sparse set of orbitals More...
 
struct  StringReplacements
 Single replacements as a sparse matrix. More...
 
class  StringToSortedVector
 A (sparse) set of Slater determinants represented as a sparse matrix. More...
 
class  Supercell
 3-d box of unit cells with an odd number of cells in each dimension More...
 
class  SyscallFailed
 
class  SystemException
 
class  Task
 
class  ToleranceExceeded
 
struct  TupleOfOrbitalTuplesGenerator
 
class  UCSRDiagonalMatrixPreconditioner
 
class  UCSRMatrix
 
class  Uncomputable
 
class  Unit
 The Unit class is used to perform unit conversions. More...
 
class  UnitCell
 
class  UnitFactory
 
class  Wavefunction
 
class  WavefunctionProperty
 WavefunctionProperty computes a Taylor expansion of a molecular property using a visiting Wavefunction . More...
 
class  WavefunctionWorld
 It provides an execution context (madness::World) and a set of atoms. More...
 

Typedefs

template<typename Tile >
using tile_builder_fnref_t = std::variant< TiledArray::function_ref< Tile(const TiledArray::Range::index &, const TiledArray::Range &)>, TiledArray::function_ref< madness::Future< Tile >(const TiledArray::Range::index &, const TiledArray::Range &)> >
 
template<typename Policy = MPQC_DEFAULT_TA_POLICY_CLASS>
using DirectArrayD = DirectArray< MPQC_DEFAULT_REAL_TA_TILE_CLASS, Policy >
 
template<typename T >
using RowMatrix = ::Eigen::Matrix< T, ::Eigen::Dynamic, ::Eigen::Dynamic, ::Eigen::RowMajor >
 
using RowMatrixXd = RowMatrix< double >
 
template<typename T >
using EigenVector = ::Eigen::Matrix< T, ::Eigen::Dynamic, 1 >
 
template<typename T >
using RowEigenVector = ::Eigen::Matrix< T, 1, ::Eigen::Dynamic >
 
using MatrixZ = RowMatrix< std::complex< double > >
 
using RowMatrixXz = MatrixZ
 
using VectorZ = EigenVector< std::complex< double > >
 
using VectorD = EigenVector< double >
 
using Vector3d = Eigen::Vector3d
 
using Vector3i = Eigen::Vector3i
 
template<typename Tile = MPQC_DEFAULT_REAL_TA_TILE_CLASS, typename Policy = MPQC_DEFAULT_TA_POLICY_CLASS>
using DistArray = TiledArray::DistArray< Tile, Policy >
 
template<typename Policy = MPQC_DEFAULT_TA_POLICY_CLASS>
using DistArrayD = mpqc::DistArray< MPQC_DEFAULT_REAL_TA_TILE_CLASS, Policy >
 
template<typename Policy = MPQC_DEFAULT_TA_POLICY_CLASS>
using DistArrayZ = mpqc::DistArray< MPQC_DEFAULT_COMPLEX_TA_TILE_CLASS, Policy >
 
template<typename Tile = MPQC_DEFAULT_REAL_TA_TILE_CLASS, typename Policy = MPQC_DEFAULT_TA_POLICY_CLASS>
using DistArrayVector = TiledArray::DistArrayVector< Tile, Policy >
 
using DArray = TiledArray::DistArray< MPQC_DEFAULT_REAL_TA_TILE_CLASS, MPQC_DEFAULT_TA_POLICY_CLASS >
 
using ZArray = TiledArray::DistArray< MPQC_DEFAULT_REAL_TA_TILE_CLASS, MPQC_DEFAULT_TA_POLICY_CLASS >
 
using Array = DArray
 
using orbital_index_type = short
 
using PopulatedSparseOrbitalRange = PartitionedSparseOrbitalRange< 2 >
 
template<typename T >
using Describable = std::is_base_of< DescribedClass, std::decay_t< T > >
 
using time_point = std::chrono::high_resolution_clock::time_point
 
using TAPolicy = TiledArray::SparsePolicy
 

Enumerations

enum  NSpinCases { NSpinCases1 = 2, NSpinCases2 = 3 }
 
enum  NPureSpinCases { NPureSpinCases2 = 2 }
 
enum  SpinCase1 { SpinCase1::Any = -1, SpinCase1::Alpha = 0, SpinCase1::Beta = 1, SpinCase1::Invalid = 2 }
 
enum  { SpinCase1Any = static_cast<int>(SpinCase1::Any), SpinCase1Alpha = static_cast<int>(SpinCase1::Alpha), SpinCase1Beta = static_cast<int>(SpinCase1::Beta), SpinCase1Invalid = static_cast<int>(SpinCase1::Invalid) }
 
enum  SpinCase2 {
  SpinCase2::Any = -1, SpinCase2::AlphaBeta = 0, SpinCase2::AlphaAlpha = 1, SpinCase2::BetaBeta = 2,
  SpinCase2::Invalid = 3
}
 
enum  PureSpinCase2 { PureSpinCase2::Any = -1, PureSpinCase2::Singlet = 0, PureSpinCase2::Triplet = 1, PureSpinCase2::Invalid = 2 }
 

Functions

void announce (madness::World &world)
 
void initialize (int &argc, char **argv, madness::World &top_world, std::shared_ptr< GetLongOpt > opt=std::shared_ptr< GetLongOpt >())
 Static MPQC initializer. Must be called from the main thread of every process in the default MADWorld World before doing any MPQC-specific computation (e.g. before creating mpqc::MPQC objects). More...
 
void initialize (int &argc, char **argv, std::shared_ptr< GetLongOpt > opt=std::shared_ptr< GetLongOpt >())
 The version of initialize that initializes MADNESS and TiledArray and uses MPI_COMM_WORLD as the context. More...
 
void finalize ()
 Finalize MPQC. More...
 
void initialize_fpe ()
 Initializes the floating point exceptions. More...
 
std::tuple< std::shared_ptr< mpqc::KeyVal >, KeyVal::InputFormatmake_keyval (madness::World &world, const std::string &filename)
 Constructs a KeyVal object on every rank of world by reading file filename on rank 0. More...
 
std::shared_ptr< GetLongOptmake_options ()
 Creates a default options parser object for an MPQC executable. More...
 
std::tuple< std::string, std::optional< std::string > > process_options (madness::World &world, const std::shared_ptr< GetLongOpt > &options)
 Processes command-line options parsed by the options parser. More...
 
std::string to_string (MPQCInit::InputFormat f)
 
std::ostream & operator<< (std::ostream &os, Atom const &a)
 
Vector3d const & center (Atom const &a)
 Returns the center of the atom. More...
 
double mass (Atom const &a)
 Returns the mass of the atom. More...
 
int64_t total_atomic_number (Atom const &a)
 
size_t natoms (Atom const &a)
 
Vector3d const & center_of_mass (Atom const &a)
 Returns the center of the atom. More...
 
void update (AtomBasedCluster &abc, const std::vector< Atom > &atoms, size_t &pos)
 
double mass (AtomBasedClusterable const &ac)
 
int64_t total_atomic_number (AtomBasedClusterable const &ac)
 
size_t natoms (AtomBasedClusterable const &ac)
 
Vector3d const & center (AtomBasedClusterable const &ac)
 
Vector3d const & center_of_mass (AtomBasedClusterable const &ac)
 
std::vector< Atomcollapse_to_atoms (AtomBasedClusterable const &ac)
 
void update (AtomBasedClusterable &ac, const std::vector< Atom > &atoms, size_t &pos)
 
void apply (AtomBasedClusterable &ac, const std::function< void(Atom &)> &op)
 
std::vector< Atomcollapse_to_atoms (const Atom &a)
 collapse_to_atoms ends the recursive loop of the templated version of the function More...
 
void update (Atom &a, const std::vector< Atom > &atoms, size_t &pos)
 
void apply (Atom &a, const std::function< void(Atom &)> &op)
 
template<typename T >
std::vector< Atomcollapse_to_atoms (T const &t)
 collapse to atoms takes an iteratable list of types which provide an atoms function. More...
 
Molecule attach_hydrogens_and_kmeans (std::vector< AtomBasedClusterable > const &clusterables, size_t nclusters)
 
Molecule kmeans (std::vector< AtomBasedClusterable > const &clusterables, size_t nclusters)
 
std::vector< libint2::Atom > to_libint_atom (std::vector< Atom > const &atoms)
 
template<std::size_t N>
void increment (MolecularCoordinates *coords, std::array< size_t, N > coord_idxs, std::array< double, N > step)
 
template<std::size_t N>
void decrement (MolecularCoordinates *coords, std::array< size_t, N > coord_idxs, std::array< double, N > step)
 
std::ostream & operator<< (std::ostream &os, const MolecularCoordinates &coord)
 
bool zero_intersect (const std::pair< int64_t, int64_t > A, const std::pair< int64_t, int64_t > &B)
 return true if intervals A and B do not intersect More...
 
bool zero_intersect (const Supercell &A, const Supercell &B)
 
std::pair< int64_t, int64_t > intersect (const std::pair< int64_t, int64_t > A, const std::pair< int64_t, int64_t > &B)
 return the intersect of intervals A and B More...
 
Supercell intersect (const Supercell &A, const Supercell &B)
 return the intersect of Supercells A and B More...
 
std::ostream & operator<< (std::ostream &, UnitCell const &)
 print UnitCell to a std::ostream More...
 
std::shared_ptr< UnitCellto_reciprocal (const std::shared_ptr< const UnitCell > &unitcell)
 computes the unit cell of the reciprocal lattice More...
 
std::shared_ptr< UnitCellto_monkhorstpack_fbz_mesh (const std::shared_ptr< const UnitCell > &rcell, const Vector3i &mesh_size)
 computes the unit cell of the Monkhorst-Pack mesh for the first Brllouin zone of the reciprocal lattice More...
 
std::ostream & operator<< (std::ostream &, Molecule const &)
 Make Molecules printable. More...
 
std::vector< Formulapermutations (const Formula &formula)
 
template<typename Tile , typename Policy >
void d_ab_inplace (TA::DistArray< Tile, Policy > &abij, const Eigen::VectorXd &faa, const double e_i, const double e_j)
 
template<typename Tile , typename Policy >
double norm_vec_tensor (const std::vector< TA::DistArray< Tile, Policy >> &vec_array)
 
 MPQC_LCAO_WAVEFUNCTION_WORLD_PROPERTY (madness::World &world)
 
 MPQC_LCAO_WAVEFUNCTION_WORLD_PROPERTY (const std::shared_ptr< const Molecule > &atoms)
 
 MPQC_LCAO_WAVEFUNCTION_WORLD_PROPERTY (const std::shared_ptr< lcao::gaussian::OrbitalBasisRegistry > &basis_registry)
 
std::string to_string (const Unit &unit)
 
template<typename Tile , typename Policy >
std::vector< typename Tile::numeric_type > array_max_n (const TA::DistArray< Tile, Policy > &a, std::size_t n)
 
template<typename Tile , typename Policy >
std::vector< typename Tile::numeric_type > array_abs_max_n (const TA::DistArray< Tile, Policy > &a, std::size_t n)
 
template<typename Tile , typename Policy >
std::vector< std::pair< typename Tile::numeric_type, TA::Range::index > > array_max_n_index (const TA::DistArray< Tile, Policy > &a, std::size_t n)
 
template<typename Tile , typename Policy >
std::vector< std::pair< typename Tile::numeric_type, TA::Range::index > > array_abs_max_n_index (const TA::DistArray< Tile, Policy > &a, std::size_t n)
 
template<typename Tile , typename Policy >
TA::DistArray< Tile, Policy > fuse_arrays (const std::vector< TA::DistArray< Tile, Policy >> &arrays)
 
template<typename Tile , typename Policy >
TA::DistArray< Tile, Policy > split_fused_array (const TA::DistArray< Tile, Policy > &fused_array, std::size_t i, const TA::TiledRange &split_trange)
 
constexpr SpinCase1 flip (SpinCase1 sc)
 
template<bool UnitCoefficient = false, bool Lower = false, typename StringSet , typename ResultIterator , typename Coefficient , typename Tensor1 , typename Tensor2 , typename AReplacementFactory = decltype(make_str_repl<1, 1, typename StringSet::value_type>), typename BReplacementFactory = decltype(make_str_repl<1, 1, typename StringSet::value_type>), typename AAReplacementFactory = decltype(make_str_repl<2, 2, typename StringSet::value_type>), typename BBReplacementFactory = decltype(make_str_repl<2, 2, typename StringSet::value_type>), typename ABReplacementFactory = decltype(make_sd_repl<1, 1, 1, 1, StringSet>)>
auto compute_Hc_I (const SlaterDeterminantSet< StringSet > &dets, const typename SlaterDeterminantSet< StringSet >::const_iterator sditer, ResultIterator &sigma, const Coefficient c_I, const Tensor1 &h1, const Tensor2 &h2_1122, bool debug=false, const AReplacementFactory &make_A_replacement=make_str_repl< 1, 1, typename StringSet::value_type >, const BReplacementFactory &make_B_replacement=make_str_repl< 1, 1, typename StringSet::value_type >, const AAReplacementFactory &make_AA_replacement=make_str_repl< 2, 2, typename StringSet::value_type >, const BBReplacementFactory &make_BB_replacement=make_str_repl< 2, 2, typename StringSet::value_type >, const ABReplacementFactory &make_AB_replacement=make_sd_repl< 1, 1, 1, 1, StringSet >)
 
template<typename StringSet , typename ResultIterator , typename CoefficientIterator , typename Tensor1 , typename Tensor2 >
void compute_Hc (const SlaterDeterminantSet< StringSet > &dets, ResultIterator &sigma, const CoefficientIterator &C, const Tensor1 &h1, const Tensor2 &h2_1122, bool debug=false)
 
template<typename StringSet , typename Tensor1 , typename Tensor2 >
auto compute_H (const SlaterDeterminantSet< StringSet > &dets, const Tensor1 &h1, const Tensor2 &h2_1122, bool debug=false)
 
template<typename StringSet , typename Tensor1 , typename Tensor2 , typename AReplacementFactory = decltype(make_str_repl<1, 1, typename StringSet::value_type>), typename BReplacementFactory = decltype(make_str_repl<1, 1, typename StringSet::value_type>), typename AAReplacementFactory = decltype(make_str_repl<2, 2, typename StringSet::value_type>), typename BBReplacementFactory = decltype(make_str_repl<2, 2, typename StringSet::value_type>), typename ABReplacementFactory = decltype(make_sd_repl<1, 1, 1, 1, StringSet>)>
UCSRMatrix< double > compute_sparse_H (const SlaterDeterminantSet< StringSet > &dets, const Tensor1 &h1, const Tensor2 &h2_1122, bool debug=false, const AReplacementFactory &make_A_replacement=make_str_repl< 1, 1, typename StringSet::value_type >, const BReplacementFactory &make_B_replacement=make_str_repl< 1, 1, typename StringSet::value_type >, const AAReplacementFactory &make_AA_replacement=make_str_repl< 2, 2, typename StringSet::value_type >, const BBReplacementFactory &make_BB_replacement=make_str_repl< 2, 2, typename StringSet::value_type >, const ABReplacementFactory &make_AB_replacement=make_sd_repl< 1, 1, 1, 1, StringSet >)
 
template<bool LowerH = true, typename String , typename Tensor1 , typename Tensor2 >
void compute_sparse_H_v2 (UCSRMatrix< double > &H, const SlaterDeterminantSparseMatrix< String, SpinCase1::Alpha > &dets, const Tensor1 &h1, const Tensor2 &h2_1122, bool debug=false)
 
template<class CharT , class Traits , std::size_t RR, std::size_t RC, bool Abs, typename T >
std::basic_ostream< CharT, Traits > & operator<< (std::basic_ostream< CharT, Traits > &os, const HBSparseTensor< RR, RC, Abs, T > &hbtensor)
 
template<std::size_t RC, std::size_t RA, typename FString , bool AbsHBData = true, typename T = double>
auto make_hb_str_repl (const FString &str, std::shared_ptr< HBSparseTensor< RC, RA, AbsHBData, T >> hb_data, typename HBSparseTensor< RC, RA, AbsHBData, T >::magnitude_type eps)
 makes an HB generator of string replacements More...
 
template<std::size_t RCa, std::size_t RCb, std::size_t RAa, std::size_t RAb, typename StringOrStringRange , bool AbsHBData = true, typename T = double>
auto make_hb_sd_repl (const SlaterDeterminant< StringOrStringRange > &sd, std::shared_ptr< HBSparseTensor< RCa+RCb, RAa+RAb, AbsHBData, T >> hb_data, typename HBSparseTensor< RCa+RCb, RAa+RAb, AbsHBData, T >::magnitude_type eps)
 makes an HB generator of determinant replacements More...
 
template<std::size_t RCa, std::size_t RCb, std::size_t RAa, std::size_t RAb, typename String , bool AbsHBData = true, typename T = double, typename ReplFilter = const decltype(detail::always_true)&>
void generate_hb_sd_repl (SlaterDeterminantSparseMatrix< String > &hb_sdset, std::shared_ptr< HBSparseTensor< RCa+RCb, RAa+RAb, AbsHBData, T >> hb_data, const SlaterDeterminant< String > &sd, typename HBSparseTensor< RCa+RCb, RAa+RAb, AbsHBData, T >::magnitude_type eps, bool debug=false, ReplFilter &&repl_filter=detail::always_true)
 HB-generates determinants and inserts into a SlaterDeterminantSet. More...
 
template<std::size_t RCa, std::size_t RCb, std::size_t RAa, std::size_t RAb, typename StringSet , bool AbsHBData = true, typename T = double>
void generate_hb_sd_repl (SlaterDeterminantSet< StringSet > &hb_sdset, std::shared_ptr< HBSparseTensor< RCa+RCb, RAa+RAb, AbsHBData, T >> hb_data, const SlaterDeterminant< StringSet > &sd, typename HBSparseTensor< RCa+RCb, RAa+RAb, AbsHBData, T >::magnitude_type eps, bool debug=false)
 HB-generates determinants and inserts into a SlaterDeterminantSet. More...
 
template<std::size_t RC, std::size_t RA = RC, typename FString >
auto make_str_repl (const FString &str)
 
template<typename FermionStringSet >
FermionStringSet make_stringset (size_t m, size_t n)
 
template<typename SDSet >
SDSet make_sdset (size_t m, size_t na, size_t nb)
 
template<std::size_t RCa, std::size_t RCb, std::size_t RAa = RCa, std::size_t RAb = RCb, typename StringRange >
auto make_sd_repl (const SlaterDeterminant< StringRange > &det)
 
template<std::size_t RCa, std::size_t RCb, std::size_t RAa = RCa, std::size_t RAb = RCb, typename StringRange >
auto make_det_inrepl (const SlaterDeterminant< StringRange > &det, const StringRange &astrange, const StringRange &bstrange)
 
template<std::size_t Nc, std::size_t Na = Nc, bool Lower = false, typename String >
auto make_str_repl_list (const FermionStringSparseSet< String > &sset)
 
template<std::size_t Nc, std::size_t Na = Nc, bool Lower = false, typename String >
auto make_str_repl_list (const FermionStringSparseSet< String > &sset_from, const FermionStringSparseSet< String > &sset_to)
 
template<std::size_t I, std::size_t NP>
auto get (const PartitionedSparseOrbitalRange< NP > &porbs)
 adapts PartitionedSparseOrbitalRange for structured bindings More...
 
auto occ_range (const PopulatedSparseOrbitalRange &poporbs)
 
auto nocc (const PopulatedSparseOrbitalRange &poporbs)
 
auto uocc_range (const PopulatedSparseOrbitalRange &poporbs)
 
auto nuocc (const PopulatedSparseOrbitalRange &poporbs)
 
template<typename String , typename = std::enable_if_t<meta::is_string_v<String>>>
const auto & value (const String &str)
 
template<typename T1 , typename T2 >
const auto & alpha (const std::pair< T1, T2 > &p)
 
template<typename T1 , typename T2 >
auto & alpha (std::pair< T1, T2 > &p)
 
template<typename T1 , typename T2 >
const auto & alpha (const std::tuple< T1, T2 > &p)
 
template<typename T1 , typename T2 >
auto & alpha (std::tuple< T1, T2 > &p)
 
template<typename T1 , typename T2 >
const auto & beta (const std::pair< T1, T2 > &p)
 
template<typename T1 , typename T2 >
auto & beta (std::pair< T1, T2 > &p)
 
template<typename T1 , typename T2 >
const auto & beta (const std::tuple< T1, T2 > &p)
 
template<typename T1 , typename T2 >
auto & beta (std::tuple< T1, T2 > &p)
 
template<typename StringRange , typename = std::enable_if_t<!meta::is_string_v<StringRange>>>
const auto & alpha (const SlaterDeterminant< StringRange > &det)
 
template<typename StringRange , typename = std::enable_if_t<!meta::is_string_v<StringRange>>>
auto & alpha (SlaterDeterminant< StringRange > &det)
 
template<typename StringRange , typename = std::enable_if_t<!meta::is_string_v<StringRange>>>
const auto & alphastr (const SlaterDeterminant< StringRange > &det)
 
template<typename String , typename = std::enable_if_t<meta::is_string_v<String>>>
auto alpha (const SlaterDeterminant< String > &det)
 
template<typename String , typename = std::enable_if_t<meta::is_string_v<String>>>
auto alpha (SlaterDeterminant< String > &det)
 
template<typename String >
std::enable_if_t< meta::is_string_v< String >, const String & > alphastr (const SlaterDeterminant< String > &det)
 
template<typename StringRange , typename = std::enable_if_t<!meta::is_string_v<StringRange>>>
decltype(auto) beta (const SlaterDeterminant< StringRange > &det)
 
template<typename StringRange , typename = std::enable_if_t<!meta::is_string_v<StringRange>>>
auto & beta (SlaterDeterminant< StringRange > &det)
 
template<typename StringRange , typename = std::enable_if_t<!meta::is_string_v<StringRange>>>
const auto & betastr (const SlaterDeterminant< StringRange > &det)
 
template<typename String , typename = std::enable_if_t<meta::is_string_v<String>>>
auto beta (const SlaterDeterminant< String > &det)
 
template<typename String , typename = std::enable_if_t<meta::is_string_v<String>>>
auto beta (SlaterDeterminant< String > &det)
 
template<typename String >
std::enable_if_t< meta::is_string_v< String >, const String & > betastr (const SlaterDeterminant< String > &det)
 
template<typename RangeA , typename RangeB , typename = std::enable_if_t<TiledArray::detail::is_range_v<RangeA> && TiledArray::detail::is_range_v<RangeB>>>
decltype(auto) begin (const std::pair< RangeA, RangeB > &rng)
 
template<typename RangeA , typename RangeB , typename = std::enable_if_t<TiledArray::detail::is_range_v<RangeA> && TiledArray::detail::is_range_v<RangeB>>>
decltype(auto) end (const std::pair< RangeA, RangeB > &rng)
 
auto occ_range (const std::pair< PopulatedSparseOrbitalRange, PopulatedSparseOrbitalRange > &poporbs)
 
auto uocc_range (const std::pair< PopulatedSparseOrbitalRange, PopulatedSparseOrbitalRange > &poporbs)
 
template<class CharT , class Traits , typename StringRange >
std::basic_ostream< CharT, Traits > & operator<< (std::basic_ostream< CharT, Traits > &os, const SlaterDeterminant< StringRange > &sd)
 
template<SpinCase1 Spin = SpinCase1::Alpha, typename StringRange , template< typename... > class HashMap>
auto make_sdmat (const SlaterDeterminantSet< StringRange, HashMap > &sdset, bool recompute_ordinals=false)
 
constexpr unsigned int nspincases1 (bool spin_polarized)
 Returns the number of unique spin cases (1 or 2) More...
 
constexpr unsigned int nspincases2 (bool spin_polarized)
 Returns the number of unique combinations of 2 spin cases (1 or 3) More...
 
constexpr unsigned int npurespincases2 ()
 Returns the number of pure 2 spin cases. More...
 
constexpr SpinCase1 case1 (SpinCase2 S)
 returns the first spin case of the 2-spin S More...
 
constexpr SpinCase1 case2 (SpinCase2 S)
 returns the second spin case of the 2-spin S More...
 
constexpr SpinCase2 case12 (SpinCase1 S1, SpinCase1 S2)
 combines 2 spins to give 1 2-spin More...
 
constexpr SpinCase1 other (SpinCase1 S)
 given 1-spin return the other 1-spin More...
 
constexpr const char * to_string (SpinCase1 S)
 
constexpr const char * to_string (SpinCase2 S)
 
constexpr const char * to_string (PureSpinCase2 S)
 
SpinCase1 to_spincase1 (std::string key)
 
PureSpinCase2 to_purespincase2 (std::string key)
 
SpinCase2 to_spincase2 (std::string key)
 
std::string prepend_spincase (SpinCase2 S, const std::string &R, bool lowercase=false)
 Prepend string representation of S to R and return. More...
 
std::string prepend_spincase (PureSpinCase2 S, const std::string &R, bool lowercase=false)
 Prepend string representation of S to R and return. More...
 
std::string prepend_spincase (SpinCase1 S, const std::string &R, bool lowercase=false)
 Prepend string representation of S to R and return. More...
 
template<class CharT , class Traits , size_t Ns>
std::basic_ostream< CharT, Traits > & operator<< (std::basic_ostream< CharT, Traits > &os, const FermionOccupationNBitString< Ns > &x)
 
FermionOccupationDBitString operator+ (const FermionOccupationDBitString &s1, const FermionOccupationDBitString &s2)
 
template<class CharT , class Traits >
std::basic_ostream< CharT, Traits > & operator<< (std::basic_ostream< CharT, Traits > &os, const FermionOccupationDBitString &x)
 
template<class CharT , class Traits >
std::basic_ostream< CharT, Traits > & operator<< (std::basic_ostream< CharT, Traits > &os, const FermionOccupationBlockString &x)
 
FermionOccupationBlockString operator+ (const FermionOccupationBlockString &s1, const FermionOccupationBlockString &s2)
 
template<typename String >
auto make_string (orbital_index_type size, std::initializer_list< orbital_index_type > block_sizes)
 string factory More...
 
template<typename String >
auto make_string (orbital_index_type size, std::initializer_list< orbital_index_type > specified_states, bool default_occupancy)
 string factory More...
 
template<typename StringSet >
StringSet make_product (const StringSet &strset1, const StringSet &strset2, std::function< void(const typename StringSet::const_iterator &, const typename StringSet::const_iterator &, const typename StringSet::const_iterator &)> op={})
 constructs tensor product of two StringSet objects More...
 
template<typename Iterator , typename = std::enable_if_t<TiledArray::detail::is_pair_v< std::decay_t<typename std::iterator_traits< std::remove_reference_t<Iterator>>::value_type>>>>
const auto & value (Iterator &&it)
 
template<typename Iterator , typename = std::enable_if_t<TiledArray::detail::is_pair_v< std::decay_t<typename std::iterator_traits< std::remove_reference_t<Iterator>>::value_type>>>>
auto ordinal (Iterator &&it)
 
template<size_t NC1, size_t NA1, size_t NC2, size_t NA2>
bool operator== (const PrimitiveOperator< NC1, NA1 > &op1, const PrimitiveOperator< NC2, NA2 > &op2)
 
template<class CharT , class Traits , size_t Nc, size_t Na>
std::basic_ostream< CharT, Traits > & operator<< (std::basic_ostream< CharT, Traits > &os, const PrimitiveOperator< Nc, Na > &o)
 
template<typename CreIndices , typename AnnIndices >
auto make_oper_from_tuple (const std::tuple< const CreIndices &, const AnnIndices & > &creidx_annidx)
 
template<typename CreIndicesTuple , typename AnnIndicesTuple >
auto make_opers_from_tuples (const std::tuple< const CreIndicesTuple &, const AnnIndicesTuple & > creidxt_annidxt)
 
template<std::size_t Nc, std::size_t Na, typename FString >
PrimitiveOperator< Nc, Na > make_op (const FString &to, const FString &from)
 computes operator converting string from to to More...
 
template<typename T >
std::size_t hash_value (const T &obj)
 
template<typename T >
auto chop (const UCSRMatrix< T > &m, T eps=1000 *std::numeric_limits< T >::epsilon())
 
template<typename Ch , typename Tr , typename T >
std::basic_ostream< Ch, Tr > & operator<< (std::basic_ostream< Ch, Tr > &os, const UCSRMatrix< T > &mat)
 
template<bool SelfAdjoint = false, typename T >
auto to_eigen_dense (const UCSRMatrix< T > &spmat)
 
template<typename Tile , typename Policy >
void minimize_storage (TA::DistArray< Tile, Policy > &A)
 
template<typename Tile , typename Policy >
void minimize_storage (TA::DistArrayVector< Tile, Policy > &A)
 
template<typename Tile >
void minimize_storage (TA::DistArray< Tile, TA::SparsePolicy > &A, double truncate_threshold)
 
template<typename Tile >
void minimize_storage (TA::DistArrayVector< Tile, TA::SparsePolicy > &A, double truncate_threshold)
 
void minimize_storage (TA::DistArray< TA::Tile< math::DecomposedTensor< double >>, TA::SparsePolicy > &A, double clr_threshold)
 
template std::ios & indent< char > (std::ios &)
 
template std::ios & decindent< char > (std::ios &)
 
template std::ios & incindent< char > (std::ios &)
 
template std::ios & skipnextindent< char > (std::ios &)
 
template std::wios & indent< wchar_t > (std::wios &)
 
template std::wios & decindent< wchar_t > (std::wios &)
 
template std::wios & incindent< wchar_t > (std::wios &)
 
template std::wios & skipnextindent< wchar_t > (std::wios &)
 
std::basic_ostream< char > & operator<< (std::basic_ostream< char > &o, const mpqcprintf< char > &s)
 
std::basic_ostream< wchar_t > & operator<< (std::basic_ostream< wchar_t > &o, const mpqcprintf< wchar_t > &s)
 
template<typename Char >
std::basic_ios< Char > & indent (std::basic_ios< Char > &o)
 
template<typename Char >
std::basic_ios< Char > & decindent (std::basic_ios< Char > &o)
 
template<typename Char >
std::basic_ios< Char > & incindent (std::basic_ios< Char > &o)
 
template<typename Char >
std::basic_ios< Char > & skipnextindent (std::basic_ios< Char > &o)
 
template<typename Char >
std::basic_ostream< Char > & operator<< (std::basic_ostream< Char > &, const mpqcprintf< Char > &)
 
template<typename Char , typename... Args>
mpqcprintf< Char > printf (const Char *fmt, Args &&... args)
 
template<typename T >
std::shared_ptr< T > make_shared_from_list (std::initializer_list< typename detail::init_list_ctor< T >::type > list)
 wrapper for std::shared_ptr that can accept a std::initializer_list More...
 
madness::World * get_comm_world ()
 
void set_default_world (madness::World *world)
 
void set_local_world (madness::World *world)
 
madness::World * get_default_world ()
 
madness::World * get_local_world ()
 
KeyVal operator+ (const KeyVal &first, const KeyVal &second)
 
KeyVal make_kv ()
 make an empty KeyVal More...
 
template<typename... Args>
std::enable_if_t< sizeof...(Args) % 2==1, KeyValmake_kv (Args &&... args)
 report an error if make_kv receives an odd number of arguments More...
 
template<typename Key , typename Value , typename... RestOfArgs>
KeyVal make_kv (Key &&key, Value &&value, RestOfArgs &&... rest_of_args)
 
template<typename... Args>
auto make_list (Args &&... args)
 
void launch_gdb_xterm ()
 Use this to launch GNU debugger in xterm. More...
 
void launch_lldb_xterm ()
 Use this to launch LLVM debugger in xterm. More...
 
madness::World * get_default_world (const KeyVal &kv)
 accesses the default madness::World object More...
 
template<typename String >
void print_meminfo (String &&tag, const madness::World &world= *get_default_world(), const std::string filename_prefix=FormIO::default_basename()+std::string(".meminfo"))
 print aggregate memory stats to file <basename>.meminfo.<comm_world_rank>, tagged with \п tag More...
 
template<typename Property , typename Provider >
bool provides (const std::shared_ptr< Provider > &provider)
 
template<typename Property , typename Object >
std::enable_if_t< detail::has_provider< Property >::value, typename Property::Provider * > as_provider_of (Object *ptr)
 
template<typename Property , typename Object >
std::enable_if_t< detail::has_provider< Property >::value, typename Property::Provider * > as_provider_of (const std::shared_ptr< Object > &ptr)
 
template<typename Property , typename Provider , typename... EvaluateArgs>
std::enable_if_t<(detail::has_provider< Property >::value &&!std::is_const< Property >::value), void > evaluate (Property &property, Provider &provider, EvaluateArgs... eval_args)
 Evaluates property using provider. More...
 
template<typename Property , typename Provider >
std::enable_if_t< detail::has_provider< Property >::value, Property & > operator<< (Property &property, Provider &provider)
 Evaluates property using provider. More...
 
time_point now ()
 
std::chrono::system_clock::time_point system_now ()
 
double duration_in_s (time_point const &t0, time_point const &t1)
 
int64_t duration_in_ns (time_point const &t0, time_point const &t1)
 
time_point fenced_now (madness::World &world)
 
time_point now (madness::World &world, bool fence)
 
std::vector< Atomcollapse_to_atoms (AtomBasedCluster const &abc)
 
void apply (AtomBasedCluster &abc, const std::function< void(Atom &)> &op)
 
std::ostream & operator<< (std::ostream &, AtomBasedCluster const &)
 print the cluster by printing each of its elements More...
 
Vector3d const & center (AtomBasedCluster const &c)
 Returns the Center of mass of the cluster. More...
 
double mass (AtomBasedCluster const &c)
 print the cluster by printing each of its elements More...
 
int64_t total_atomic_number (AtomBasedCluster const &c)
 print the cluster by printing each of its elements More...
 
Vector3d const & center_of_mass (AtomBasedCluster const &c)
 print the cluster by printing each of its elements More...
 
size_t natoms (AtomBasedCluster const &c)
 print the cluster by printing each of its elements More...
 
void set_center (AtomBasedCluster &c, Vector3d const &point)
 print the cluster by printing each of its elements More...
 
void remove_clusterables (AtomBasedCluster &c)
 print the cluster by printing each of its elements More...
 
template<typename T >
void attach_clusterable (AtomBasedCluster &c, T t)
 print the cluster by printing each of its elements More...
 
void update_center (AtomBasedCluster &c)
 print the cluster by printing each of its elements More...
 
std::ostream & operator<< (std::ostream &, Cluster const &)
 print the cluster by printing each of its elements More...
 
Vector3d const & center (Cluster const &c)
 print the cluster by printing each of its elements More...
 

Variables

bool mpqc_is_linked = true
 
ForceLink< Molecule, UnitCellchemistry_molecule_forcelink
 
const std::unordered_map< wchar_t, std::string > greek_to_english_name
 
ForceLink< Energy, ExcitationEnergy, GFRealPole, MolecularFiniteDifferenceDerivative< 1, double > > chemistry_qc_properties_forcelink
 
ForceLink< WavefunctionWorld, Wavefunctionchemistry_qc_wfn_forcelink
 
const auto log = Log{}
 

Typedef Documentation

◆ Array

using mpqc::Array = typedef DArray

◆ DArray

using mpqc::DArray = typedef TiledArray::DistArray<MPQC_DEFAULT_REAL_TA_TILE_CLASS, MPQC_DEFAULT_TA_POLICY_CLASS>

◆ Describable

template<typename T >
using mpqc::Describable = typedef std::is_base_of<DescribedClass, std::decay_t<T> >

◆ DirectArrayD

template<typename Policy = MPQC_DEFAULT_TA_POLICY_CLASS>
using mpqc::DirectArrayD = typedef DirectArray<MPQC_DEFAULT_REAL_TA_TILE_CLASS, Policy>

◆ DistArray

template<typename Tile = MPQC_DEFAULT_REAL_TA_TILE_CLASS, typename Policy = MPQC_DEFAULT_TA_POLICY_CLASS>
using mpqc::DistArray = typedef TiledArray::DistArray<Tile, Policy>

◆ DistArrayD

template<typename Policy = MPQC_DEFAULT_TA_POLICY_CLASS>
using mpqc::DistArrayD = typedef mpqc::DistArray<MPQC_DEFAULT_REAL_TA_TILE_CLASS, Policy>

◆ DistArrayVector

template<typename Tile = MPQC_DEFAULT_REAL_TA_TILE_CLASS, typename Policy = MPQC_DEFAULT_TA_POLICY_CLASS>
using mpqc::DistArrayVector = typedef TiledArray::DistArrayVector<Tile, Policy>

◆ DistArrayZ

template<typename Policy = MPQC_DEFAULT_TA_POLICY_CLASS>
using mpqc::DistArrayZ = typedef mpqc::DistArray<MPQC_DEFAULT_COMPLEX_TA_TILE_CLASS, Policy>

◆ EigenVector

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

◆ MatrixZ

using mpqc::MatrixZ = typedef RowMatrix<std::complex<double> >

◆ orbital_index_type

typedef short mpqc::orbital_index_type

this type indexes orbitals

Note
this is signed to be nullable

◆ PopulatedSparseOrbitalRange

◆ RowEigenVector

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

◆ RowMatrix

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

◆ RowMatrixXd

using mpqc::RowMatrixXd = typedef RowMatrix<double>

◆ RowMatrixXz

using mpqc::RowMatrixXz = typedef MatrixZ

◆ TAPolicy

using mpqc::TAPolicy = typedef TiledArray::SparsePolicy

◆ tile_builder_fnref_t

template<typename Tile >
using mpqc::tile_builder_fnref_t = typedef std::variant< TiledArray::function_ref<Tile(const TiledArray::Range::index &, const TiledArray::Range &)>, TiledArray::function_ref<madness::Future<Tile>( const TiledArray::Range::index &, const TiledArray::Range &)> >

a type-erasing reference to {synchronous,asynchronous} Tile builder

Template Parameters
Tilea TiledArray tile type

◆ time_point

using mpqc::time_point = typedef std::chrono::high_resolution_clock::time_point

◆ Vector3d

using mpqc::Vector3d = typedef Eigen::Vector3d

◆ Vector3i

using mpqc::Vector3i = typedef Eigen::Vector3i

◆ VectorD

using mpqc::VectorD = typedef EigenVector<double>

◆ VectorZ

using mpqc::VectorZ = typedef EigenVector<std::complex<double> >

◆ ZArray

using mpqc::ZArray = typedef TiledArray::DistArray<MPQC_DEFAULT_REAL_TA_TILE_CLASS, MPQC_DEFAULT_TA_POLICY_CLASS>

Function Documentation

◆ announce()

void mpqc::announce ( madness::World &  world)

◆ apply() [1/2]

void mpqc::apply ( Atom a,
const std::function< void(Atom &)> &  op 
)
inline

◆ apply() [2/2]

void mpqc::apply ( AtomBasedCluster abc,
const std::function< void(Atom &)> &  op 
)

applies a function to each atom in AtomBasedCluster

Parameters
abcthe AtomBasedCluster object
opa function taking non-const lvalue ref to Atom and returning void, stored in std::function

◆ array_abs_max_n()

template<typename Tile , typename Policy >
std::vector<typename Tile::numeric_type> mpqc::array_abs_max_n ( const TA::DistArray< Tile, Policy > &  a,
std::size_t  n 
)

◆ array_abs_max_n_index()

template<typename Tile , typename Policy >
std::vector<std::pair<typename Tile::numeric_type, TA::Range::index> > mpqc::array_abs_max_n_index ( const TA::DistArray< Tile, Policy > &  a,
std::size_t  n 
)

◆ array_max_n()

template<typename Tile , typename Policy >
std::vector<typename Tile::numeric_type> mpqc::array_max_n ( const TA::DistArray< Tile, Policy > &  a,
std::size_t  n 
)

◆ array_max_n_index()

template<typename Tile , typename Policy >
std::vector<std::pair<typename Tile::numeric_type, TA::Range::index> > mpqc::array_max_n_index ( const TA::DistArray< Tile, Policy > &  a,
std::size_t  n 
)

◆ as_provider_of() [1/2]

template<typename Property , typename Object >
std::enable_if_t<detail::has_provider<Property>::value, typename Property::Provider*> mpqc::as_provider_of ( const std::shared_ptr< Object > &  ptr)
Returns
object pointer cast to Property::Provider

◆ as_provider_of() [2/2]

template<typename Property , typename Object >
std::enable_if_t<detail::has_provider<Property>::value, typename Property::Provider*> mpqc::as_provider_of ( Object *  ptr)
Returns
object pointer cast to Property::Provider

◆ attach_clusterable()

template<typename T >
void mpqc::attach_clusterable ( AtomBasedCluster c,
t 
)
inline

print the cluster by printing each of its elements

◆ attach_hydrogens_and_kmeans()

Molecule mpqc::attach_hydrogens_and_kmeans ( std::vector< AtomBasedClusterable > const &  clusterables,
size_t  nclusters 
)

◆ center() [1/3]

Vector3d const& mpqc::center ( Atom const &  a)
inline

Returns the center of the atom.

◆ center() [2/3]

Vector3d const& mpqc::center ( AtomBasedCluster const &  c)
inline

Returns the Center of mass of the cluster.

This function exists to allow interfacing with generic clustering code. Overloading center to return the center of mass is in some sense specializing center for atoms.

◆ center() [3/3]

Vector3d const& mpqc::center ( Cluster const &  c)
inline

print the cluster by printing each of its elements

◆ center_of_mass() [1/2]

Vector3d const& mpqc::center_of_mass ( Atom const &  a)
inline

Returns the center of the atom.

Center of mass is part of the clustering interface.

◆ center_of_mass() [2/2]

Vector3d const& mpqc::center_of_mass ( AtomBasedCluster const &  c)
inline

print the cluster by printing each of its elements

◆ chop()

template<typename T >
auto mpqc::chop ( const UCSRMatrix< T > &  m,
eps = 1000 * std::numeric_limits<T>::epsilon() 
)

◆ collapse_to_atoms() [1/3]

std::vector< Atom > mpqc::collapse_to_atoms ( AtomBasedCluster const &  abc)

converts an AtomBasedCluster to a vector of Atom's

Parameters
abcthe AtomBasedCluster object
Returns
the atoms sequence

◆ collapse_to_atoms() [2/3]

std::vector<Atom> mpqc::collapse_to_atoms ( const Atom a)
inline

collapse_to_atoms ends the recursive loop of the templated version of the function

The return type is a vector because that simplifies the interface of the overload which does captures all types that are not Atoms.

◆ collapse_to_atoms() [3/3]

template<typename T >
std::vector<Atom> mpqc::collapse_to_atoms ( T const &  t)

collapse to atoms takes an iteratable list of types which provide an atoms function.

The idea is that something like a vector of clusters of clusters will call call atoms on each cluster of clusters, which will internally call collapse_to_atoms on it's type recursing until the type being passed to collapse_to_atoms is just an atom.

◆ compute_H()

template<typename StringSet , typename Tensor1 , typename Tensor2 >
auto mpqc::compute_H ( const SlaterDeterminantSet< StringSet > &  dets,
const Tensor1 &  h1,
const Tensor2 &  h2_1122,
bool  debug = false 
)

◆ compute_Hc()

template<typename StringSet , typename ResultIterator , typename CoefficientIterator , typename Tensor1 , typename Tensor2 >
void mpqc::compute_Hc ( const SlaterDeterminantSet< StringSet > &  dets,
ResultIterator &  sigma,
const CoefficientIterator &  C,
const Tensor1 &  h1,
const Tensor2 &  h2_1122,
bool  debug = false 
)

◆ compute_Hc_I()

template<bool UnitCoefficient = false, bool Lower = false, typename StringSet , typename ResultIterator , typename Coefficient , typename Tensor1 , typename Tensor2 , typename AReplacementFactory = decltype(make_str_repl<1, 1, typename StringSet::value_type>), typename BReplacementFactory = decltype(make_str_repl<1, 1, typename StringSet::value_type>), typename AAReplacementFactory = decltype(make_str_repl<2, 2, typename StringSet::value_type>), typename BBReplacementFactory = decltype(make_str_repl<2, 2, typename StringSet::value_type>), typename ABReplacementFactory = decltype(make_sd_repl<1, 1, 1, 1, StringSet>)>
auto mpqc::compute_Hc_I ( const SlaterDeterminantSet< StringSet > &  dets,
const typename SlaterDeterminantSet< StringSet >::const_iterator  sditer,
ResultIterator &  sigma,
const Coefficient  c_I,
const Tensor1 &  h1,
const Tensor2 &  h2_1122,
bool  debug = false,
const AReplacementFactory &  make_A_replacement = make_str_repl<1, 1, typename StringSet::value_type>,
const BReplacementFactory &  make_B_replacement = make_str_repl<1, 1, typename StringSet::value_type>,
const AAReplacementFactory &  make_AA_replacement = make_str_repl<2, 2, typename StringSet::value_type>,
const BBReplacementFactory &  make_BB_replacement = make_str_repl<2, 2, typename StringSet::value_type>,
const ABReplacementFactory &  make_AB_replacement = make_sd_repl<1, 1, 1, 1, StringSet> 
)
Template Parameters
UnitCoefficientif true, will ignore the value of c_I
Lowerif true, will only compute contributions from the lower traingle of the Hamiltonian
StringSet
ResultIterator
Coefficient
Tensor1
Tensor2
AReplacementFactory
BReplacementFactory
AAReplacementFactory
BBReplacementFactory
ABReplacementFactory
Parameters
dets
sditer
sigma
c_I
h1
h2_1122
debug
make_A_replacement
make_B_replacement
make_AA_replacement
make_BB_replacement
make_AB_replacement
Returns

◆ compute_sparse_H()

template<typename StringSet , typename Tensor1 , typename Tensor2 , typename AReplacementFactory = decltype(make_str_repl<1, 1, typename StringSet::value_type>), typename BReplacementFactory = decltype(make_str_repl<1, 1, typename StringSet::value_type>), typename AAReplacementFactory = decltype(make_str_repl<2, 2, typename StringSet::value_type>), typename BBReplacementFactory = decltype(make_str_repl<2, 2, typename StringSet::value_type>), typename ABReplacementFactory = decltype(make_sd_repl<1, 1, 1, 1, StringSet>)>
UCSRMatrix<double> mpqc::compute_sparse_H ( const SlaterDeterminantSet< StringSet > &  dets,
const Tensor1 &  h1,
const Tensor2 &  h2_1122,
bool  debug = false,
const AReplacementFactory &  make_A_replacement = make_str_repl<1, 1, typename StringSet::value_type>,
const BReplacementFactory &  make_B_replacement = make_str_repl<1, 1, typename StringSet::value_type>,
const AAReplacementFactory &  make_AA_replacement = make_str_repl<2, 2, typename StringSet::value_type>,
const BBReplacementFactory &  make_BB_replacement = make_str_repl<2, 2, typename StringSet::value_type>,
const ABReplacementFactory &  make_AB_replacement = make_sd_repl<1, 1, 1, 1, StringSet> 
)

◆ compute_sparse_H_v2()

template<bool LowerH = true, typename String , typename Tensor1 , typename Tensor2 >
void mpqc::compute_sparse_H_v2 ( UCSRMatrix< double > &  H,
const SlaterDeterminantSparseMatrix< String, SpinCase1::Alpha > &  dets,
const Tensor1 &  h1,
const Tensor2 &  h2_1122,
bool  debug = false 
)
Warning
this only computes alpha-beta replacements so far!

◆ d_ab_inplace()

template<typename Tile , typename Policy >
void mpqc::d_ab_inplace ( TA::DistArray< Tile, Policy > &  abij,
const Eigen::VectorXd &  faa,
const double  e_i,
const double  e_j 
)

◆ decindent()

template<typename Char >
std::basic_ios< Char > & mpqc::decindent ( std::basic_ios< Char > &  o)

◆ decindent< char >()

template std::ios& mpqc::decindent< char > ( std::ios &  )

◆ decindent< wchar_t >()

template std::wios& mpqc::decindent< wchar_t > ( std::wios &  )

◆ decrement()

template<std::size_t N>
void mpqc::decrement ( MolecularCoordinates coords,
std::array< size_t, N >  coord_idxs,
std::array< double, N >  step 
)

◆ duration_in_ns()

int64_t mpqc::duration_in_ns ( time_point const &  t0,
time_point const &  t1 
)
inline

◆ duration_in_s()

double mpqc::duration_in_s ( time_point const &  t0,
time_point const &  t1 
)
inline

◆ evaluate()

template<typename Property , typename Provider , typename... EvaluateArgs>
std::enable_if_t<(detail::has_provider<Property>::value && !std::is_const<Property>::value), void> mpqc::evaluate ( Property property,
Provider &  provider,
EvaluateArgs...  eval_args 
)

Evaluates property using provider.

◆ fenced_now()

time_point mpqc::fenced_now ( madness::World &  world)
inline

◆ fuse_arrays()

template<typename Tile , typename Policy >
TA::DistArray<Tile, Policy> mpqc::fuse_arrays ( const std::vector< TA::DistArray< Tile, Policy >> &  arrays)

Fuse a vector of N dimension Array of the same TiledRange into a N+1 dimension Array. The new dimension will be the leading dimension, and will be blocked by 1. All the arrays must have the same TiledRange object.

Parameters
arraysa vector of DistArray object with the same TiledRanges

copy the data from tile

write to blocks of fused_array

◆ generate_hb_sd_repl() [1/2]

template<std::size_t RCa, std::size_t RCb, std::size_t RAa, std::size_t RAb, typename StringSet , bool AbsHBData = true, typename T = double>
void mpqc::generate_hb_sd_repl ( SlaterDeterminantSet< StringSet > &  hb_sdset,
std::shared_ptr< HBSparseTensor< RCa+RCb, RAa+RAb, AbsHBData, T >>  hb_data,
const SlaterDeterminant< StringSet > &  sd,
typename HBSparseTensor< RCa+RCb, RAa+RAb, AbsHBData, T >::magnitude_type  eps,
bool  debug = false 
)

HB-generates determinants and inserts into a SlaterDeterminantSet.

Template Parameters
RCacreation rank for first case
RAaannihilation rank first case
RCbcreation rank for second case
RAbannihilation rank for second case
StringSeta string set type
AbsHBDatawhether uses a HBData object containing absolute values (magnitudes) or not
Tthe values contained by the HBData object
Parameters
[in,out]hb_sdseta set of Slater determinants into which to insert generated determinants
[in]hb_dataHB data tensor
[in]sdthe Slater determinant from which to generate HB replacements; does not have to be present in hb_sdset
[in]epsthe HB screening threshold; replacement i is generated if hb_data[i] >= eps
debugif true, will print extra debugging output to mpqc::ExEnv::out0()

◆ generate_hb_sd_repl() [2/2]

template<std::size_t RCa, std::size_t RCb, std::size_t RAa, std::size_t RAb, typename String , bool AbsHBData = true, typename T = double, typename ReplFilter = const decltype(detail::always_true)&>
void mpqc::generate_hb_sd_repl ( SlaterDeterminantSparseMatrix< String > &  hb_sdset,
std::shared_ptr< HBSparseTensor< RCa+RCb, RAa+RAb, AbsHBData, T >>  hb_data,
const SlaterDeterminant< String > &  sd,
typename HBSparseTensor< RCa+RCb, RAa+RAb, AbsHBData, T >::magnitude_type  eps,
bool  debug = false,
ReplFilter &&  repl_filter = detail::always_true 
)

HB-generates determinants and inserts into a SlaterDeterminantSet.

Template Parameters
RCacreation rank for first case
RAaannihilation rank first case
RCbcreation rank for second case
RAbannihilation rank for second case
Stringstring type
ReplFilterreplacement filter, a callable with signature bool(const std::tuple<OperA, OperB>&)
AbsHBDatawhether uses a HBData object containing absolute values (magnitudes) or not
Tthe values contained by the HBData object
Parameters
[in,out]hb_sdseta set of Slater determinants into which to insert generated determinants
[in]hb_dataHB data tensor
[in]sdthe Slater determinant from which to generate HB replacements; does not have to be present in hb_sdset
[in]epsthe HB screening threshold; replacement i is generated if hb_data[i] >= eps
[in]debugif true, will print extra debugging output to mpqc::ExEnv::out0()
[in]repl_filteran object that filters out produced replacements; the default is to not filter out any replacements

◆ get()

template<std::size_t I, std::size_t NP>
auto mpqc::get ( const PartitionedSparseOrbitalRange< NP > &  porbs)

adapts PartitionedSparseOrbitalRange for structured bindings

Template Parameters
Iindex of the element to access
NPnumber of partitions in the PartitionedSparseOrbitalRange
Precondition
I<=NP get<0> returns the full range, get<1> returns the size of first partition.
Warning
The size of the last partition is not accessible via this function.

◆ get_comm_world()

madness::World * mpqc::get_comm_world ( )
Returns
the pointer to the World object associated with the MPI COMM_WORLD communicator

◆ get_default_world() [1/2]

madness::World * mpqc::get_default_world ( )
Returns
the pointer to the default madness::World object, either set by the most recent call to set_default_world(), or corresponding to the MPI COMM_WORLD communicator

◆ get_default_world() [2/2]

madness::World * mpqc::get_default_world ( const KeyVal kv)

accesses the default madness::World object

First reads the "$:world" from kv , if not given will use get_default_world()

Parameters
[in]KeyVal
Returns
reference to the default madness::World object, either set by the most recent call to set_default_world(), or corresponding to the MPI COMM_WORLD communicator

◆ get_local_world()

madness::World * mpqc::get_local_world ( )
Returns
the pointer to the local madness::World object, either set by the most recent call to set_local_world(), or nullptr

◆ hash_value()

template<typename T >
std::size_t mpqc::hash_value ( const T &  obj)

adapts class T such that T::hash_value() returns std::size_t for the use with boost.hash

◆ incindent()

template<typename Char >
std::basic_ios< Char > & mpqc::incindent ( std::basic_ios< Char > &  o)

◆ incindent< char >()

template std::ios& mpqc::incindent< char > ( std::ios &  )

◆ incindent< wchar_t >()

template std::wios& mpqc::incindent< wchar_t > ( std::wios &  )

◆ increment()

template<std::size_t N>
void mpqc::increment ( MolecularCoordinates coords,
std::array< size_t, N >  coord_idxs,
std::array< double, N >  step 
)

◆ indent()

template<typename Char >
std::basic_ios< Char > & mpqc::indent ( std::basic_ios< Char > &  o)

◆ indent< char >()

template std::ios& mpqc::indent< char > ( std::ios &  )

◆ indent< wchar_t >()

template std::wios& mpqc::indent< wchar_t > ( std::wios &  )

◆ intersect() [1/2]

std::pair<int64_t, int64_t> mpqc::intersect ( const std::pair< int64_t, int64_t >  A,
const std::pair< int64_t, int64_t > &  B 
)
inline

return the intersect of intervals A and B

◆ intersect() [2/2]

Supercell mpqc::intersect ( const Supercell A,
const Supercell B 
)
inline

return the intersect of Supercells A and B

◆ kmeans()

Molecule mpqc::kmeans ( std::vector< AtomBasedClusterable > const &  clusterables,
size_t  nclusters 
)

◆ launch_gdb_xterm()

void mpqc::launch_gdb_xterm ( )

Use this to launch GNU debugger in xterm.

◆ launch_lldb_xterm()

void mpqc::launch_lldb_xterm ( )

Use this to launch LLVM debugger in xterm.

◆ make_det_inrepl()

template<std::size_t RCa, std::size_t RCb, std::size_t RAa = RCa, std::size_t RAb = RCb, typename StringRange >
auto mpqc::make_det_inrepl ( const SlaterDeterminant< StringRange > &  det,
const StringRange &  astrange,
const StringRange &  bstrange 
)

◆ make_hb_sd_repl()

template<std::size_t RCa, std::size_t RCb, std::size_t RAa, std::size_t RAb, typename StringOrStringRange , bool AbsHBData = true, typename T = double>
auto mpqc::make_hb_sd_repl ( const SlaterDeterminant< StringOrStringRange > &  sd,
std::shared_ptr< HBSparseTensor< RCa+RCb, RAa+RAb, AbsHBData, T >>  hb_data,
typename HBSparseTensor< RCa+RCb, RAa+RAb, AbsHBData, T >::magnitude_type  eps 
)

makes an HB generator of determinant replacements

Template Parameters
RCacreation rank for alpha spin
RCbcreation rank for beta spin
RAaannihilation rank for alpha spin
RAbannihilation rank for beta spin
StringOrStringRangestring or string range type
AbsHBDatawhether uses a HBData object containing absolute values (magnitudes) or not
Tthe values contained by the HBData object
Parameters
[in]sdthe Slater determinant from which to generate HB replacements
[in]hb_dataHB data tensor
[in]epsthe HB screening threshold; replacement i is generated if hb_data[i] >= eps

◆ make_hb_str_repl()

template<std::size_t RC, std::size_t RA, typename FString , bool AbsHBData = true, typename T = double>
auto mpqc::make_hb_str_repl ( const FString &  str,
std::shared_ptr< HBSparseTensor< RC, RA, AbsHBData, T >>  hb_data,
typename HBSparseTensor< RC, RA, AbsHBData, T >::magnitude_type  eps 
)

makes an HB generator of string replacements

Template Parameters
RCcreation rank
RAannihilation rank
FStringstring type
AbsHBDatawhether uses a HBData object containing absolute values (magnitudes) or not
Tthe values contained by the HBData object
Parameters
[in]strthe string from which to generate HB replacements
[in]hb_dataHB data tensor
[in]epsthe HB screening threshold; replacement i is generated if hb_data[i] >= eps

◆ make_kv() [1/3]

KeyVal mpqc::make_kv ( )
inline

make an empty KeyVal

◆ make_kv() [2/3]

template<typename... Args>
std::enable_if_t<sizeof...(Args) % 2 == 1, KeyVal> mpqc::make_kv ( Args &&...  args)
inline

report an error if make_kv receives an odd number of arguments

◆ make_kv() [3/3]

template<typename Key , typename Value , typename... RestOfArgs>
KeyVal mpqc::make_kv ( Key &&  key,
Value &&  value,
RestOfArgs &&...  rest_of_args 
)
inline

makes an empty KeyVal, then assigns {key,value} pair and the rest of kv pairs in rest_of_args

◆ make_list()

template<typename... Args>
auto mpqc::make_list ( Args &&...  args)

make a std::array for working around std::initializer_list's non-deducibility in many contexts, e.g. in make_kv

◆ make_op()

template<std::size_t Nc, std::size_t Na, typename FString >
PrimitiveOperator<Nc, Na> mpqc::make_op ( const FString &  to,
const FString &  from 
)

computes operator converting string from to to

◆ make_oper_from_tuple()

template<typename CreIndices , typename AnnIndices >
auto mpqc::make_oper_from_tuple ( const std::tuple< const CreIndices &, const AnnIndices & > &  creidx_annidx)

creates a PrimitiveOperator from a tuple of refs to std::array of created and annihilated indices

Template Parameters
CreIndicesstd::array of integers
AnnIndicesstd::array of integers
Parameters
[in]creidx_annidxstd::tuple of refs to arrays of created and annihilated indices

◆ make_opers_from_tuples()

template<typename CreIndicesTuple , typename AnnIndicesTuple >
auto mpqc::make_opers_from_tuples ( const std::tuple< const CreIndicesTuple &, const AnnIndicesTuple & >  creidxt_annidxt)

creates a 2-tuple of PrimitiveOperator objects from a tuple of refs to 2-tuple of std::array of created and annihilated indices

Template Parameters
CreIndicesTuplestd::tuple<const std::array&, const std::array&>
AnnIndicesTuplestd::tuple<const std::array&, const std::array&>
Parameters
[in]creidxt_annidxtcreated and annihilated indices

◆ make_product()

template<typename StringSet >
StringSet mpqc::make_product ( const StringSet &  strset1,
const StringSet &  strset2,
std::function< void(const typename StringSet::const_iterator &, const typename StringSet::const_iterator &, const typename StringSet::const_iterator &)>  op = {} 
)

constructs tensor product of two StringSet objects

Template Parameters
StringSeta string set type
Parameters
[in]strset1first StringSet object
[in]strset1second StringSet object
opthis callable will be called, if nonnull, for each product string with the iterators pointing to the string from strset1 and strset2 sets used to construct the product, as well as the iterator to the product string in the result
Returns
the tensor product of this with other
Postcondition
result.size()==first.size()*second.size()

◆ make_sd_repl()

template<std::size_t RCa, std::size_t RCb, std::size_t RAa = RCa, std::size_t RAb = RCb, typename StringRange >
auto mpqc::make_sd_repl ( const SlaterDeterminant< StringRange > &  det)

◆ make_sdmat()

template<SpinCase1 Spin = SpinCase1::Alpha, typename StringRange , template< typename... > class HashMap>
auto mpqc::make_sdmat ( const SlaterDeterminantSet< StringRange, HashMap > &  sdset,
bool  recompute_ordinals = false 
)

◆ make_sdset()

template<typename SDSet >
SDSet mpqc::make_sdset ( size_t  m,
size_t  na,
size_t  nb 
)

creates a set of determinants by distributing na particles of spin-up and nb particles of spin-down

Parameters
mnumber of states
nanumber of spin-up particles
nbnumber of spin-down particles

◆ make_shared_from_list()

template<typename T >
std::shared_ptr<T> mpqc::make_shared_from_list ( std::initializer_list< typename detail::init_list_ctor< T >::type >  list)

wrapper for std::shared_ptr that can accept a std::initializer_list

Template Parameters
Ta type
Parameters
[in]lista std::initializer_list from which to construct an object of type T
Returns
std::shard_ptr<T>
See also
detail::init_list_ctor

◆ make_str_repl()

template<std::size_t RC, std::size_t RA = RC, typename FString >
auto mpqc::make_str_repl ( const FString &  str)

◆ make_str_repl_list() [1/2]

template<std::size_t Nc, std::size_t Na = Nc, bool Lower = false, typename String >
auto mpqc::make_str_repl_list ( const FermionStringSparseSet< String > &  sset)

generates a set of string replacements between strings within the given (sparse) set of strings

◆ make_str_repl_list() [2/2]

template<std::size_t Nc, std::size_t Na = Nc, bool Lower = false, typename String >
auto mpqc::make_str_repl_list ( const FermionStringSparseSet< String > &  sset_from,
const FermionStringSparseSet< String > &  sset_to 
)

generates a set of string replacements from strings within one given (sparse) set of strings to another

◆ make_string() [1/2]

template<typename String >
auto mpqc::make_string ( orbital_index_type  size,
std::initializer_list< orbital_index_type block_sizes 
)

string factory

Template Parameters
Stringa string type
Parameters
sizestring size
block_sizessizes of blocks with alternating occupancies
Synopsis
make_string<String>(20, {5, 4, 3}) creates a string with first 5 states occupied, next 4 unoccupied, next 3 occupied, and the rest unoccupied.

◆ make_string() [2/2]

template<typename String >
auto mpqc::make_string ( orbital_index_type  size,
std::initializer_list< orbital_index_type specified_states,
bool  default_occupancy 
)

string factory

Template Parameters
Stringa string type
Parameters
sizestring size
specified_statesspecifies the states whose occupancy is opposite of default_occupancy
default_occupancyspecifies the default occupancies of all states not mentioned in specified_states
Synopsis
make_string<String>(20, {0, 1, 5}, false) creates a string with states 0, 1, and 5 occupied, the rest unoccupied.

◆ make_stringset()

template<typename FermionStringSet >
FermionStringSet mpqc::make_stringset ( size_t  m,
size_t  n 
)

Makes a set of strings by distributing n particles in m states

Parameters
mnumber of states
nnumber of particles

◆ mass() [1/2]

double mpqc::mass ( Atom const &  a)
inline

Returns the mass of the atom.

◆ mass() [2/2]

double mpqc::mass ( AtomBasedCluster const &  c)
inline

print the cluster by printing each of its elements

◆ minimize_storage() [1/5]

void mpqc::minimize_storage ( TA::DistArray< TA::Tile< math::DecomposedTensor< double >>, TA::SparsePolicy > &  A,
double  clr_threshold 
)
inline

◆ minimize_storage() [2/5]

template<typename Tile , typename Policy >
void mpqc::minimize_storage ( TA::DistArray< Tile, Policy > &  A)
inline

◆ minimize_storage() [3/5]

template<typename Tile >
void mpqc::minimize_storage ( TA::DistArray< Tile, TA::SparsePolicy > &  A,
double  truncate_threshold 
)
inline

◆ minimize_storage() [4/5]

template<typename Tile , typename Policy >
void mpqc::minimize_storage ( TA::DistArrayVector< Tile, Policy > &  A)
inline

◆ minimize_storage() [5/5]

template<typename Tile >
void mpqc::minimize_storage ( TA::DistArrayVector< Tile, TA::SparsePolicy > &  A,
double  truncate_threshold 
)
inline

◆ MPQC_LCAO_WAVEFUNCTION_WORLD_PROPERTY() [1/3]

mpqc::MPQC_LCAO_WAVEFUNCTION_WORLD_PROPERTY ( const std::shared_ptr< const Molecule > &  atoms)

◆ MPQC_LCAO_WAVEFUNCTION_WORLD_PROPERTY() [2/3]

mpqc::MPQC_LCAO_WAVEFUNCTION_WORLD_PROPERTY ( const std::shared_ptr< lcao::gaussian::OrbitalBasisRegistry > &  basis_registry)

◆ MPQC_LCAO_WAVEFUNCTION_WORLD_PROPERTY() [3/3]

mpqc::MPQC_LCAO_WAVEFUNCTION_WORLD_PROPERTY ( madness::World &  world)

◆ natoms() [1/2]

size_t mpqc::natoms ( Atom const &  a)
inline

of atoms in an atom, always returns 1

◆ natoms() [2/2]

size_t mpqc::natoms ( AtomBasedCluster const &  c)
inline

print the cluster by printing each of its elements

◆ nocc()

auto mpqc::nocc ( const PopulatedSparseOrbitalRange poporbs)
inline

◆ norm_vec_tensor()

template<typename Tile , typename Policy >
double mpqc::norm_vec_tensor ( const std::vector< TA::DistArray< Tile, Policy >> &  vec_array)

◆ now() [1/2]

time_point mpqc::now ( )
inline

◆ now() [2/2]

time_point mpqc::now ( madness::World &  world,
bool  fence 
)
inline

return fence() time based on bool fence

◆ nuocc()

auto mpqc::nuocc ( const PopulatedSparseOrbitalRange poporbs)
inline

◆ occ_range()

auto mpqc::occ_range ( const PopulatedSparseOrbitalRange poporbs)
inline

◆ operator+() [1/3]

FermionOccupationBlockString mpqc::operator+ ( const FermionOccupationBlockString s1,
const FermionOccupationBlockString s2 
)
inline

concatenates two strings

Note
asserts that the current states of s1 and s2 are nonnull
Postcondition
result.size()==s1.size()+s2.size()

◆ operator+() [2/3]

FermionOccupationDBitString mpqc::operator+ ( const FermionOccupationDBitString s1,
const FermionOccupationDBitString s2 
)
inline

concatenates this string and another string

Parameters
s1a string
s2another string
Returns
the result of concatenating 2 strings
Note
asserts that the current states of s1 and s2 are nonnull

◆ operator+() [3/3]

KeyVal mpqc::operator+ ( const KeyVal first,
const KeyVal second 
)

union of two KeyVal objects

Note
obsoletes sc::AggregateKeyVal from MPQC3

◆ operator<<() [1/16]

template<typename Property , typename Provider >
std::enable_if_t<detail::has_provider<Property>::value, Property&> mpqc::operator<< ( Property property,
Provider &  provider 
)

Evaluates property using provider.

◆ operator<<() [2/16]

template<typename Ch , typename Tr , typename T >
std::basic_ostream<Ch, Tr>& mpqc::operator<< ( std::basic_ostream< Ch, Tr > &  os,
const UCSRMatrix< T > &  mat 
)

◆ operator<<() [3/16]

template<typename Char >
std::basic_ostream<Char>& mpqc::operator<< ( std::basic_ostream< Char > &  ,
const mpqcprintf< Char > &   
)

◆ operator<<() [4/16]

std::basic_ostream< char > & mpqc::operator<< ( std::basic_ostream< char > &  o,
const mpqcprintf< char > &  s 
)

◆ operator<<() [5/16]

template<class CharT , class Traits >
std::basic_ostream<CharT, Traits>& mpqc::operator<< ( std::basic_ostream< CharT, Traits > &  os,
const FermionOccupationBlockString x 
)

◆ operator<<() [6/16]

template<class CharT , class Traits >
std::basic_ostream<CharT, Traits>& mpqc::operator<< ( std::basic_ostream< CharT, Traits > &  os,
const FermionOccupationDBitString x 
)

◆ operator<<() [7/16]

template<class CharT , class Traits , size_t Ns>
std::basic_ostream<CharT, Traits>& mpqc::operator<< ( std::basic_ostream< CharT, Traits > &  os,
const FermionOccupationNBitString< Ns > &  x 
)

concatenates two strings

Note
asserts that the current states of s1 and s2 are nonnull
Postcondition
result.size()==s1.size()+s2.size()

◆ operator<<() [8/16]

template<class CharT , class Traits , std::size_t RR, std::size_t RC, bool Abs, typename T >
std::basic_ostream<CharT, Traits>& mpqc::operator<< ( std::basic_ostream< CharT, Traits > &  os,
const HBSparseTensor< RR, RC, Abs, T > &  hbtensor 
)

◆ operator<<() [9/16]

template<class CharT , class Traits , size_t Nc, size_t Na>
std::basic_ostream<CharT, Traits>& mpqc::operator<< ( std::basic_ostream< CharT, Traits > &  os,
const PrimitiveOperator< Nc, Na > &  o 
)

◆ operator<<() [10/16]

template<class CharT , class Traits , typename StringRange >
std::basic_ostream<CharT, Traits>& mpqc::operator<< ( std::basic_ostream< CharT, Traits > &  os,
const SlaterDeterminant< StringRange > &  sd 
)

◆ operator<<() [11/16]

std::basic_ostream< wchar_t > & mpqc::operator<< ( std::basic_ostream< wchar_t > &  o,
const mpqcprintf< wchar_t > &  s 
)

◆ operator<<() [12/16]

std::ostream & mpqc::operator<< ( std::ostream &  os,
AtomBasedCluster const &  c 
)

print the cluster by printing each of its elements

◆ operator<<() [13/16]

std::ostream & mpqc::operator<< ( std::ostream &  os,
Cluster const &  c 
)

print the cluster by printing each of its elements

◆ operator<<() [14/16]

std::ostream & mpqc::operator<< ( std::ostream &  os,
UnitCell const &  unitcell 
)

print UnitCell to a std::ostream

◆ operator<<() [15/16]

std::ostream & mpqc::operator<< ( std::ostream &  os,
Atom const &  a 
)

◆ operator<<() [16/16]

std::ostream & mpqc::operator<< ( std::ostream &  os,
const MolecularCoordinates coord 
)

◆ operator==()

template<size_t NC1, size_t NA1, size_t NC2, size_t NA2>
bool mpqc::operator== ( const PrimitiveOperator< NC1, NA1 > &  op1,
const PrimitiveOperator< NC2, NA2 > &  op2 
)

◆ ordinal()

template<typename Iterator , typename = std::enable_if_t<TiledArray::detail::is_pair_v< std::decay_t<typename std::iterator_traits< std::remove_reference_t<Iterator>>::value_type>>>>
auto mpqc::ordinal ( Iterator &&  it)

syntactic sugar to extract ordinal from string set iterator {value,ordinal}

Template Parameters
Iteratora string set iterator; dereferences to a {value,ordinal} pair
Parameters
[in]ita string set iterator

◆ permutations()

std::vector< Formula > mpqc::permutations ( const Formula formula)

◆ print_meminfo()

template<typename String >
void mpqc::print_meminfo ( String &&  tag,
const madness::World &  world = *get_default_world(),
const std::string  filename_prefix = FormIO::default_basename() + std::string(".meminfo") 
)

print aggregate memory stats to file <basename>.meminfo.<comm_world_rank>, tagged with \п tag

Parameters
[in]tagrecord tag as any string type, e.g. const char[] , std::string , or std::wstring
[in]worldthe World object to use for this collective call
[in]filename_prefixfile name prefix; the default value is <basename>.meminfo
Note
To make print_meminfo() usable execute ::madness::print_meminfo_enable() prior to calling this.
must set global locale properly with std::locale::global() if tag has nontrivial encoding
Warning
this DOES fence using world

◆ printf()

template<typename Char , typename... Args>
mpqcprintf<Char> mpqc::printf ( const Char *  fmt,
Args &&...  args 
)

◆ provides()

template<typename Property , typename Provider >
bool mpqc::provides ( const std::shared_ptr< Provider > &  provider)
Returns
true if Provider can provide Property

◆ remove_clusterables()

void mpqc::remove_clusterables ( AtomBasedCluster c)
inline

print the cluster by printing each of its elements

◆ set_center()

void mpqc::set_center ( AtomBasedCluster c,
Vector3d const &  point 
)
inline

print the cluster by printing each of its elements

◆ set_default_world()

void mpqc::set_default_world ( madness::World *  world)
Parameters
[in]worldthe new default madness::World object

◆ set_local_world()

void mpqc::set_local_world ( madness::World *  world)
Parameters
[in]worldthe new local madness::World object

◆ skipnextindent()

template<typename Char >
std::basic_ios< Char > & mpqc::skipnextindent ( std::basic_ios< Char > &  o)

◆ skipnextindent< char >()

template std::ios& mpqc::skipnextindent< char > ( std::ios &  )

◆ skipnextindent< wchar_t >()

template std::wios& mpqc::skipnextindent< wchar_t > ( std::wios &  )

◆ split_fused_array()

template<typename Tile , typename Policy >
TA::DistArray<Tile, Policy> mpqc::split_fused_array ( const TA::DistArray< Tile, Policy > &  fused_array,
std::size_t  i,
const TA::TiledRange &  split_trange 
)

split the ith mode of a fused Array into a Array object

Parameters
fused_arraya DistArray object with leading dimension blocked by 1
ithe number of mode to split
trangeTiledRange of the split Array object

copy the data from tile

write to blocks of fused_array

◆ system_now()

std::chrono::system_clock::time_point mpqc::system_now ( )
inline

◆ to_eigen_dense()

template<bool SelfAdjoint = false, typename T >
auto mpqc::to_eigen_dense ( const UCSRMatrix< T > &  spmat)

◆ to_libint_atom()

std::vector<libint2::Atom> mpqc::to_libint_atom ( std::vector< Atom > const &  )

Function takes mpqc::molecule::Atom vector and converts it to a vector of libint atoms.

◆ to_monkhorstpack_fbz_mesh()

std::shared_ptr< UnitCell > mpqc::to_monkhorstpack_fbz_mesh ( const std::shared_ptr< const UnitCell > &  rcell,
const Vector3i mesh_size 
)

computes the unit cell of the Monkhorst-Pack mesh for the first Brllouin zone of the reciprocal lattice

Parameters
rcellthe unit cell of the reciprocal lattice
mesh_sizethe size of the mesh in each dimension (each element must be a positive odd number)
Returns
shared_ptr to the unit cell of the MP mesh for FBZ
Note
for nonperiodic dimensions the primitive vectors are zero

◆ to_reciprocal()

std::shared_ptr< UnitCell > mpqc::to_reciprocal ( const std::shared_ptr< const UnitCell > &  unitcell)

computes the unit cell of the reciprocal lattice

Parameters
unitcellthe unit cell of the direct lattice
Returns
shared_ptr to the UnitCell of the reciprocal lattice
Note
for nonperiodic dimensions the primitive vectors are zero

◆ to_string()

std::string mpqc::to_string ( const Unit unit)
inline

◆ total_atomic_number() [1/2]

int64_t mpqc::total_atomic_number ( Atom const &  a)
inline
Returns
the atomic number of the Atom object

◆ total_atomic_number() [2/2]

int64_t mpqc::total_atomic_number ( AtomBasedCluster const &  c)
inline

print the cluster by printing each of its elements

◆ uocc_range()

auto mpqc::uocc_range ( const PopulatedSparseOrbitalRange poporbs)
inline

◆ update() [1/2]

void mpqc::update ( Atom a,
const std::vector< Atom > &  atoms,
size_t &  pos 
)
inline

◆ update() [2/2]

void mpqc::update ( AtomBasedCluster abc,
const std::vector< Atom > &  atoms,
size_t &  pos 
)

updates an AtomBasedCluster using a vector of Atom's

Parameters
abcthe AtomBasedCluster object
atomsthe vector of Atom objects
posthe index of the next Atom object in Atoms to be used

◆ update_center()

void mpqc::update_center ( AtomBasedCluster c)
inline

print the cluster by printing each of its elements

◆ value() [1/2]

template<typename String , typename = std::enable_if_t<meta::is_string_v<String>>>
const auto& mpqc::value ( const String &  str)

◆ value() [2/2]

template<typename Iterator , typename = std::enable_if_t<TiledArray::detail::is_pair_v< std::decay_t<typename std::iterator_traits< std::remove_reference_t<Iterator>>::value_type>>>>
const auto& mpqc::value ( Iterator &&  it)

syntactic sugar to extract value from string set iterator {value,ordinal}

Template Parameters
Iteratora string set iterator; dereferences to a {value,ordinal} pair
Parameters
[in]ita string set iterator

◆ zero_intersect() [1/2]

bool mpqc::zero_intersect ( const std::pair< int64_t, int64_t >  A,
const std::pair< int64_t, int64_t > &  B 
)
inline

return true if intervals A and B do not intersect

◆ zero_intersect() [2/2]

bool mpqc::zero_intersect ( const Supercell A,
const Supercell B 
)
inline

Variable Documentation

◆ chemistry_molecule_forcelink

ForceLink<Molecule, UnitCell> mpqc::chemistry_molecule_forcelink

◆ chemistry_qc_properties_forcelink

ForceLink<Energy, ExcitationEnergy, GFRealPole, MolecularFiniteDifferenceDerivative<1, double> > mpqc::chemistry_qc_properties_forcelink

◆ chemistry_qc_wfn_forcelink

ForceLink<WavefunctionWorld, Wavefunction> mpqc::chemistry_qc_wfn_forcelink

◆ greek_to_english_name

const std::unordered_map<wchar_t, std::string> mpqc::greek_to_english_name
Initial value:
= {
{L'Α', "ALPHA"}, {L'Β', "BETA"}, {L'Γ', "GAMMA"}, {L'Δ', "DELTA"},
{L'Ε', "EPSILON"}, {L'Ζ', "ZETA"}, {L'Η', "ETA"}, {L'Θ', "THETA"},
{L'Ι', "IOTA"}, {L'Κ', "KAPPA"}, {L'Λ', "LAMBDA"}, {L'Μ', "MU"},
{L'Ν', "NU"}, {L'Ξ', "XI"}, {L'Ο', "OMICRON"}, {L'Π', "PI"},
{L'Ρ', "RHO"}, {L'Σ', "SIGMA"}, {L'Τ', "TAU"}, {L'Υ', "UPSILON"},
{L'Φ', "PHI"}, {L'Χ', "CHI"}, {L'Ψ', "PSI"}, {L'Ω', "OMEGA"},
{L'α', "alpha"}, {L'β', "beta"}, {L'γ', "gamma"}, {L'δ', "delta"},
{L'ε', "epsilon"}, {L'ζ', "zeta"}, {L'η', "eta"}, {L'θ', "theta"},
{L'ι', "iota"}, {L'κ', "kappa"}, {L'λ', "lambda"}, {L'μ', "mu"},
{L'ν', "nu"}, {L'ξ', "xi"}, {L'ο', "omicron"}, {L'π', "pi"},
{L'ρ', "rho"}, {L'σ', "sigma"}, {L'τ', "tau"}, {L'υ', "upsilon"},
{L'φ', "phi"}, {L'χ', "chi"}, {L'ψ', "psi"}, {L'ω', "omega"}}

◆ log

const auto mpqc::log = Log{}

◆ mpqc_is_linked

bool mpqc::mpqc_is_linked = true
const Vector3d & L() const
Definition: group.cpp:74