MPQC 3.0.0-alpha
Loading...
Searching...
No Matches
r12int_eval.h
1//
2// r12int_eval.h
3//
4// Copyright (C) 2004 Edward Valeev
5//
6// Author: Edward Valeev <evaleev@vt.edu>
7// Maintainer: EV
8//
9// This file is part of the SC Toolkit.
10//
11// The SC Toolkit is free software; you can redistribute it and/or modify
12// it under the terms of the GNU Library General Public License as published by
13// the Free Software Foundation; either version 2, or (at your option)
14// any later version.
15//
16// The SC Toolkit is distributed in the hope that it will be useful,
17// but WITHOUT ANY WARRANTY; without even the implied warranty of
18// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
19// GNU Library General Public License for more details.
20//
21// You should have received a copy of the GNU Library General Public License
22// along with the SC Toolkit; see the file COPYING.LIB. If not, write to
23// the Free Software Foundation, 675 Mass Ave, Cambridge, MA 02139, USA.
24//
25// The U.S. Government is granted a limited license as per AL 91-7.
26//
27
28#ifndef _chemistry_qc_mbptr12_r12inteval_h
29#define _chemistry_qc_mbptr12_r12inteval_h
30
31#include <util/ref/ref.h>
32#include <chemistry/qc/mbptr12/r12wfnworld.h>
33#include <chemistry/qc/mbptr12/r12technology.h>
34#include <chemistry/qc/mbptr12/twoparticlecontraction.h>
35#include <chemistry/qc/wfn/spin.h>
36#include <chemistry/qc/mbptr12/fixedcoefficient.h>
37#include <chemistry/qc/mbptr12/twobodytensorinfo.h>
38
39namespace sc {
40
41 class TwoBodyMOIntsTransform;
42 class DistArray4;
43 class R12Amplitudes;
44 class R12EnergyIntermediates;
45
50class R12IntEval : virtual public SavableState {
51 private:
52 // change to false to use the old fock builder
53 static const bool USE_FOCKBUILD = true;
54
55
56 bool evaluated_;
57
58 // Calculation information (number of basis functions, R12 approximation, etc.)
60
61 RefSCMatrix V_[NSpinCases2];
62 RefSCMatrix X_[NSpinCases2];
63 // Note that intermediate B is symmetric but is stored as a full matrix to simplify the code
64 // that computes asymmetric form of B
65 RefSCMatrix B_[NSpinCases2];
66 RefSCMatrix BB_[NSpinCases2]; // The difference between B intermediate of approximation B and A'
67 RefSCMatrix A_[NSpinCases2];
69 Ref<CuspConsistentGeminalCoefficient> cuspconsistentgeminalcoefficient_[NSpinCases2];
70
71 RefSCVector emp2pair_[NSpinCases2];
72 RefSCDimension dim_oo_[NSpinCases2];
73 RefSCDimension dim_vv_[NSpinCases2];
75 RefSCDimension dim_aa_[NSpinCases2];
77 RefSCDimension dim_f12_[NSpinCases2];
79 RefSCDimension dim_GG_[NSpinCases2];
81 RefSCDimension dim_gg_[NSpinCases2];
82
83 Ref<R12Amplitudes> Amps_; // First-order amplitudes of various contributions to the pair functions
84 double emp2_obs_singles_;
85 double emp2_cabs_singles_;
86 RefSCMatrix T1_cabs_[NSpinCases1];
87 int debug_;
88
89 mutable RefSymmSCMatrix ordm_[NSpinCases1];
90
92 void spinadapt_mospace_labels(SpinCase1 spin, std::string& id, std::string& name) const;
93
96 void f_bra_ket(SpinCase1 spin,
97 bool make_F,
98 bool make_hJ,
99 bool make_K,
103 const Ref<OrbitalSpace>& extspace,
104 const Ref<OrbitalSpace>& intspace
105 );
106 Ref<OrbitalSpace> hj_i_p_[NSpinCases1];
107 Ref<OrbitalSpace> hj_i_A_[NSpinCases1];
108 Ref<OrbitalSpace> hj_i_P_[NSpinCases1];
109 Ref<OrbitalSpace> hj_i_m_[NSpinCases1];
110 Ref<OrbitalSpace> hj_i_a_[NSpinCases1];
111 Ref<OrbitalSpace> hj_m_m_[NSpinCases1];
112 Ref<OrbitalSpace> hj_m_p_[NSpinCases1];
113 Ref<OrbitalSpace> hj_a_A_[NSpinCases1];
114 Ref<OrbitalSpace> hj_p_p_[NSpinCases1];
115 Ref<OrbitalSpace> hj_p_A_[NSpinCases1];
116 Ref<OrbitalSpace> hj_p_P_[NSpinCases1];
117 Ref<OrbitalSpace> hj_P_P_[NSpinCases1];
118 Ref<OrbitalSpace> hj_p_m_[NSpinCases1];
119 Ref<OrbitalSpace> hj_p_a_[NSpinCases1];
120 Ref<OrbitalSpace> K_i_p_[NSpinCases1];
121 Ref<OrbitalSpace> K_i_m_[NSpinCases1];
122 Ref<OrbitalSpace> K_i_a_[NSpinCases1];
123 Ref<OrbitalSpace> K_i_A_[NSpinCases1];
124 Ref<OrbitalSpace> K_i_P_[NSpinCases1];
125 Ref<OrbitalSpace> K_m_a_[NSpinCases1];
126 Ref<OrbitalSpace> K_a_a_[NSpinCases1];
127 Ref<OrbitalSpace> K_a_p_[NSpinCases1];
128 Ref<OrbitalSpace> K_a_P_[NSpinCases1];
129 Ref<OrbitalSpace> K_p_p_[NSpinCases1];
130 Ref<OrbitalSpace> K_p_m_[NSpinCases1];
131 Ref<OrbitalSpace> K_p_a_[NSpinCases1];
132 Ref<OrbitalSpace> K_p_A_[NSpinCases1];
133 Ref<OrbitalSpace> K_p_P_[NSpinCases1];
134 Ref<OrbitalSpace> K_A_P_[NSpinCases1];
135 Ref<OrbitalSpace> K_P_P_[NSpinCases1];
136 Ref<OrbitalSpace> F_P_P_[NSpinCases1];
137 Ref<OrbitalSpace> F_p_A_[NSpinCases1];
138 Ref<OrbitalSpace> F_p_p_[NSpinCases1];
139 Ref<OrbitalSpace> F_p_m_[NSpinCases1];
140 Ref<OrbitalSpace> F_p_a_[NSpinCases1];
141 Ref<OrbitalSpace> F_m_m_[NSpinCases1];
142 Ref<OrbitalSpace> F_m_a_[NSpinCases1];
143 Ref<OrbitalSpace> F_m_p_[NSpinCases1];
144 Ref<OrbitalSpace> F_m_P_[NSpinCases1];
145 Ref<OrbitalSpace> F_gg_P_[NSpinCases1];
146 Ref<OrbitalSpace> F_m_A_[NSpinCases1];
147 Ref<OrbitalSpace> F_i_A_[NSpinCases1];
148 Ref<OrbitalSpace> F_i_m_[NSpinCases1];
149 Ref<OrbitalSpace> F_i_a_[NSpinCases1];
150 Ref<OrbitalSpace> F_i_p_[NSpinCases1];
151 Ref<OrbitalSpace> F_i_P_[NSpinCases1];
152 Ref<OrbitalSpace> F_a_a_[NSpinCases1];
153 Ref<OrbitalSpace> F_a_A_[NSpinCases1];
154 Ref<OrbitalSpace> h_P_P_[NSpinCases1];
155 Ref<OrbitalSpace> J_i_p_[NSpinCases1];
156 Ref<OrbitalSpace> J_i_P_[NSpinCases1];
157 Ref<OrbitalSpace> J_P_P_[NSpinCases1];
158 Ref<OrbitalSpace> F_A_A_[NSpinCases1];
159 Ref<OrbitalSpace> F_p_P_[NSpinCases1];
160 Ref<OrbitalSpace> gamma_p_p_[NSpinCases1];
161 Ref<OrbitalSpace> gamma_m_m_[NSpinCases1];
162 Ref<OrbitalSpace> gammaFgamma_p_p_[NSpinCases1];
163 Ref<OrbitalSpace> Fgamma_p_P_[NSpinCases1];
164
166 Ref<OrbitalSpace> cabs_space_fockcanonical(SpinCase1 s,
167 double scale_H,
168 double scale_J,
169 double scale_K);
170
172 void init_intermeds_();
174 void init_intermeds_r12_();
176 void init_intermeds_g12_();
178 void r2_contrib_to_X_new_();
180 RefSCMatrix compute_I_(const Ref<OrbitalSpace>& space1,
181 const Ref<OrbitalSpace>& space2,
182 const Ref<OrbitalSpace>& space3,
183 const Ref<OrbitalSpace>& space4);
185 RefSCMatrix compute_r2_(const Ref<OrbitalSpace>& space1,
186 const Ref<OrbitalSpace>& space2,
187 const Ref<OrbitalSpace>& space3,
188 const Ref<OrbitalSpace>& space4);
189
190#if 0
194 RefSCMatrix Delta_DKH_(const Ref<OrbitalSpace>& bra_space,
195 const Ref<OrbitalSpace>& ket_space,
196 SpinCase1 S = Alpha);
197 // Computes T, V and the mass-velocity term in the momentum basis.
198 // It's a modified version of Wavefunction::core_hamiltonian_dk
199 RefSymmSCMatrix hcore_plus_massvelocity_(const Ref<GaussianBasisSet> &bas,
200 const Ref<GaussianBasisSet> &p_bas,
201 bool include_T = true,
202 bool include_V = true,
203 bool include_MV = true);
204#endif
205
206 // Computes the skalar Pauli-Hamiltonian (T + V + mass_velocity + Darwin),
207 // with the mass-velocity term evaluated in the momentum basis.
208 RefSymmSCMatrix pauli(const Ref<GaussianBasisSet> &bas,
209 const Ref<GaussianBasisSet> &pbas = 0,
210 const bool momentum=false);
211 RefSymmSCMatrix pauli_realspace(const Ref<GaussianBasisSet> &bas);
212 RefSymmSCMatrix pauli_momentumspace(const Ref<GaussianBasisSet> &bas,
213 const Ref<GaussianBasisSet> &p_bas);
214
216 RefSCMatrix coulomb_(const SpinCase1 &spin, const Ref<OrbitalSpace>& bra_space,
217 const Ref<OrbitalSpace>& ket_space);
218 RefSCMatrix coulomb_(const Ref<OrbitalSpace>& occ_space, const Ref<OrbitalSpace>& bra_space,
219 const Ref<OrbitalSpace>& ket_space);
221 RefSCMatrix exchange_(const SpinCase1 &spin, const Ref<OrbitalSpace>& bra_space,
222 const Ref<OrbitalSpace>& ket_space);
223 RefSCMatrix exchange_(const Ref<OrbitalSpace>& occ_space, const Ref<OrbitalSpace>& bra_space,
224 const Ref<OrbitalSpace>& ket_space);
225
227 void checkpoint_() const;
228
230 void contrib_to_VXB_a_();
231 void contrib_to_VXB_c_ansatz1_();
232 void contrib_to_VX_GenRefansatz2_();
233 void contrib_to_VX_GenRefansatz2_spinfree_();
235 void contrib_to_VXB_a_vbsneqobs_();
237 void contrib_to_B_DKH_a_();
238
240 void compute_mp2_pair_energies_(RefSCVector& emp2pair,
241 SpinCase2 S,
242 const Ref<OrbitalSpace>& space1,
243 const Ref<OrbitalSpace>& space2,
244 const Ref<OrbitalSpace>& space3,
245 const Ref<OrbitalSpace>& space4,
246 const std::string& tform_key);
247
248 void compute_A(); //< computes A ("coupling") matrix for all spincases
249
255 void compute_A_direct_(RefSCMatrix& A,
256 const Ref<OrbitalSpace>& space1,
257 const Ref<OrbitalSpace>& space2,
258 const Ref<OrbitalSpace>& space3,
259 const Ref<OrbitalSpace>& space4,
260 const Ref<OrbitalSpace>& rispace2,
261 const Ref<OrbitalSpace>& rispace4,
262 bool antisymmetrize);
263
264 // make compute_tbint_tensor_ public for now
265 public:
290 template <typename DataProcess, bool CorrFactorInBra, bool CorrFactorInKet>
292 TwoBodyOper::type tbint_type,
293 const Ref<OrbitalSpace>& space1,
294 const Ref<OrbitalSpace>& space2,
295 const Ref<OrbitalSpace>& space3,
296 const Ref<OrbitalSpace>& space4,
297 bool antisymmetrize,
298 const std::vector<std::string>& tform_keys);
299
300 private:
340 template <typename DataProcessBra,
341 typename DataProcessKet,
342 typename DataProcessBraKet,
343 bool CorrFactorInBra,
344 bool CorrFactorInKet,
345 bool CorrFactorInInt>
346 DEPRECATED void contract_tbint_tensor(
347 RefSCMatrix& T,
348 TwoBodyOper::type tbint_type_bra,
349 TwoBodyOper::type tbint_type_ket,
350 const Ref<OrbitalSpace>& space1_bra,
351 const Ref<OrbitalSpace>& space2_bra,
352 const Ref<OrbitalSpace>& space1_intb,
353 const Ref<OrbitalSpace>& space2_intb,
354 const Ref<OrbitalSpace>& space1_ket,
355 const Ref<OrbitalSpace>& space2_ket,
356 const Ref<OrbitalSpace>& space1_intk,
357 const Ref<OrbitalSpace>& space2_intk,
358 const Ref<mbptr12::TwoParticleContraction>& tpcontract,
359 bool antisymmetrize,
360 const std::vector<std::string>& tformkeys_bra,
361 const std::vector<std::string>& tformkeys_ket);
362
370 template <bool CorrFactorInBra,
371 bool CorrFactorInKet>
372 void DEPRECATED contract_tbint_tensor(
373 RefSCMatrix& T,
374 TwoBodyOper::type tbint_type_bra,
375 TwoBodyOper::type tbint_type_ket,
376 double scale,
377 const Ref<OrbitalSpace>& space1_bra,
378 const Ref<OrbitalSpace>& space2_bra,
379 const Ref<OrbitalSpace>& space1_intb,
380 const Ref<OrbitalSpace>& space2_intb,
381 const Ref<OrbitalSpace>& space1_ket,
382 const Ref<OrbitalSpace>& space2_ket,
383 const Ref<OrbitalSpace>& space1_intk,
384 const Ref<OrbitalSpace>& space2_intk,
385 bool antisymmetrize,
386 const std::vector<std::string>& tformkeys_bra,
387 const std::vector<std::string>& tformkeys_ket);
388
394 template <bool CorrFactorInBra,
395 bool CorrFactorInKet>
396 void
397 contract_tbint_tensor(
398 std::vector< Ref<DistArray4> >& results,
399 TwoBodyOper::type tbint_type_bra,
400 TwoBodyOper::type tbint_type_ket,
401 double scale,
402 const Ref<OrbitalSpace>& space1_bra,
403 const Ref<OrbitalSpace>& space2_bra,
404 const Ref<OrbitalSpace>& space1_intb,
405 const Ref<OrbitalSpace>& space2_intb,
406 const Ref<OrbitalSpace>& space1_ket,
407 const Ref<OrbitalSpace>& space2_ket,
408 const Ref<OrbitalSpace>& space1_intk,
409 const Ref<OrbitalSpace>& space2_intk,
410 bool antisymmetrize,
411 const std::vector<std::string>& tformkeys_bra,
412 const std::vector<std::string>& tformkeys_ket);
413
414 // for now make it public
415 public:
451 template <typename DataProcessBra,
452 typename DataProcessKet,
453 bool CorrFactorInBra,
454 bool CorrFactorInKet,
455 bool CorrFactorInInt>
457 SpinCase2 pairspin,
458 TwoBodyTensorInfo tbtensor_type_bra,
459 TwoBodyTensorInfo tbtensor_type_ket,
460 const Ref<OrbitalSpace>& space1_bra,
461 const Ref<OrbitalSpace>& space1_intb,
462 const Ref<OrbitalSpace>& space2_intb,
463 const Ref<OrbitalSpace>& space3_intb,
464 const Ref<OrbitalSpace>& space1_ket,
465 const Ref<OrbitalSpace>& space1_intk,
466 const Ref<OrbitalSpace>& space2_intk,
467 const Ref<OrbitalSpace>& space3_intk,
468 const Ref<mbptr12::TwoParticleContraction>& tpcontract,
469 const std::vector<std::string>& tformkeys_bra,
470 const std::vector<std::string>& tformkeys_ket);
471
472 private:
484 void compute_X_(RefSCMatrix& X,
485 SpinCase2 sc2,
486 const Ref<OrbitalSpace>& bra1,
487 const Ref<OrbitalSpace>& bra2,
488 const Ref<OrbitalSpace>& ket1,
489 const Ref<OrbitalSpace>& ket2,
490 bool F2_only = false);
509 void compute_FxF_(RefSCMatrix& FxF,
510 SpinCase2 sc2,
511 const Ref<OrbitalSpace>& bra1,
512 const Ref<OrbitalSpace>& bra2,
513 const Ref<OrbitalSpace>& ket1,
514 const Ref<OrbitalSpace>& ket2,
515 const Ref<OrbitalSpace>& int1,
516 const Ref<OrbitalSpace>& int2,
517 const Ref<OrbitalSpace>& intk1,
518 const Ref<OrbitalSpace>& intk2,
519 const Ref<OrbitalSpace>& intkx1,
520 const Ref<OrbitalSpace>& intkx2
521 );
522
524 void AT2_contrib_to_V_();
525
527 void AF12_contrib_to_B_();
528
530 void compute_B_gbc_();
531
534 void compute_B_fX_();
535
537 void compute_BB_();
538
540 void compute_BC_();
541 void compute_BC_ansatz1_();
542 void compute_BC_GenRefansatz2_();
543 void compute_BC_GenRefansatz2_spinfree();
544
545
547 void compute_BCp_();
548
550 void compute_BApp_();
551
553 void compute_B_DKH_();
554
556 const RefSCVector& compute_emp2(SpinCase2 s);
561 double compute_emp2_obs_singles(bool obs_singles);
562#if 0
566 double compute_emp2_cabs_singles();
567#endif
574 double compute_emp2_cabs_singles_noncanonical(bool vir_cabs_coupling);
578 double compute_emp2_cabs_singles_noncanonical_ccsd(const RefSCMatrix& T1_ia_alpha,
579 const RefSCMatrix& T1_ia_beta);
580
583 void compute_R_vbsneqobs_(const Ref<TwoBodyMOIntsTransform>& ipjq_tform,
584 RefSCMatrix& Raa, RefSCMatrix& Rab);
585
587 void compute_amps_();
588
594 void globally_sum_scmatrix_(RefSCMatrix& A, bool to_all_tasks = false, bool to_average = false);
595
601 void globally_sum_scvector_(RefSCVector& A, bool to_all_tasks = false, bool to_average = false);
602
607 void globally_sum_intermeds_(bool to_all_tasks = false);
608
609public:
613 ~R12IntEval();
614
616 virtual void obsolete();
617
618 void debug(int d) { debug_ = d; }
619
620 int dk() const { return r12world()->refwfn()->dk(); }
621
622#if 1
623 const Ref<R12Technology::CorrelationFactor>& corrfactor() const { return r12world()->r12tech()->corrfactor(); }
624 R12Technology::ABSMethod abs_method() const { return r12world()->r12tech()->abs_method(); }
625 const Ref<R12Technology::R12Ansatz>& ansatz() const { return r12world()->r12tech()->ansatz(); }
626 bool spin_polarized() const { return r12world()->refwfn()->spin_polarized(); }
627 bool gbc() const { return r12world()->r12tech()->gbc(); }
628 bool ebc() const { return r12world()->r12tech()->ebc(); }
629 bool coupling() const { return r12world()->r12tech()->coupling(); }
630 bool compute_1rdm() const { return r12world()->r12tech()->compute_1rdm(); }
631 bool coupling_1rdm_f12b() const { return r12world()->r12tech()->coupling_1rdm_f12b(); }
632 R12Technology::StandardApproximation stdapprox() const { return r12world()->r12tech()->stdapprox(); }
633 bool omit_P() const { return r12world()->r12tech()->omit_P(); }
634 const Ref<MOIntsTransformFactory>& tfactory() const { return r12world()->world()->tfactory(); };
635 const Ref<MOIntsRuntime>& moints_runtime() const { return r12world()->world()->moints_runtime(); };
636 const Ref<TwoBodyFourCenterMOIntsRuntime>& moints_runtime4() const { return r12world()->world()->moints_runtime()->runtime_4c(); };
637 const Ref<FockBuildRuntime>& fockbuild_runtime() const { return r12world()->world()->fockbuild_runtime(); };
638
639#endif
640
642 const Ref<R12WavefunctionWorld>& r12world() const { return r12world_; }
643
645 RefSCDimension dim_oo(SpinCase2 S) const { return dim_oo_[S]; }
647 RefSCDimension dim_vv(SpinCase2 S) const { return dim_vv_[S]; }
649 RefSCDimension dim_aa(SpinCase2 S) const { return dim_aa_[S]; }
651 RefSCDimension dim_f12(SpinCase2 S) const { return dim_f12_[S]; }
653 RefSCDimension dim_GG(SpinCase2 S) const { return dim_GG_[S]; }
655 RefSCDimension dim_gg(SpinCase2 S) const { return dim_gg_[S]; }
656
658 int nspincases1() const {
659 if (r12world()->spinadapted()) return 1;
660 return ::sc::nspincases1(r12world()->refwfn()->spin_polarized());
661 }
663 int nspincases2() const {
664 if (r12world()->spinadapted()) return 1;
665 return ::sc::nspincases2(r12world()->refwfn()->spin_polarized());
666 }
667
669 virtual void compute();
670
672 bool bc() const;
673
676
678 const RefSCMatrix& V(SpinCase2 S);
680 const RefSCMatrix& V();
682 RefSymmSCMatrix X(SpinCase2 S);
686 RefSymmSCMatrix B(SpinCase2 S);
690 RefSymmSCMatrix BB(SpinCase2 S);
692 const RefSCMatrix& A(SpinCase2 S);
694 const RefSCMatrix& T2(SpinCase2 S);
696 const RefSCMatrix& F12(SpinCase2 S);
697
699 RefSCMatrix V(SpinCase2 spincase2,
700 const Ref<OrbitalSpace>& p,
701 const Ref<OrbitalSpace>& q);
703 std::vector<Ref<DistArray4> > V_distarray4(SpinCase2 spincase2,
704 const Ref<OrbitalSpace>& p,
705 const Ref<OrbitalSpace>& q);
707 std::vector<Ref<DistArray4> > U_distarray4(SpinCase2 spincase2,
708 const Ref<OrbitalSpace>& p,
709 const Ref<OrbitalSpace>& q);
711 RefSymmSCMatrix P(SpinCase2 S);
713 double C_CuspConsistent(int i,int j,int k,int l,SpinCase2 pairspin);
714
722 double emp2_cabs_singles(bool vir_cabs_coupling = true);
724 double emp2_cabs_singles(const RefSCMatrix& T1_ia_alpha,
725 const RefSCMatrix& T1_ia_beta);
727 const RefSCVector& emp2(SpinCase2 S);
728
735 const RefSCMatrix& T1_cabs(SpinCase1 spin) const;
736
738 const Ref<OrbitalSpace>& occ_act(SpinCase1 S) const;
740 const Ref<OrbitalSpace>& occ(SpinCase1 S) const;
742 const Ref<OrbitalSpace>& vir_act(SpinCase1 S) const;
744 const Ref<OrbitalSpace>& vir(SpinCase1 S) const;
746 const Ref<OrbitalSpace>& orbs(SpinCase1 S) const;
748 const Ref<OrbitalSpace>& GGspace(SpinCase1 S) const;
750 const Ref<OrbitalSpace>& ggspace(SpinCase1 S) const;
756 RefSymmSCMatrix ordm(SpinCase1 S) const;
761 RefSymmSCMatrix ordm_occ_av() const;
763 const Ref<OrbitalSpace>& hj_x_P(SpinCase1 S);
765 const Ref<OrbitalSpace>& hj_x_A(SpinCase1 S);
767 const Ref<OrbitalSpace>& hj_x_p(SpinCase1 S);
769 const Ref<OrbitalSpace>& hj_x_m(SpinCase1 S);
771 const Ref<OrbitalSpace>& hj_x_a(SpinCase1 S);
773 const Ref<OrbitalSpace>& hj_i_P(SpinCase1 S);
775 const Ref<OrbitalSpace>& hj_i_A(SpinCase1 S);
777 const Ref<OrbitalSpace>& hj_i_p(SpinCase1 S);
779 const Ref<OrbitalSpace>& hj_i_m(SpinCase1 S);
781 const Ref<OrbitalSpace>& hj_i_a(SpinCase1 S);
783 const Ref<OrbitalSpace>& hj_m_m(SpinCase1 S);
785 const Ref<OrbitalSpace>& hj_m_p(SpinCase1 S);
787 const Ref<OrbitalSpace>& hj_a_A(SpinCase1 S);
789 const Ref<OrbitalSpace>& hj_p_P(SpinCase1 S);
791 const Ref<OrbitalSpace>& hj_p_A(SpinCase1 S);
793 const Ref<OrbitalSpace>& hj_p_p(SpinCase1 S);
795 const Ref<OrbitalSpace>& hj_p_m(SpinCase1 S);
797 const Ref<OrbitalSpace>& hj_p_a(SpinCase1 S);
799 const Ref<OrbitalSpace>& hj_P_P(SpinCase1 S);
801 const Ref<OrbitalSpace>& K_x_P(SpinCase1 S);
803 const Ref<OrbitalSpace>& K_x_A(SpinCase1 S);
805 const Ref<OrbitalSpace>& K_x_p(SpinCase1 S);
807 const Ref<OrbitalSpace>& K_x_m(SpinCase1 S);
809 const Ref<OrbitalSpace>& K_x_a(SpinCase1 S);
811 const Ref<OrbitalSpace>& K_i_P(SpinCase1 S);
813 const Ref<OrbitalSpace>& K_i_A(SpinCase1 S);
815 const Ref<OrbitalSpace>& K_i_p(SpinCase1 S);
817 const Ref<OrbitalSpace>& K_i_m(SpinCase1 S);
819 const Ref<OrbitalSpace>& K_i_a(SpinCase1 S);
821 const Ref<OrbitalSpace>& K_m_a(SpinCase1 S);
823 const Ref<OrbitalSpace>& K_a_a(SpinCase1 S);
825 const Ref<OrbitalSpace>& K_a_p(SpinCase1 S);
827 const Ref<OrbitalSpace>& K_a_P(SpinCase1 S);
829 const Ref<OrbitalSpace>& K_p_P(SpinCase1 S);
831 const Ref<OrbitalSpace>& K_p_A(SpinCase1 S);
833 const Ref<OrbitalSpace>& K_p_p(SpinCase1 S);
835 const Ref<OrbitalSpace>& K_p_m(SpinCase1 S);
837 const Ref<OrbitalSpace>& K_p_a(SpinCase1 S);
839 const Ref<OrbitalSpace>& K_A_P(SpinCase1 S);
841 const Ref<OrbitalSpace>& K_P_P(SpinCase1 S);
843 const Ref<OrbitalSpace>& F_x_P(SpinCase1 S);
845 const Ref<OrbitalSpace>& F_x_A(SpinCase1 S);
847 const Ref<OrbitalSpace>& F_x_p(SpinCase1 S);
849 const Ref<OrbitalSpace>& F_x_m(SpinCase1 S);
851 const Ref<OrbitalSpace>& F_x_a(SpinCase1 S);
853 const Ref<OrbitalSpace>& F_i_P(SpinCase1 S);
855 const Ref<OrbitalSpace>& F_i_A(SpinCase1 S);
857 const Ref<OrbitalSpace>& F_i_p(SpinCase1 S);
859 const Ref<OrbitalSpace>& F_i_m(SpinCase1 S);
861 const Ref<OrbitalSpace>& F_i_a(SpinCase1 S);
863 const Ref<OrbitalSpace>& F_m_m(SpinCase1 S);
865 const Ref<OrbitalSpace>& F_m_a(SpinCase1 S);
867 const Ref<OrbitalSpace>& F_m_p(SpinCase1 S);
869 const Ref<OrbitalSpace>& F_m_P(SpinCase1 S);
871 const Ref<OrbitalSpace>& F_gg_P(SpinCase1 S);
873 const Ref<OrbitalSpace>& F_m_A(SpinCase1 S);
875 const Ref<OrbitalSpace>& F_a_a(SpinCase1 S);
877 const Ref<OrbitalSpace>& F_a_A(SpinCase1 S);
879 const Ref<OrbitalSpace>& F_A_A(SpinCase1 S);
881 const Ref<OrbitalSpace>& F_p_P(SpinCase1 S);
883 const Ref<OrbitalSpace>& F_p_A(SpinCase1 S);
885 const Ref<OrbitalSpace>& F_p_p(SpinCase1 S);
887 const Ref<OrbitalSpace>& F_p_m(SpinCase1 S);
889 const Ref<OrbitalSpace>& F_p_a(SpinCase1 S);
891 const Ref<OrbitalSpace>& F_P_P(SpinCase1 S);
893 const Ref<OrbitalSpace>& h_P_P(SpinCase1 S);
895 const Ref<OrbitalSpace>& J_x_p(SpinCase1 S);
897 const Ref<OrbitalSpace>& J_i_p(SpinCase1 S);
899 const Ref<OrbitalSpace>& J_i_P(SpinCase1 S);
901 const Ref<OrbitalSpace>& J_P_P(SpinCase1 S);
902
904 const Ref<OrbitalSpace>& gamma_p_p(SpinCase1 S);
907 const Ref<OrbitalSpace>& gamma_m_m_av();
912 const Ref<OrbitalSpace>& Fgamma_p_P(SpinCase1 S);
915// const Ref<OrbitalSpace>& Fgamma_P_p(SpinCase1 S,const RefSymmSCMatrix &gamma,Ref<OrbitalSpace>& FGamma);
916// const Ref<OrbitalSpace>& gammaF_p_P(SpinCase1 S,const RefSymmSCMatrix &gamma,Ref<OrbitalSpace>& GammaF);
919 Ref<OrbitalSpace> obtensor_p_A(const RefSCMatrix &obtensor,SpinCase1 S);
920
922 std::string transform_label(const Ref<OrbitalSpace>& space1,
923 const Ref<OrbitalSpace>& space2,
924 const Ref<OrbitalSpace>& space3,
925 const Ref<OrbitalSpace>& space4,
926 const std::string& operator_label = std::string()) const;
928 std::string transform_label(const Ref<OrbitalSpace>& space1,
929 const Ref<OrbitalSpace>& space2,
930 const Ref<OrbitalSpace>& space3,
931 const Ref<OrbitalSpace>& space4,
932 unsigned int f12,
933 const std::string& operator_label = std::string()) const;
935 std::string transform_label(const Ref<OrbitalSpace>& space1,
936 const Ref<OrbitalSpace>& space2,
937 const Ref<OrbitalSpace>& space3,
938 const Ref<OrbitalSpace>& space4,
939 unsigned int f12_left,
940 unsigned int f12_right,
941 const std::string& operator_label = std::string()) const;
942
948 const Ref<OrbitalSpace>& space1,
949 const Ref<OrbitalSpace>& space2,
950 const Ref<OrbitalSpace>& space3,
951 const Ref<OrbitalSpace>& space4,
952 bool antisymmetrize,
953 const std::string& tform_key);
960 const Ref<OrbitalSpace>& space1,
961 const Ref<OrbitalSpace>& space2,
962 const Ref<OrbitalSpace>& space3,
963 const Ref<OrbitalSpace>& space4,
964 bool antisymmetrize,
965 const std::vector<std::string>& transform_keys);
966
972 const Ref<OrbitalSpace>& ket_space, SpinCase1 S = Alpha,
973 double scale_J = 1.0, double scale_K = 1.0, double scale_H = 1.0,
974 // leave at default to use R12Technology's pauli flag, else set to 0 or 1
975 int override_pauli = -1);
976
977
979 RefSCMatrix V_cabs(SpinCase2 spincase2,
980 const Ref<OrbitalSpace>& p,
981 const Ref<OrbitalSpace>& q);
982
985 const Ref<OrbitalSpace>& q);
986
987#if defined(MPQC_NEW_FEATURES)
991 void V_diag_ta();
992
993 // \param orbital index of the orbital for self-energy calculation, -1 = HOMO (default), +1 = LUMO, 0 = do nothing
994 void gf2_r12(int orbital = -1);
995
996 // TiledArray-based MP2-F12 one-electron properties
997 void compute_TA_mp2f12_1rdm();
998#endif
999
1000 void compute_ccr12_1rdm(const RefSCMatrix& T1, const Ref<DistArray4> (&T2)[NSpinCases2]);
1001
1004 return this->r12world()->world()->tfactory()->orbital_registry();
1005 }
1008 return this->r12world()->world()->tfactory()->ao_registry();
1009 }
1010
1011};
1012
1013std::vector< Ref<DistArray4> >
1014A_distarray4(SpinCase2 spincase2, const Ref<R12IntEval>& r12eval);
1015
1016// compute orbital Z-vector from F12 contribution
1017RefSCMatrix Onerdm_X_F12(SpinCase1 spin, const Ref<R12IntEval>& r12eval, int debug);
1018
1019// compute orbital relaxation contributions from CABS Singles
1020RefSCMatrix Onerdm_X_CABS_Singles(SpinCase1 spin,
1021 const Ref<R12IntEval>& r12eval,
1022 const Ref<R12EnergyIntermediates>& r12intermediates,
1023 int debug);
1024
1025} // end of namespace sc
1026
1027#endif
1028
1029// Local Variables:
1030// mode: c++
1031// c-file-style: "CLJ"
1032// End:
1033
1034
R12IntEval is the top-level class which computes intermediates occuring in R12 theories.
Definition r12int_eval.h:50
const Ref< OrbitalSpace > & F_m_a(SpinCase1 S)
Form <a|F|m> space.
RefSCMatrix V_cabs(SpinCase2 spincase2, const Ref< OrbitalSpace > &p, const Ref< OrbitalSpace > &q)
Compute V intermediates using CABS/CABS+ approach.
const Ref< OrbitalSpace > & J_i_P(SpinCase1 S)
Form <P|J|i> space.
RefSymmSCMatrix X()
Returns spin-free intermediate X.
const Ref< OrbitalSpace > & F_A_A(SpinCase1 S)
Form <A|F|A> space.
const Ref< OrbitalSpace > & hj_i_m(SpinCase1 S)
Form <m|h+J|i> space.
const RefSCMatrix & A(SpinCase2 S)
Returns S block of intermediate A.
const Ref< OrbitalSpace > & hj_x_m(SpinCase1 S)
Form <m|h+J|x> space.
RefSCMatrix V(SpinCase2 spincase2, const Ref< OrbitalSpace > &p, const Ref< OrbitalSpace > &q)
Compute .
const Ref< OrbitalSpace > & ggspace(SpinCase1 S) const
Returns the space for spin case S from which geminal-generating substitutions are allowed.
const Ref< OrbitalSpace > & K_x_P(SpinCase1 S)
Form <P|K|x> space.
const Ref< OrbitalSpace > & F_x_m(SpinCase1 S)
Form <m|F|x> space.
int nspincases1() const
Returns the number of unique spin cases.
Definition r12int_eval.h:658
RefSymmSCMatrix B()
Returns spin-free intermediate B.
const Ref< OrbitalSpace > & Fgamma_p_P()
Form <P|Fgamma|p> space using spin-average rdm1.
const Ref< OrbitalSpace > & F_gg_P(SpinCase1 S)
Form <P|F|gg> space.
const Ref< OrbitalSpace > & F_m_m(SpinCase1 S)
Form <m|F|m> space.
const Ref< OrbitalSpace > & K_p_a(SpinCase1 S)
Form <a|K|p> space.
R12IntEval(const Ref< R12WavefunctionWorld > &r12w)
Constructs R12IntEval.
const Ref< OrbitalSpace > & F_i_A(SpinCase1 S)
Form <A|F|i> space.
const Ref< OrbitalSpace > & vir(SpinCase1 S) const
Returns the vir space for spin case S.
const Ref< OrbitalSpace > & F_i_p(SpinCase1 S)
Form <p|F|i> space.
const RefSCMatrix & V(SpinCase2 S)
Returns S block of intermediate V.
const Ref< OrbitalSpace > & hj_p_p(SpinCase1 S)
Form <p|h+J|p> space.
const Ref< OrbitalSpace > & J_x_p(SpinCase1 S)
Form <p|J|x> space.
RefSCDimension dim_aa(SpinCase2 S) const
Dimension for any/any pairs of spin case S.
Definition r12int_eval.h:649
const Ref< OrbitalSpace > & F_p_P(SpinCase1 S)
Form <P|F|p> space.
const Ref< OrbitalSpace > & K_x_m(SpinCase1 S)
Form <i|K|x> space.
const Ref< OrbitalSpace > & K_a_p(SpinCase1 S)
Form <p|K|a> space.
const Ref< OrbitalSpace > & occ_act(SpinCase1 S) const
Returns the act occ space for spin case S.
const Ref< OrbitalSpace > & hj_x_P(SpinCase1 S)
Form <P|h+J|x> space.
DEPRECATED void compute_tbint_tensor(RefSCMatrix &T, TwoBodyOper::type tbint_type, const Ref< OrbitalSpace > &space1, const Ref< OrbitalSpace > &space2, const Ref< OrbitalSpace > &space3, const Ref< OrbitalSpace > &space4, bool antisymmetrize, const std::vector< std::string > &tform_keys)
compute_tbint_tensor computes a 2-body tensor T using integrals of type tbint_type.
void compute_F12_(RefSCMatrix &F12, const Ref< OrbitalSpace > &space1, const Ref< OrbitalSpace > &space2, const Ref< OrbitalSpace > &space3, const Ref< OrbitalSpace > &space4, bool antisymmetrize, const std::vector< std::string > &transform_keys)
Compute F12 integrals in basis <space1, space3 | f12 | space2, space4>.
const RefSCMatrix & T1_cabs(SpinCase1 spin) const
returns CABS singles amplitudes.
const Ref< OrbitalSpace > & gamma_p_p(SpinCase1 S)
Form <p|gamma|p> space.
const Ref< OrbitalSpaceRegistry > & orbital_registry() const
returns the OrbitalSpaceRegistry object
Definition r12int_eval.h:1003
const Ref< OrbitalSpace > & K_a_a(SpinCase1 S)
Form <a|K|a> space.
std::string transform_label(const Ref< OrbitalSpace > &space1, const Ref< OrbitalSpace > &space2, const Ref< OrbitalSpace > &space3, const Ref< OrbitalSpace > &space4, unsigned int f12_left, unsigned int f12_right, const std::string &operator_label=std::string()) const
version of transform_label() applicable when left and right correlation factors differ
const Ref< OrbitalSpace > & gamma_p_p_av()
Form spin-average <p|gamma|p> space using average (instead of total) rdm.
double emp2_cabs_singles(const RefSCMatrix &T1_ia_alpha, const RefSCMatrix &T1_ia_beta)
Returns the CABS singles MP2 energy, with fixed conventional T1 amplitudes.
RefSCDimension dim_f12(SpinCase2 S) const
Dimension for geminal functions of spin case S = # of correlation factors x dim_GG.
Definition r12int_eval.h:651
const RefSCMatrix & T2(SpinCase2 S)
Returns S block of intermediate T2.
const Ref< OrbitalSpace > & cabs_space_hcanonical(SpinCase1 s)
compute CABS space for spin s canonicalized by diagonalization of core hamitonian
std::vector< Ref< DistArray4 > > U_distarray4(SpinCase2 spincase2, const Ref< OrbitalSpace > &p, const Ref< OrbitalSpace > &q)
Compute .
const Ref< OrbitalSpace > & K_p_P(SpinCase1 S)
Form <P|K|p> space.
const Ref< OrbitalSpace > & hj_p_A(SpinCase1 S)
Form <A|h+J|p> space.
const Ref< OrbitalSpace > & F_p_p(SpinCase1 S)
Form <p|F|p> space.
const Ref< OrbitalSpace > & J_P_P(SpinCase1 S)
Form <P|J|P> space.
const Ref< OrbitalSpace > & vir_act(SpinCase1 S) const
Returns the act vir space for spin case S.
Ref< OrbitalSpace > obtensor_p_A(const RefSCMatrix &obtensor, SpinCase1 S)
Form <A|obtensor|p> space obtensor should have the dimension ncabs*norbs.
RefSymmSCMatrix BB(SpinCase2 S)
Returns S block of the difference between intermediate B of approximations B and A'.
RefSCDimension dim_gg(SpinCase2 S) const
Dimension of orbital product space from which geminal substitutions are allowed.
Definition r12int_eval.h:655
bool bc() const
does Brillouin condition hold?
const Ref< OrbitalSpace > & cabs_space_canonical(SpinCase1 s)
compute canonical CABS space for spin s
const Ref< OrbitalSpace > & gammaFgamma_p_p(SpinCase1 S)
Form <p|gammaFgamma|p> space.
const Ref< OrbitalSpace > & hj_i_P(SpinCase1 S)
Form <P|h+J|i> space.
const Ref< OrbitalSpace > & F_x_a(SpinCase1 S)
Form <a|F|x> space.
RefSCDimension dim_vv(SpinCase2 S) const
Dimension for active-vir/active-vir pairs of spin case S.
Definition r12int_eval.h:647
const Ref< OrbitalSpace > & hj_i_a(SpinCase1 S)
Form <a|h+J|i> space.
const Ref< OrbitalSpace > & hj_x_a(SpinCase1 S)
Form <a|h+J|x> space.
const Ref< OrbitalSpace > & hj_m_m(SpinCase1 S)
Form <m|h+J|m> space.
const Ref< OrbitalSpace > & F_m_p(SpinCase1 S)
Form <p|F|m> space.
const Ref< OrbitalSpace > & F_p_A(SpinCase1 S)
Form <A|F|p> space.
Ref< R12Amplitudes > amps()
Returns amplitudes of pair correlation functions.
int nspincases2() const
Returns the number of unique combinations of 2 spin cases.
Definition r12int_eval.h:663
const Ref< OrbitalSpace > & K_a_P(SpinCase1 S)
Form <P|K|a> space.
const Ref< OrbitalSpace > & F_m_A(SpinCase1 S)
Form <A|F|m> space.
const RefSCMatrix & F12(SpinCase2 S)
Returns S block of intermediate F12.
RefSymmSCMatrix P(SpinCase2 S)
Compute P = RgR.
const Ref< OrbitalSpace > & occ(SpinCase1 S) const
Returns the occ space for spin case S.
const Ref< OrbitalSpace > & F_x_A(SpinCase1 S)
Form <A|F|x> space.
const Ref< OrbitalSpace > & K_i_P(SpinCase1 S)
Form <P|K|i> space.
RefSymmSCMatrix B(SpinCase2 S)
Returns S block of intermediate B.
double C_CuspConsistent(int i, int j, int k, int l, SpinCase2 pairspin)
Returns the cusp consistent coefficient .
RefSymmSCMatrix ordm_av() const
Returns the average 1-RDM in the `‘MO’' basis (i.e. that provided by orbs() ): (alpha-rdm1 + beta-rdm...
const Ref< OrbitalSpace > & F_i_m(SpinCase1 S)
Form <m|F|i> space.
const Ref< OrbitalSpace > & K_x_a(SpinCase1 S)
Form <a|K|x> space.
const Ref< OrbitalSpace > & F_m_P(SpinCase1 S)
Form <P|F|m> space.
const Ref< R12WavefunctionWorld > & r12world() const
the R12World in which this object lives
Definition r12int_eval.h:642
const Ref< OrbitalSpace > & Fgamma_p_P(SpinCase1 S)
Form <P|Fgamma|p> space.
const Ref< OrbitalSpace > & K_m_a(SpinCase1 S)
Form <a|K|m> space.
const Ref< OrbitalSpace > & J_i_p(SpinCase1 S)
Form <p|J|i> space.
const RefSCVector & emp2(SpinCase2 S)
Returns alpha-alpha MP2 pair energies.
const Ref< OrbitalSpace > & hj_P_P(SpinCase1 S)
Form <P|h+J|P> space.
const Ref< OrbitalSpace > & F_x_p(SpinCase1 S)
Form <p|F|x> space.
const Ref< OrbitalSpace > & K_x_p(SpinCase1 S)
Form <p|K|x> space.
const RefSCMatrix & V()
Returns spin-free intermediate V.
RefSCDimension dim_oo(SpinCase2 S) const
Dimension for active-occ/active-occ pairs of spin case S.
Definition r12int_eval.h:645
RefSCMatrix fock(const Ref< OrbitalSpace > &bra_space, const Ref< OrbitalSpace > &ket_space, SpinCase1 S=Alpha, double scale_J=1.0, double scale_K=1.0, double scale_H=1.0, int override_pauli=-1)
Compute the Fock matrix between bra_ and ket_ spaces of spin S.
const Ref< OrbitalSpace > & hj_p_a(SpinCase1 S)
Form <a|h+J|p> space.
const Ref< OrbitalSpace > & hj_i_p(SpinCase1 S)
Form <p|h+J|i> space.
const Ref< OrbitalSpace > & h_P_P(SpinCase1 S)
Form <P|h|P> space.
std::string transform_label(const Ref< OrbitalSpace > &space1, const Ref< OrbitalSpace > &space2, const Ref< OrbitalSpace > &space3, const Ref< OrbitalSpace > &space4, unsigned int f12, const std::string &operator_label=std::string()) const
Generates canonical id for transform. f12 is the index of the correlation function.
const Ref< OrbitalSpace > & GGspace(SpinCase1 S) const
Returns the geminal-generating orbital space for spin case S.
std::vector< Ref< DistArray4 > > V_distarray4(SpinCase2 spincase2, const Ref< OrbitalSpace > &p, const Ref< OrbitalSpace > &q)
Compute .
RefSCMatrix V_genref_spinfree(const Ref< OrbitalSpace > &p, const Ref< OrbitalSpace > &q)
Compute V intermediates in SF-[2]R12.
RefSCDimension dim_GG(SpinCase2 S) const
Dimension of orbital product space that multiply the correlation factor to produce geminal functions.
Definition r12int_eval.h:653
void save_data_state(StateOut &)
Save the base classes (with save_data_state) and the members in the same order that the StateIn CTOR ...
const Ref< OrbitalSpace > & F_p_a(SpinCase1 S)
Form <a|F|p> space.
double emp2_obs_singles()
Returns the OBS singles correction to the MP2 energy.
std::string transform_label(const Ref< OrbitalSpace > &space1, const Ref< OrbitalSpace > &space2, const Ref< OrbitalSpace > &space3, const Ref< OrbitalSpace > &space4, const std::string &operator_label=std::string()) const
Generates canonical id for transform. no correlation function included.
const Ref< OrbitalSpace > & F_a_a(SpinCase1 S)
Form <a|F|a> space.
RefSymmSCMatrix X(SpinCase2 S)
Returns S block of intermediate X.
const Ref< OrbitalSpace > & hj_x_A(SpinCase1 S)
Form <A|h+J|x> space.
void compute_T2_(RefSCMatrix &T2, const Ref< OrbitalSpace > &space1, const Ref< OrbitalSpace > &space2, const Ref< OrbitalSpace > &space3, const Ref< OrbitalSpace > &space4, bool antisymmetrize, const std::string &tform_key)
Compute T2 amplitude in basis <space1, space3 | space2, space4>.
double emp2_cabs_singles(bool vir_cabs_coupling=true)
Returns the CABS singles correction to the MP2 energy.
const Ref< OrbitalSpace > & hj_x_p(SpinCase1 S)
Form <p|h+J|x> space.
const Ref< OrbitalSpace > & K_p_A(SpinCase1 S)
Form <A|K|p> space.
RefSymmSCMatrix ordm() const
Returns the total 1-RDM in the `‘MO’' basis (i.e. that provided by orbs() )
const Ref< OrbitalSpace > & hj_i_A(SpinCase1 S)
Form <A|h+J|i> space.
const Ref< OrbitalSpace > & K_P_P(SpinCase1 S)
Form <P|K|P> space.
const Ref< OrbitalSpace > & F_i_P(SpinCase1 S)
Form <P|F|i> space.
const Ref< OrbitalSpace > & orbs(SpinCase1 S) const
Returns the OBS space for spin case S.
const Ref< OrbitalSpace > & K_x_A(SpinCase1 S)
Form <A|K|x> space.
const Ref< OrbitalSpace > & hj_p_m(SpinCase1 S)
Form <m|h+J|p> space.
void contract_tbint_tensors_to_obtensor(RefSCMatrix &T, SpinCase2 pairspin, TwoBodyTensorInfo tbtensor_type_bra, TwoBodyTensorInfo tbtensor_type_ket, const Ref< OrbitalSpace > &space1_bra, const Ref< OrbitalSpace > &space1_intb, const Ref< OrbitalSpace > &space2_intb, const Ref< OrbitalSpace > &space3_intb, const Ref< OrbitalSpace > &space1_ket, const Ref< OrbitalSpace > &space1_intk, const Ref< OrbitalSpace > &space2_intk, const Ref< OrbitalSpace > &space3_intk, const Ref< mbptr12::TwoParticleContraction > &tpcontract, const std::vector< std::string > &tformkeys_bra, const std::vector< std::string > &tformkeys_ket)
<space1bra space1_intb |Tbra| space2_intb space3_intb> * <space2_intk space3_intk |Tket| space1ket sp...
Definition contract_tbint_tensors_to_obtensor.h:46
const Ref< OrbitalSpace > & K_i_A(SpinCase1 S)
Form <A|K|i> space.
const Ref< OrbitalSpace > & K_i_p(SpinCase1 S)
Form <p|K|i> space.
const Ref< OrbitalSpace > & F_P_P(SpinCase1 S)
Form <P|F|P> space.
const Ref< OrbitalSpace > & hj_m_p(SpinCase1 S)
Form <p|h+J|m> space.
const Ref< OrbitalSpace > & K_i_a(SpinCase1 S)
Form <a|K|i> space.
const Ref< OrbitalSpace > & F_x_P(SpinCase1 S)
Form <P|F|x> space.
const Ref< OrbitalSpace > & K_p_m(SpinCase1 S)
Form <m|K|p> space.
const Ref< OrbitalSpace > & F_i_a(SpinCase1 S)
Form <a|F|i> space.
const Ref< OrbitalSpace > & hj_p_P(SpinCase1 S)
Form <P|h+J|p> space.
const Ref< OrbitalSpace > & K_i_m(SpinCase1 S)
Form <m|K|i> space.
const Ref< OrbitalSpace > & F_a_A(SpinCase1 S)
Form <A|F|a> space.
virtual void compute()
This function causes the intermediate matrices to be computed.
const Ref< OrbitalSpace > & K_p_p(SpinCase1 S)
Form <p|K|p> space.
const Ref< OrbitalSpace > & K_A_P(SpinCase1 S)
Form <P|K|A> space.
const Ref< AOSpaceRegistry > & ao_registry() const
returns the AOSpaceRegistry object
Definition r12int_eval.h:1007
const Ref< OrbitalSpace > & hj_a_A(SpinCase1 S)
Form <A|h+J|a> space.
const Ref< OrbitalSpace > & F_p_m(SpinCase1 S)
Form <m|F|p> space.
RefSymmSCMatrix ordm(SpinCase1 S) const
Returns the 1-RDM for spin S in the `‘MO’' basis (i.e. that provided by orbs(S) )
The RefSCDimension class is a smart pointer to an SCDimension specialization.
Definition dim.h:152
The RefSCMatrix class is a smart pointer to an SCMatrix specialization.
Definition matrix.h:135
The RefSCVector class is a smart pointer to an SCVector specialization.
Definition matrix.h:55
The RefSymmSCMatrix class is a smart pointer to an SCSymmSCMatrix specialization.
Definition matrix.h:265
A template class that maintains references counts.
Definition ref.h:361
Base class for objects that can save/restore state.
Definition state.h:45
Restores fundamental and user-defined types from images created with StateOut.
Definition statein.h:79
Serializes fundamental and user-defined types.
Definition stateout.h:71
Provides information about the type of a two body tensor.
Definition twobodytensorinfo.h:37
Definition reftestx.h:32
Contains all MPQC code up to version 3.
Definition mpqcin.h:14
void antisymmetrize(RefSCMatrix &Aanti, const RefSCMatrix &A, const Ref< OrbitalSpace > &bra, const Ref< OrbitalSpace > &ket, bool accumulate=false)
Antisymmetrizes 4-index quantity <ij|A|kl> -> <ij|A|kl> - <ij|A|lk> and saves to Aanti.
type
types of known two-body operators
Definition operator.h:318

Generated at Wed Sep 25 2024 02:45:29 for MPQC 3.0.0-alpha using the documentation package Doxygen 1.12.0.