MPQC  3.0.0-alpha
sc::SingleReference_R12Intermediates< T > Class Template Reference

SingleReference_R12Intermediates computes R12/F12 intermediates using MPQC3 runtime. More...

#include <chemistry/qc/mbptr12/sr_r12intermediates.h>

Public Types

typedef TA::Array< T, 4 > TArray4
 standard 4-index tensor
 
typedef TA::Array< T, 4, DA4_Tile< T > > TArray4d
 4-index tensor with lazy tiles
 
typedef TA::Array< T, 2 > TArray2
 standard 2-index tensor
 
typedef TA::Array< T, 2, DA4_Tile< T > > TArray2d
 2-index tensor with lazy tiles
 
typedef TA::Array< TA::Tensor< T >, 2 > TArray22
 2-index tensor with lazy tiles More...
 
typedef TA::Array< TA::Tensor< T >, 2, DA4_Tile34< T > > TArray22d
 2-index tensor of lazy 2-index tensors
 
typedef TA::Array< T, 4, LazyTensor< T, 4, expressions::TGeminalGenerator< T > > > TArray4Tg
 geminal T tensor
 

Public Member Functions

 SingleReference_R12Intermediates (madness::World &world, const Ref< R12WavefunctionWorld > &r12world)
 Constructs an SingleReference_R12Intermediates object. More...
 
const Ref< R12WavefunctionWorld > & r12world () const
 
std::pair< TArray2, TArray2V_diag ()
 computes diagonal (spin-restricted, for now) V intermediate More...
 
std::pair< TArray2, TArray2X_diag ()
 computes diagonal (spin-restricted, for now) X intermediate More...
 
std::pair< TArray2, TArray2B_diag ()
 computes diagonal (spin-restricted, for now) B intermediate More...
 
void gf2_r12 (int orbital)
 Computes second-order Green's function IPs and EAs \parame orbital the index of the orbital, -1 = HOMO, +1 = LUMO.
 
TArray2 rdm1 ()
 returns the 1-particle reduced density matrix More...
 
void compute_multipole ()
 
TArray2 XaiAddToXam (const TA::Array< double, 2 > &Xam, const TA::Array< double, 2 > &Xai)
 
TArray2 BPk_Qk (const char *p, const char *q, const double C_0, const double C_1)
 
TArray4 Bpr_qs (const char *p, const char *q)
 
TArray4 VPq_Rs (const char *p, const char *q, const char *r, const char *s, const double C_0, const double C_1)
 
TArray2 VRk_Sk (const char *r, const char *s, const double C_0, const double C_1)
 
TArray2 Xam_CabsSingles (const TArray2 &TmA, const TArray2 &Tma)
 
TArray2 Xam_mp2 (const TArray4 &T2_ijab, const TArray2 &Dij, const TArray2 &Dab)
 
TArray2 Xam_Cmp2f12 (const double C_0, const double C_1, const TArray4 &T2_ijab, const TArray4 &A_ijab, const TArray2 &Dij, const TArray2 &Dab, const TArray2 &RT_apb)
 
void compute_Df12_XB (const double C_0, const double C_1, TArray2 &D_f12_ij, TArray2 &D_f12_ab, TArray2 &D_f12_apbp, TArray2 &D_f12_apb)
 
TArray2 Xam_Df12_XB (const double C_0, const double C_1, const TArray2 &Df12_ij, const TArray2 &Df12_ab, const TArray2 &Df12_apbp, const TArray2 &Df12_apb)
 
TArray2 Xam_V (const double C_0, const double C_1)
 
TArray2 Xam_X (const double C_0, const double C_1)
 
TArray2 Xam_B (const double C_0, const double C_1)
 
TArray2 Xiip_VBX (const double C_0, const double C_1)
 
TArray2 Xiip_CVT (const double C_0, const double C_1, const TArray2 &T1, const TArray4 &T2)
 
TArray2 Xam_CT2_ccsd (const double C_0, const double C_1, const TArray4 &T2, const TArray2 &RT2_aPb)
 
TArray2 Xam_VT_ccsd (const double C_0, const double C_1, const TArray2 &T1, const TArray4 &T2)
 
void compute_T_cc2 (TArray2 &T1, TArray4 &T2)
 
void compute_lambda_cc2 (const TArray2 &t1, const TArray4 &t2, TArray2 &L1, TArray4 &L2)
 
void compute_lambda_cc2_2 (const TArray2 &t1, const TArray4 &t2, TArray2 &L1, TArray4 &L2)
 
void compute_Delta_cc (const TA::TiledRange &TR_ai, const TA::TiledRange &TR_abij, TArray2 &D_ai, TArray4 &D_abij)
 
void compute_intermediates_TFW_ccsd (const TArray2 &t1, const TArray4 &t2, const TArray4d &g_abij, const TArray4d &g_aikl, const TArray4d &g_aibj, const TArray4d &g_aijb, const TArray4d &g_abci, const TArray4d &g_ijkl, const TArray4d &g_abcd, TArray2 &TFkc, TArray2 &TFac, TArray2 &TFki, TArray4 &TW_KbCj_ab, TArray4 &TW_KbcJ_ba, TArray4 &TW_AbCd_ab, TArray4 &TW_KlIj_ab)
 
void compute_intermediates_CW_ccsd (const TArray2 &t1, const TArray4 &t2, TArray2 &CFkc, TArray2 &CFac, TArray2 &CFki, TArray4 &CW_KbEj_ab, TArray4 &CW_KbeJ_ba, TArray4 &CW_AbCd_ab, TArray4 &CW_AbCi_ab, TArray4 &CW_KlIj_ab, TArray4 &CW_KbIj_ab, TArray4 &CW_KlIc_ab, TArray4 &CW_KliC_ab, TArray4 &CW_AkCd_ab)
 
void compute_T_ccsd (TArray2 &t1, TArray4 &t2, const std::string method)
 
void compute_lambda_ccsd (const TArray2 &t1, const TArray4 &t2, TArray2 &L1, TArray4 &L2, const std::string method)
 
void compute_cc2_1rdm_amp (const TArray2 &T1_cc2, const TArray4 &T2_cc2, const TArray2 &L1_cc2, const TArray4 &L2_cc2, TArray2 &Dij_cc2, TArray2 &Dab_cc2, TArray2 &Dia_cc2, TArray2 &Dai_cc2)
 
void compute_Gamma_ijab_ccsd (const TArray2 &T1, const TArray4 &T2, const TArray4 &tau_ab, const TArray2 &L1, const TArray4 &L2, TArray4 &Gamma_IjAb_ab)
 
void compute_ccsd_1rdm_amp (const TArray2 &T1, const TArray4 &T2, const TArray2 &L1, const TArray4 &L2, TArray2 &Dij, TArray2 &Dab, TArray2 &Dia, TArray2 &Dai)
 
void compute_Gamma2_ccsd (const TArray2 &T1, const TArray4 &T2, const TArray2 &L1, const TArray4 &L2, const TArray4 &tau_aa, const TArray4 &tau_ab, TArray4 &Gamma_IjKa_ab, TArray4 &Gamma_AbCi_ab, TArray4 &Gamma_iBjA_ab, TArray4 &Gamma_iBJa_ba, TArray4 &Gamma_AbCd_ab, TArray4 &Gamma_IjKl_ab, TArray4 &Gamma_IjAb_ab)
 
void compute_Xam_ccsd (const TArray2 &T1, const TArray4 &T2, const TArray2 &L1, const TArray4 &L2, TArray2 &Xam_tot, TArray2 &Xiip)
 
void compute_multipole_F12b_coupling ()
 
TArray4 VPQ_RS (const char *p, const char *q, const char *r, const char *s)
 
TArray2 Xam_CT2L2_f12b (const double C_0, const double C_1, const TArray4 &T2, const TArray2 &RT2_aPb, const TArray4 &L2, const TArray2 &RL2_aPb)
 
TArray2 Xam_VTL_f12b (const double C_0, const double C_1, const TArray2 &T1, const TArray4 &T2, const TArray2 &L1, const TArray4 &L2)
 
TArray2 Xam_VT1T1_f12b (const double C_0, const double C_1, const TArray2 &T1)
 
TArray2 Xam_VT1T1_f12b_test (const double C_0, const double C_1, const TArray4 &tauT1_ab)
 
TArray2 Xam_VL2T1_f12b (const double C_0, const double C_1, const TArray2 &T1, const TArray4 &L2)
 
TArray2 Xam_VL2T1_f12b_test (const double C_0, const double C_1, const TArray2 &T1, const TArray4 &L2)
 
TArray2 Xam_VL2T1T1_f12b (const double C_0, const double C_1, const TArray2 &T1, const TArray4 &L2)
 
TArray2 Xam_VL2T1T1_f12b_test (const double C_0, const double C_1, const TArray2 &T1, const TArray4 &L2)
 
TArray2 Xiip_CVT_f12b (const double C_0, const double C_1, const TArray2 &T1, const TArray4 &T2, const TArray2 &L1, const TArray4 &L2)
 
TArray2 Xiip_VT1T1_f12b (const double C_0, const double C_1, const TArray2 &T1, const TArray4 &T2, const TArray4 &L2)
 
TArray4 rdm2 ()
 returns the 2-particle density matrix More...
 
TArray4 V_spinfree (bool symmetrize_p1_p2=false)
 computes spin-free V intermediate More...
 
TArray4 X_spinfree (bool symmetrize_p1_p2=false)
 computes spin-free X intermediate More...
 
TArray4 B_spinfree (bool symmetrize_p1_p2=false)
 computes spin-free B intermediate More...
 
void set_T1 (const RefSCMatrix &t1)
 provides T1 amplitude tensor More...
 
void set_T1_cabs (const RefSCMatrix &t1_cabs)
 provides T1 CABS amplitude tensor More...
 
void set_L1 (const RefSCMatrix &l1)
 provides Lambda1 amplitude tensor More...
 
void set_T2 (const Ref< DistArray4 >(&t2)[NSpinCases2])
 provides T2 amplitudes More...
 
void set_L2 (const Ref< DistArray4 >(&l2)[NSpinCases2])
 provides Lambda2 amplitudes More...
 
void set_rdm2 (const Ref< SpinFreeRDM< Two > > &rdm2)
 provides (spin-free) RDM2 More...
 
TArray4dijxy (const std::string &key)
 
TArray22dij_xy (const std::string &key)
 
TArray2xy (const std::string &key)
 
TArray2sieve (const TArray2 &input, const std::string &output_annotation)
 sieves x|o1|y -> x'|o1|y' does not throw only if each tile of the result depends only on 1 tile of the input
 
TA::expressions::TsrExpr< const TArray4d, true > _4 (const std::string &key)
 Given a descriptive key, creates a rank-4 Array of integrals, or other related quantities The syntax of key is similar to that used by ParsedTwoBodyInt and TwoBodyMOIntsRuntime, but with te_type embedded into key. More...
 
TA::expressions::TsrExpr< const TArray2, true > _2 (const std::string &key)
 Given a descriptive key, creates a rank-2 Array of integrals, or other related quantities The syntax of key is similar to that used by ParsedOneBodyInt, but with oe_type embedded into key. More...
 
TA::expressions::TsrExpr< const TArray4Tg, true > _Tg (const std::string &key)
 like _4, produces geminal T tensor
 

Detailed Description

template<typename T>
class sc::SingleReference_R12Intermediates< T >

SingleReference_R12Intermediates computes R12/F12 intermediates using MPQC3 runtime.

Template Parameters
Tthe numeric type supporting all tensors

Member Typedef Documentation

◆ TArray22

template<typename T>
typedef TA::Array<TA::Tensor<T>, 2> sc::SingleReference_R12Intermediates< T >::TArray22

2-index tensor with lazy tiles

2-index tensor of 2-index tensors

Constructor & Destructor Documentation

◆ SingleReference_R12Intermediates()

template<typename T>
sc::SingleReference_R12Intermediates< T >::SingleReference_R12Intermediates ( madness::World &  world,
const Ref< R12WavefunctionWorld > &  r12world 
)
inline

Constructs an SingleReference_R12Intermediates object.

There can be many.

Parameters
worldMADNESS world in which the computed Arrays will live
r12worldR12WavefunctionWorld that provides computation of integrals and OrbitalSpace objects

Member Function Documentation

◆ _2()

template<typename T >
TA::expressions::TsrExpr< const typename SingleReference_R12Intermediates< T >::TArray2, true > sc::SingleReference_R12Intermediates< T >::_2 ( const std::string &  key)

Given a descriptive key, creates a rank-2 Array of integrals, or other related quantities The syntax of key is similar to that used by ParsedOneBodyInt, but with oe_type embedded into key.

The following oe_types are understood:

  • mu_i : electric dipole integral $ - \bf{r}_i $ (i = x,y,z)
  • q_ij : electric quadrupole integral $ - \bf{r}_i \bf{r}_j $ (ij = xx, xy, xz, yy, yz, zz)
  • gamma : $ \Gamma $, 1-RDM of the reference wave function
  • T1 : 1-body excitation amplitudes, e.g. of the coupled-cluster method, provided by set_T1() and/or set_T1_cabs()
  • L1 : 1-body de-excitation amplitudes, e.g. of the coupled-cluster method, provided by set_L1()
  • I : identity matrix

Indices are interpreted using to_space().

Example1: <i|mu_x|p> will create an array of <act_occ| -x |obs> integrals Example2: <j_F(p')|q_xy|a'> will create an array of <cbs| -xy |cabs> integrals with the second index transformed using Fock matrix between act_occ and cbs spaces.

Note that the indices are used to create the resulting tensor expression, hence these keys can be used in composing expressions.

See also
_4()
xy()

◆ _4()

template<typename T >
TA::expressions::TsrExpr< const typename SingleReference_R12Intermediates< T >::TArray4d, true > sc::SingleReference_R12Intermediates< T >::_4 ( const std::string &  key)

Given a descriptive key, creates a rank-4 Array of integrals, or other related quantities The syntax of key is similar to that used by ParsedTwoBodyInt and TwoBodyMOIntsRuntime, but with te_type embedded into key.

The following te_types are understood:

  • g : $ r_{12}^{-1} $
  • r : $ f(r_{12}) $
  • gr : $ r_{12}^{-1} f(r_{12}) $
  • rTr : $ [f(r_{12}), [ \hat{T}_1 , f(r_{12})]] $
  • gamma : $ \Gamma $, 2-RDM, provided by set_rdm2()
  • T2 : 2-body excitation amplitudes, e.g. of the coupled-cluster method, provided by set_T2()
  • L2 : 2-body de-excitation amplitudes, e.g. of the coupled-cluster method, provided by set_L2()

Indices are interpreted using to_space().

Example1: <i j|g|p a'> will create an array of <act_occ act_occ| r_{12}^{-1} |obs cabs> integrals Example2: <i j_F(p')|g|p a'> will create an array of <act_occ cbs| r_{12}^{-1} |obs cabs> integrals with the second index transformed using Fock matrix between act_occ and cbs spaces.

Note that the indices are used to create the resulting tensor expression, hence these keys can be used in composing expressions. For example, the full (non-diagonal) V intermediate of the R12 method without the CABS terms can be computed as follows: TArray4d V_ij_kl = _4("<i j|gr|k l>") - _4("<i j|g|p q>") * _4("<k l|r|p q>");

See also
ijxy()

◆ B_diag()

template<typename T >
std::pair< typename SingleReference_R12Intermediates< T >::TArray2, typename SingleReference_R12Intermediates< T >::TArray2 > sc::SingleReference_R12Intermediates< T >::B_diag ( )

computes diagonal (spin-restricted, for now) B intermediate

Returns
$ B^{ij}_{ij} $ and $ B^{ij}_{ji} $, respectively

this is incomplete at the moment

◆ B_spinfree()

template<typename T >
SingleReference_R12Intermediates< T >::TArray4 sc::SingleReference_R12Intermediates< T >::B_spinfree ( bool  symmetrize_p1_p2 = false)

computes spin-free B intermediate

Parameters
symmetrize_p1_p2if yes, make sure this holds for the result: B^{ij}_{kl} = B^{ji}_{kl}
Returns
$ B^{ij}_{kl} $

◆ compute_multipole()

template<typename T >
void sc::SingleReference_R12Intermediates< T >::compute_multipole ( )

< evaluator of preconditioner

References sc::Molecule::natom(), and sc::ExEnv::out0().

◆ compute_multipole_F12b_coupling()

template<typename T >
void sc::SingleReference_R12Intermediates< T >::compute_multipole_F12b_coupling ( )

< evaluator of preconditioner

References sc::ExEnv::out0().

◆ ij_xy()

◆ ijxy()

◆ rdm1()

returns the 1-particle reduced density matrix

Returns
$ \gamma^{p}_{q} $, respectively

this is just an example of how to compute the density

◆ rdm2()

template<typename T>
TArray4 sc::SingleReference_R12Intermediates< T >::rdm2 ( )

returns the 2-particle density matrix

Returns
$ \gamma^{pq}_{rs} $, respectively

Referenced by sc::SingleReference_R12Intermediates< T >::set_rdm2().

◆ set_L1()

template<typename T>
void sc::SingleReference_R12Intermediates< T >::set_L1 ( const RefSCMatrix l1)
inline

provides Lambda1 amplitude tensor

Parameters
l1act_occ by act_vir matrix

◆ set_L2()

template<typename T>
void sc::SingleReference_R12Intermediates< T >::set_L2 ( const Ref< DistArray4 >(&)  l2[NSpinCases2])
inline

provides Lambda2 amplitudes

Parameters
l2array of Lambda2 amplitudes, for AlphaBeta, AlphaAlpha, and (optionally) BetaBeta

◆ set_rdm2()

template<typename T>
void sc::SingleReference_R12Intermediates< T >::set_rdm2 ( const Ref< SpinFreeRDM< Two > > &  rdm2)
inline

provides (spin-free) RDM2

Parameters
rdm2a SpinFreeRDM<Two> object

References sc::SingleReference_R12Intermediates< T >::rdm2().

◆ set_T1()

template<typename T>
void sc::SingleReference_R12Intermediates< T >::set_T1 ( const RefSCMatrix t1)
inline

provides T1 amplitude tensor

Parameters
t1act_occ by act_vir matrix

◆ set_T1_cabs()

template<typename T>
void sc::SingleReference_R12Intermediates< T >::set_T1_cabs ( const RefSCMatrix t1_cabs)
inline

provides T1 CABS amplitude tensor

Parameters
t1act_occ by allvir (or CABS) matrix

◆ set_T2()

template<typename T>
void sc::SingleReference_R12Intermediates< T >::set_T2 ( const Ref< DistArray4 >(&)  t2[NSpinCases2])
inline

provides T2 amplitudes

Parameters
t2array of T2 amplitudes, for AlphaBeta, AlphaAlpha, and (optionally) BetaBeta

◆ V_diag()

template<typename T >
std::pair< typename SingleReference_R12Intermediates< T >::TArray2, typename SingleReference_R12Intermediates< T >::TArray2 > sc::SingleReference_R12Intermediates< T >::V_diag ( )

computes diagonal (spin-restricted, for now) V intermediate

Returns
$ V^{ij}_{ij} $ and $ V^{ij}_{ji} $, respectively

◆ V_spinfree()

template<typename T >
SingleReference_R12Intermediates< T >::TArray4 sc::SingleReference_R12Intermediates< T >::V_spinfree ( bool  symmetrize_p1_p2 = false)

computes spin-free V intermediate

Parameters
symmetrize_p1_p2if yes, make sure this holds for the result: V^{ij}_{mn} = V^{ji}_{nm}
Returns
$ V^{ij}_{mn} $, where ij refer to the geminal indices

◆ X_diag()

template<typename T >
std::pair< typename SingleReference_R12Intermediates< T >::TArray2, typename SingleReference_R12Intermediates< T >::TArray2 > sc::SingleReference_R12Intermediates< T >::X_diag ( )

computes diagonal (spin-restricted, for now) X intermediate

Returns
$ X^{ij}_{ij} $ and $ X^{ij}_{ji} $, respectively

◆ X_spinfree()

template<typename T >
SingleReference_R12Intermediates< T >::TArray4 sc::SingleReference_R12Intermediates< T >::X_spinfree ( bool  symmetrize_p1_p2 = false)

computes spin-free X intermediate

Parameters
symmetrize_p1_p2if yes, make sure this holds for the result: X^{ij}_{kl} = X^{ji}_{kl}
Returns
$ X^{ij}_{kl} $

◆ xy()

template<typename T >
SingleReference_R12Intermediates< T >::TArray2 & sc::SingleReference_R12Intermediates< T >::xy ( const std::string &  key)

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

Generated at Sun Jan 26 2020 23:24:12 for MPQC 3.0.0-alpha using the documentation package Doxygen 1.8.16.