mpqc::lcao::cc::CCSD< Tile, Policy > Class Template Reference
Collaboration diagram for mpqc::lcao::cc::CCSD< Tile, Policy >:

Documentation

template<typename Tile, typename Policy>
class mpqc::lcao::cc::CCSD< Tile, Policy >

CCSD class that computed CCSD energy

Public Types

using TArray = TA::DistArray< Tile, Policy >
 
using AOFactory = gaussian::AOFactory< Tile, Policy >
 
- Public Types inherited from mpqc::lcao::LCAOWavefunction< Tile, Policy >
using AOFactoryType = gaussian::AOFactory< Tile, Policy >
 
using ArrayType = typename AOFactoryType::TArray
 
using DirectTArray = typename AOFactoryType::DirectTArray
 
using LCAOFactoryType = LCAOFactory< Tile, Policy >
 

Public Member Functions

 CCSD ()=default
 
 CCSD (const KeyVal &kv)
 
virtual ~CCSD ()
 
void obsolete () override
 
TArray t1 () const
 
TArray t2 () const
 
bool is_df () const
 
bool is_cp_ccsdt () const
 
double cp_ccsdt_precision () const
 
double cp_ccsdt_rank () const
 
bool is_reduced_abcd_memory () const
 
void set_t1 (const TArray &t1)
 
void set_t2 (const TArray &t2)
 
bool verbose () const
 
void set_target_precision (double precision)
 
const AOFactory::DirectTArrayget_direct_ao_integral () const
 
- Public Member Functions inherited from mpqc::lcao::LCAOWavefunction< Tile, Policy >
 LCAOWavefunction (const KeyVal &kv)
 The KeyVal constructor. More...
 
virtual ~LCAOWavefunction ()=default
 
LCAOFactoryTypelcao_factory ()
 
const LCAOFactoryTypelcao_factory () const
 
AOFactoryTypeao_factory ()
 
const AOFactoryTypeao_factory () const
 
void obsolete () override
 
virtual void init_sdref (std::shared_ptr< lcao::Wavefunction > ref_wfn, double target_ref_precision)
 initializes the orbital spaces associated with a (closed-shell, for now) single-determinant reference state More...
 
const std::shared_ptr< const ::mpqc::utility::TRange1Engine > & trange1_engine () const
 
bool is_frozen_core () const
 
int charge () const
 
size_t occ_block () const
 
size_t unocc_block () const
 
auto localizer () const
 
auto cluster_occupied_orbitals () const
 

Protected Member Functions

bool can_evaluate (Energy *energy) override
 
void evaluate (Energy *energy) override
 
double compute_ccsd (TArray &t1, TArray &t2)
 
double compute_rpa_ccsd (TArray &t1, TArray &t2)
 
std::shared_ptr< const EigenVector< typename Tile::numeric_type > > orbital_energy ()
 
void set_orbital_energy (const EigenVector< typename Tile::numeric_type > &orbital_energy)
 
const TArray get_Xab ()
 get three center integral (X|ab) More...
 
const TArray get_Xij ()
 get three center integral (X|ij) More...
 
const TArray get_Xai ()
 get three center integral (X|ai) More...
 
virtual const TArray get_ijab ()
 <ij|ab> More...
 
virtual const TArray get_ijkl ()
 <ij|kl> More...
 
virtual const TArray get_abcd ()
 <ab|cd> More...
 
virtual const TArray get_iabc ()
 <ia|bc> More...
 
virtual const TArray get_iajb ()
 <ia|jb> More...
 
virtual const TArray get_ijka ()
 <ij|ka> More...
 
virtual const TArray get_fock_ia ()
 <i|f|a> More...
 
virtual const TArray get_fock_ij ()
 <i|f|j> More...
 
virtual const TArray get_fock_ab ()
 <a|f|b> More...
 
virtual const TArray get_fock_pq ()
 <p|f|q> More...
 
- Protected Member Functions inherited from mpqc::lcao::LCAOWavefunction< Tile, Policy >
std::tuple< Eigen::VectorXd, TA::DistArray< Tile, Policy > > compute_canonical_occupied_orbitals (bool use_df=true) const
 
std::unique_ptr< Eigen::VectorXd > compute_diagonal_fpq (bool df)
 computes the MO-basis Fock matrix and extracts the diagonal elements More...
 
void init_closed_shell_cabs_svd ()
 
void init_closed_shell_dualbasis_cabs_svd (std::string ri_method)
 

Protected Attributes

TArray T1_
 protected members More...
 
TArray T2_
 
KeyVal kv_
 
std::string solver_str_
 
std::shared_ptr<::mpqc::cc::Solver< TArray > > solver_
 
std::shared_ptr< Wavefunctionref_wfn_
 
AOFactory::DirectTArray direct_ao_array_
 
bool df_
 
bool reduced_abcd_memory_ = false
 
std::string method_
 
std::size_t max_iter_
 
double target_precision_
 
double computed_precision_ = std::numeric_limits<double>::max()
 
bool verbose_
 
double ccsd_corr_energy_
 
bool cp_ccsd_
 
bool cp_ccsd_t_
 
double cp3_precision_
 
double cp4_precision_
 
double cp3_rank_
 
double cp4_rank_
 
bool cp4_on_
 
bool cp_ps_
 
bool cp_robust_
 
double rank_block_factor_
 
int had_block_size_
 
bool ta_als_ = false
 
bool zero_t1_
 
std::string approximation_
 
double alpha_
 
double beta_
 
bool RPA_with_exchange_
 
std::shared_ptr< const EigenVector< typename Tile::numeric_type > > f_pq_diagonal_
 diagonal elements of the Fock matrix (not necessarily the eigenvalues) More...
 

Member Typedef Documentation

◆ AOFactory

template<typename Tile , typename Policy >
using mpqc::lcao::cc::CCSD< Tile, Policy >::AOFactory = gaussian::AOFactory<Tile, Policy>

◆ TArray

template<typename Tile , typename Policy >
using mpqc::lcao::cc::CCSD< Tile, Policy >::TArray = TA::DistArray<Tile, Policy>

Constructor & Destructor Documentation

◆ CCSD() [1/2]

template<typename Tile , typename Policy >
mpqc::lcao::cc::CCSD< Tile, Policy >::CCSD ( )
default

◆ CCSD() [2/2]

template<typename Tile , typename Policy >
mpqc::lcao::cc::CCSD< Tile, Policy >::CCSD ( const KeyVal kv)
inlineexplicit

KeyVal constructor

Parameters
kvthe KeyVal object; it will be queried for all keywords of LCAOWavefunction<Tile,Policy> and the Solver class determined by keyword solver (see below), as well as the following additional keywords:
Keyword Type Default Description
ref Wavefunction none a reference Wavefunction; currently it needs to provide Energy and satisfy requirements for LCAOWavefunction::init_sdref (i.e. provide either CanonicalOrbitalSpace or PopulatedOrbitalSpace)
method string df if df_basis is provided, direct otherwise method to compute the CCSD residual; valid choices are:
  • standard (uses 4-index MO integrals throughout)
  • direct (uses 4-index MO integrals with up to 3 unoccupied indices, and 4-center AO integrals)
  • df (approximates 4-index MO integrals using density fitting)
  • direct_df (hybrid between df and direct that avoids storing MO integrals with 3 unoccupied indices by using DF, see DOI 10.1021/acs.jpca.6b10150 for details)
  • reduced_scaling (approximates 4-index MO integrals using DF and Canonical Polyadic tensor decompositions to reduce the formal scaling of CCSD)
max_iter int 30 maxmium number of solver iterations
solver string jacobi_diis specifies the CCSD solver; valid choices are jacobi_diis (combination of Jacobi update and DIIS) and pno (simulated PNO solver; only valid if method is set to df or direct_df ); kv will also be used to construct the Solver object, hence it will be queried for the corresponding keywords.
verbose bool false if print more information in CCSD iteration
reduced_abcd_memory bool true if method==standard , avoid storing an extra abcd intermediate at the cost of increased FLOPs; if method==df , avoid storage of (ab|cd) integral in favor of lazy evaluation in batches
approximation string none approximation to the residual equations to employ (valid choices are none [no approximation, i.e. CCSD], DCSD, pCCSD, dRPA @ 2CC )
ccsd_gpu bool false if true will compute T2*Wabcd on GPU, available only when TiledArray is enabled with CUDA
cp_ccsd bool false if method==df OR method==reduced_scaling compute Xab integrals using CP decomposition
cp_ccsd_t bool false Only will enable if cp_ccsd==true. Compute CCSD(T) using CP decomposed integrals?
cp3_rank double 3.0 CP3 rank set to number of auxiliary basis functions times cp3_rank
cp4_rank double 0.6 CP4 rank set to number of auxiliary basis functions times cp4_rank
cp3_precision double 0.1 ALS threshold for CP decomposition
cp4_precision double 0.1 ALS threshold for CP decomposition
cp_ps bool false Should the particle particle ladder term be computed using the CP-PS approach?
cp_robust bool false Should the particle particle ladder term be computed using the CP-Robust approach? Isn't enabled if cp_ps==true
rank_block_factor double 1.0 Blocking of the rank dimension is rank_block_factor times the auxiliary block size
had_block_size int -1 New block size for the ab dimensions of the Xab CP decomposition. If new block size is the same as the old block size, no reblocking occurs.
cp4_on_ bool false if method==reduced_scaling Should one use the CP4 section of the code.
RPA_with_exchange bool false if true, the calculation done is ringCCD (or RPA with exchange) and if false direct-ring CCD(drCCD or dRPA)
zero_t1 bool false make T1= 0

◆ ~CCSD()

template<typename Tile , typename Policy >
virtual mpqc::lcao::cc::CCSD< Tile, Policy >::~CCSD ( )
inlinevirtual

Member Function Documentation

◆ can_evaluate()

template<typename Tile , typename Policy >
bool mpqc::lcao::cc::CCSD< Tile, Policy >::can_evaluate ( Energy energy)
inlineoverrideprotected

◆ compute_ccsd()

template<typename Tile , typename Policy >
double mpqc::lcao::cc::CCSD< Tile, Policy >::compute_ccsd ( TArray t1,
TArray t2 
)
protected

function that compute CCSD energy,

Parameters
t1T1 amplitudes, return converged T1 amplitudes
t2T2 amplitudes, return converged T2 amplitudes
Returns
double converged CCSD energy

decide if to initialize Gabcd and Giabc

◆ compute_rpa_ccsd()

template<typename Tile , typename Policy >
double mpqc::lcao::cc::CCSD< Tile, Policy >::compute_rpa_ccsd ( TArray t1,
TArray t2 
)
protected

decide if to initialize Gabcd and Giabc

◆ cp_ccsdt_precision()

template<typename Tile , typename Policy >
double mpqc::lcao::cc::CCSD< Tile, Policy >::cp_ccsdt_precision ( ) const
inline

◆ cp_ccsdt_rank()

template<typename Tile , typename Policy >
double mpqc::lcao::cc::CCSD< Tile, Policy >::cp_ccsdt_rank ( ) const
inline

◆ evaluate()

template<typename Tile , typename Policy >
void mpqc::lcao::cc::CCSD< Tile, Policy >::evaluate ( Energy energy)
overrideprotected

◆ get_abcd()

template<typename Tile , typename Policy >
virtual const TArray mpqc::lcao::cc::CCSD< Tile, Policy >::get_abcd ( )
inlineprotectedvirtual

<ab|cd>

◆ get_direct_ao_integral()

template<typename Tile , typename Policy >
const AOFactory::DirectTArray& mpqc::lcao::cc::CCSD< Tile, Policy >::get_direct_ao_integral ( ) const
inline

◆ get_fock_ab()

template<typename Tile , typename Policy >
virtual const TArray mpqc::lcao::cc::CCSD< Tile, Policy >::get_fock_ab ( )
inlineprotectedvirtual

<a|f|b>

◆ get_fock_ia()

template<typename Tile , typename Policy >
virtual const TArray mpqc::lcao::cc::CCSD< Tile, Policy >::get_fock_ia ( )
inlineprotectedvirtual

<i|f|a>

◆ get_fock_ij()

template<typename Tile , typename Policy >
virtual const TArray mpqc::lcao::cc::CCSD< Tile, Policy >::get_fock_ij ( )
inlineprotectedvirtual

<i|f|j>

◆ get_fock_pq()

template<typename Tile , typename Policy >
virtual const TArray mpqc::lcao::cc::CCSD< Tile, Policy >::get_fock_pq ( )
inlineprotectedvirtual

<p|f|q>

◆ get_iabc()

template<typename Tile , typename Policy >
virtual const TArray mpqc::lcao::cc::CCSD< Tile, Policy >::get_iabc ( )
inlineprotectedvirtual

<ia|bc>

◆ get_iajb()

template<typename Tile , typename Policy >
virtual const TArray mpqc::lcao::cc::CCSD< Tile, Policy >::get_iajb ( )
inlineprotectedvirtual

<ia|jb>

◆ get_ijab()

template<typename Tile , typename Policy >
virtual const TArray mpqc::lcao::cc::CCSD< Tile, Policy >::get_ijab ( )
inlineprotectedvirtual

<ij|ab>

◆ get_ijka()

template<typename Tile , typename Policy >
virtual const TArray mpqc::lcao::cc::CCSD< Tile, Policy >::get_ijka ( )
inlineprotectedvirtual

<ij|ka>

◆ get_ijkl()

template<typename Tile , typename Policy >
virtual const TArray mpqc::lcao::cc::CCSD< Tile, Policy >::get_ijkl ( )
inlineprotectedvirtual

<ij|kl>

◆ get_Xab()

template<typename Tile , typename Policy >
const TArray mpqc::lcao::cc::CCSD< Tile, Policy >::get_Xab ( )
inlineprotected

get three center integral (X|ab)

utility functions that computes integrals

◆ get_Xai()

template<typename Tile , typename Policy >
const TArray mpqc::lcao::cc::CCSD< Tile, Policy >::get_Xai ( )
inlineprotected

get three center integral (X|ai)

◆ get_Xij()

template<typename Tile , typename Policy >
const TArray mpqc::lcao::cc::CCSD< Tile, Policy >::get_Xij ( )
inlineprotected

get three center integral (X|ij)

◆ is_cp_ccsdt()

template<typename Tile , typename Policy >
bool mpqc::lcao::cc::CCSD< Tile, Policy >::is_cp_ccsdt ( ) const
inline

◆ is_df()

template<typename Tile , typename Policy >
bool mpqc::lcao::cc::CCSD< Tile, Policy >::is_df ( ) const
inline

◆ is_reduced_abcd_memory()

template<typename Tile , typename Policy >
bool mpqc::lcao::cc::CCSD< Tile, Policy >::is_reduced_abcd_memory ( ) const
inline

◆ obsolete()

template<typename Tile , typename Policy >
void mpqc::lcao::cc::CCSD< Tile, Policy >::obsolete ( )
inlineoverride

◆ orbital_energy()

template<typename Tile , typename Policy >
std::shared_ptr<const EigenVector<typename Tile::numeric_type> > mpqc::lcao::cc::CCSD< Tile, Policy >::orbital_energy ( )
inlineprotected

Accesses the orbital energies, i.e. the diagonal elements of the Fock matrix

Returns
the diagonal elements of the Fock matrix
Exceptions
ProgrammingErrorif orbitals were localized, hence the Fock matrix is not diagonal

◆ set_orbital_energy()

template<typename Tile , typename Policy >
void mpqc::lcao::cc::CCSD< Tile, Policy >::set_orbital_energy ( const EigenVector< typename Tile::numeric_type > &  orbital_energy)
inlineprotected

◆ set_t1()

template<typename Tile , typename Policy >
void mpqc::lcao::cc::CCSD< Tile, Policy >::set_t1 ( const TArray t1)
inline

◆ set_t2()

template<typename Tile , typename Policy >
void mpqc::lcao::cc::CCSD< Tile, Policy >::set_t2 ( const TArray t2)
inline

◆ set_target_precision()

template<typename Tile , typename Policy >
void mpqc::lcao::cc::CCSD< Tile, Policy >::set_target_precision ( double  precision)
inline

◆ t1()

template<typename Tile , typename Policy >
TArray mpqc::lcao::cc::CCSD< Tile, Policy >::t1 ( ) const
inline

◆ t2()

template<typename Tile , typename Policy >
TArray mpqc::lcao::cc::CCSD< Tile, Policy >::t2 ( ) const
inline

◆ verbose()

template<typename Tile , typename Policy >
bool mpqc::lcao::cc::CCSD< Tile, Policy >::verbose ( ) const
inline

Member Data Documentation

◆ alpha_

template<typename Tile , typename Policy >
double mpqc::lcao::cc::CCSD< Tile, Policy >::alpha_
protected

◆ approximation_

template<typename Tile , typename Policy >
std::string mpqc::lcao::cc::CCSD< Tile, Policy >::approximation_
protected

◆ beta_

template<typename Tile , typename Policy >
double mpqc::lcao::cc::CCSD< Tile, Policy >::beta_
protected

◆ ccsd_corr_energy_

template<typename Tile , typename Policy >
double mpqc::lcao::cc::CCSD< Tile, Policy >::ccsd_corr_energy_
protected

◆ computed_precision_

template<typename Tile , typename Policy >
double mpqc::lcao::cc::CCSD< Tile, Policy >::computed_precision_ = std::numeric_limits<double>::max()
protected

◆ cp3_precision_

template<typename Tile , typename Policy >
double mpqc::lcao::cc::CCSD< Tile, Policy >::cp3_precision_
protected

◆ cp3_rank_

template<typename Tile , typename Policy >
double mpqc::lcao::cc::CCSD< Tile, Policy >::cp3_rank_
protected

◆ cp4_on_

template<typename Tile , typename Policy >
bool mpqc::lcao::cc::CCSD< Tile, Policy >::cp4_on_
protected

◆ cp4_precision_

template<typename Tile , typename Policy >
double mpqc::lcao::cc::CCSD< Tile, Policy >::cp4_precision_
protected

◆ cp4_rank_

template<typename Tile , typename Policy >
double mpqc::lcao::cc::CCSD< Tile, Policy >::cp4_rank_
protected

◆ cp_ccsd_

template<typename Tile , typename Policy >
bool mpqc::lcao::cc::CCSD< Tile, Policy >::cp_ccsd_
protected

◆ cp_ccsd_t_

template<typename Tile , typename Policy >
bool mpqc::lcao::cc::CCSD< Tile, Policy >::cp_ccsd_t_
protected

◆ cp_ps_

template<typename Tile , typename Policy >
bool mpqc::lcao::cc::CCSD< Tile, Policy >::cp_ps_
protected

◆ cp_robust_

template<typename Tile , typename Policy >
bool mpqc::lcao::cc::CCSD< Tile, Policy >::cp_robust_
protected

◆ df_

template<typename Tile , typename Policy >
bool mpqc::lcao::cc::CCSD< Tile, Policy >::df_
protected

◆ direct_ao_array_

template<typename Tile , typename Policy >
AOFactory::DirectTArray mpqc::lcao::cc::CCSD< Tile, Policy >::direct_ao_array_
protected

◆ f_pq_diagonal_

template<typename Tile , typename Policy >
std::shared_ptr<const EigenVector<typename Tile::numeric_type> > mpqc::lcao::cc::CCSD< Tile, Policy >::f_pq_diagonal_
protected

diagonal elements of the Fock matrix (not necessarily the eigenvalues)

◆ had_block_size_

template<typename Tile , typename Policy >
int mpqc::lcao::cc::CCSD< Tile, Policy >::had_block_size_
protected

◆ kv_

template<typename Tile , typename Policy >
KeyVal mpqc::lcao::cc::CCSD< Tile, Policy >::kv_
protected

◆ max_iter_

template<typename Tile , typename Policy >
std::size_t mpqc::lcao::cc::CCSD< Tile, Policy >::max_iter_
protected

◆ method_

template<typename Tile , typename Policy >
std::string mpqc::lcao::cc::CCSD< Tile, Policy >::method_
protected

◆ rank_block_factor_

template<typename Tile , typename Policy >
double mpqc::lcao::cc::CCSD< Tile, Policy >::rank_block_factor_
protected

◆ reduced_abcd_memory_

template<typename Tile , typename Policy >
bool mpqc::lcao::cc::CCSD< Tile, Policy >::reduced_abcd_memory_ = false
protected

◆ ref_wfn_

template<typename Tile , typename Policy >
std::shared_ptr<Wavefunction> mpqc::lcao::cc::CCSD< Tile, Policy >::ref_wfn_
protected

◆ RPA_with_exchange_

template<typename Tile , typename Policy >
bool mpqc::lcao::cc::CCSD< Tile, Policy >::RPA_with_exchange_
protected

◆ solver_

template<typename Tile , typename Policy >
std::shared_ptr<::mpqc::cc::Solver<TArray> > mpqc::lcao::cc::CCSD< Tile, Policy >::solver_
protected

◆ solver_str_

template<typename Tile , typename Policy >
std::string mpqc::lcao::cc::CCSD< Tile, Policy >::solver_str_
protected

◆ T1_

template<typename Tile , typename Policy >
TArray mpqc::lcao::cc::CCSD< Tile, Policy >::T1_
protected

protected members

◆ T2_

template<typename Tile , typename Policy >
TArray mpqc::lcao::cc::CCSD< Tile, Policy >::T2_
protected

◆ ta_als_

template<typename Tile , typename Policy >
bool mpqc::lcao::cc::CCSD< Tile, Policy >::ta_als_ = false
protected

◆ target_precision_

template<typename Tile , typename Policy >
double mpqc::lcao::cc::CCSD< Tile, Policy >::target_precision_
protected

◆ verbose_

template<typename Tile , typename Policy >
bool mpqc::lcao::cc::CCSD< Tile, Policy >::verbose_
protected

◆ zero_t1_

template<typename Tile , typename Policy >
bool mpqc::lcao::cc::CCSD< Tile, Policy >::zero_t1_
protected

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