MPQC 3.0.0-alpha
Loading...
Searching...
No Matches
mp2r12_energy.h
1//
2// mp2r12_energy.h
3//
4// Copyright (C) 2003 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#include <util/ref/ref.h>
29#include <chemistry/qc/mbptr12/r12technology.h>
30#include <chemistry/qc/mbptr12/r12int_eval.h>
31#include <chemistry/qc/wfn/spin.h>
32#include <chemistry/qc/mbptr12/twobodygrid.h>
33#include <math/mmisc/pairiter.h>
34#include <chemistry/qc/mbptr12/mp2r12_energy_util.h>
35
36#ifndef _chemistry_qc_mbptr12_mp2r12energy_h
37#define _chemistry_qc_mbptr12_mp2r12energy_h
38
39#define MP2R12ENERGY_CAN_COMPUTE_PAIRFUNCTION 1
40
41namespace sc {
46 class R12EnergyIntermediates : virtual public SavableState {
47 private:
48 R12Technology::StandardApproximation stdapprox_;
49 Ref<R12IntEval> r12eval_;
50 bool V_computed_;
51 RefSCMatrix V_[NSpinCases2];
52 bool X_computed_;
53 RefSymmSCMatrix X_[NSpinCases2];
54 bool B_computed_;
55 RefSymmSCMatrix B_[NSpinCases2];
56 bool A_computed_;
57 RefSCMatrix A_[NSpinCases2];
58 bool T1_cc_computed_;
59 RefSCMatrix T1_cc_[NSpinCases1];
60 bool T2_cc_computed_;
61 Ref<DistArray4> T2_cc_[NSpinCases2];
62 // lambda amplitudes
63 bool L1_cc_computed_;
64 RefSCMatrix L1_cc_[NSpinCases1];
65 bool L2_cc_computed_;
66 Ref<DistArray4> L2_cc_[NSpinCases2];
67 // parameters for importing psi ccsd one-particle density
68 bool Onerdm_cc_computed_;
69 RefSCMatrix Onerdm_cc_[NSpinCases1];
70 // parameters for orbial relaxation contribution to
71 // ccsd_f12 one-particle density
72 bool Onerdm_relax_computed_;
73 RefSCMatrix Onerdm_relax_[NSpinCases1];
74 public:
75 typedef enum { V=0, X=1, B=2, A=3 } IntermediateType;
77 const R12Technology::StandardApproximation stdapp);
81 Ref<R12IntEval> r12eval() const;
82 void set_r12eval(Ref<R12IntEval> &r12eval);
83 R12Technology::StandardApproximation stdapprox() const;
84 bool V_computed() const;
85 bool X_computed() const;
86 bool B_computed() const;
87 bool A_computed() const;
88 bool T1_cc_computed() const;
89 bool T2_cc_computed() const;
90 bool L1_cc_computed() const;
91 bool L2_cc_computed() const;
92 bool Onerdm_cc_computed() const;
93 bool Onerdm_relax_computed() const;
94 const RefSCMatrix& get_V(const SpinCase2 &spincase2) const;
95 void assign_V(const SpinCase2 &spincase2, const RefSCMatrix& V);
96 const RefSymmSCMatrix& get_X(const SpinCase2 &spincase2) const;
97 void assign_X(const SpinCase2 &spincase2, const RefSymmSCMatrix& X);
98 const RefSymmSCMatrix& get_B(const SpinCase2 &spincase2) const;
99 void assign_B(const SpinCase2 &spincase2, const RefSymmSCMatrix& B);
100 const RefSCMatrix& get_A(const SpinCase2 &spincase2) const;
101 void assign_A(const SpinCase2 &spincase2, const RefSCMatrix& A);
102 const RefSCMatrix& get_T1_cc(const SpinCase1 &spincase1) const;
103 void assign_T1_cc(const SpinCase1 &spincase1, const RefSCMatrix& T1_cc);
104 const Ref<DistArray4>& get_T2_cc(const SpinCase2 &spincase2) const;
105 void assign_T2_cc(const SpinCase2 &spincase2, const Ref<DistArray4>& T2_cc);
106 const RefSCMatrix& get_L1_cc(const SpinCase1 &spincase1) const;
107 void assign_L1_cc(const SpinCase1 &spincase1, const RefSCMatrix& L1_cc);
108 const Ref<DistArray4>& get_L2_cc(const SpinCase2 &spincase2) const;
109 void assign_L2_cc(const SpinCase2 &spincase2, const Ref<DistArray4>& L2_cc);
110 const RefSCMatrix& get_1rdm_cc(const SpinCase1 &spincase1) const;
111 void assign_1rdm_cc(const SpinCase1 &spincase1, const RefSCMatrix& Onerdm_cc);
112 const RefSCMatrix& get_1rdm_relax(const SpinCase1 &spincase1) const;
113 void assign_1rdm_relax(const SpinCase1 &spincase1, const RefSCMatrix& Onerdm_relax);
114 };
115
117class MP2R12Energy : virtual public SavableState {
118 private:
123 RefSCVector compute_2body_values_(bool equiv, const Ref<OrbitalSpace>& space1, const Ref<OrbitalSpace>& space2,
124 const SCVector3& r1, const SCVector3& r2) const;
125
126 // computes f12 contributions to the mp2 pair energies
127 virtual void compute_ef12() =0;
128
129 protected:
130 Ref<R12EnergyIntermediates> r12intermediates_;
131 Ref<R12IntEval> r12eval_;
132 bool include_obs_singles_;
133 int debug_;
134 bool evaluated_;
135
136 RefSCVector ef12_[NSpinCases2], emp2f12_[NSpinCases2];
137 // The coefficients are stored xy by ij, where xy is the geminal-multiplied pair
138 RefSCMatrix C_[NSpinCases2];
139
140 // Initialize SCVectors and SCMatrices
141 void init();
142
143 public:
144
146 MP2R12Energy(const Ref<R12EnergyIntermediates>& r12intermediates,
147 bool include_obs_singles,
148 int debug);
150
152 void obsolete();
153 void print(std::ostream&o=ExEnv::out0()) const;
154
155 Ref<R12IntEval> r12eval() const;
156 const Ref<R12EnergyIntermediates>& r12intermediates() const;
157 R12Technology::StandardApproximation stdapprox() const;
158 bool include_obs_singles() const { return include_obs_singles_; }
159 void set_debug(int debug);
160 int get_debug() const;
161
163 double energy();
164 const RefSCVector& ef12(SpinCase2 S);
165 const RefSCVector& emp2f12(SpinCase2 S);
166 double emp2f12tot(SpinCase2 S);
167 double ef12tot(SpinCase2 S);
168 void print_pair_energies(bool spinadapted,
169 double cabs_singles_energy,
170 std::ostream&so=ExEnv::out0());
171
172#if MP2R12ENERGY_CAN_COMPUTE_PAIRFUNCTION
174 void compute_pair_function(unsigned int i, unsigned int j, SpinCase2 spincase2,
175 const Ref<TwoBodyGrid>& tbgrid);
176#endif
177
179 void compute();
182 RefSCMatrix C(SpinCase2 S);
185 RefSCMatrix T2(SpinCase2 S);
186};
187
194{
195 void compute_ef12();
196
197 public:
200 bool include_obs_singles,
201 int debug);
203
205};
206
212{
213 void compute_ef12();
215 void compute_ef12_10132011();
216 void activate_ints(const std::string&, const std::string&,
217 const std::string&, const std::string&,
218 const std::string&, Ref<TwoBodyFourCenterMOIntsRuntime>&,
220 // Label Y^ij_ij, Y^ij_ji, Y^ji_ij, Y^ji_ji
221 enum idx_b1b2k1k2{ij_ij, ij_ji, ji_ij, ji_ji};
222 void compute_Y(const int b1b2_k1k2, const double prefactor,
223 const unsigned int oper_idx,
224 Ref<DistArray4>& i1i2i1i2_ints, double* array_i1i2i1i2);
225 void compute_YxF(const int b1b2_k1k2, const double prefactor,
226 const unsigned int oper1_idx, const unsigned int oper2_idx,
227 const Ref<DistArray4>& i1i2x1x2_ints, const Ref<DistArray4>& i2i1x1x2_ints,
228 double* array_ijij);
229 void compute_FxT(const int b1b2_k1k2, const unsigned int f12_idx,
230 Ref<DistArray4>& F_ints, const double* Tiiaa,
231 double* V_coupling);
232 void compute_VX(const int b1b2_k1k2, std::vector<std::string>& VX_output,
233 const unsigned int oper12_idx, Ref<DistArray4>& Y12_ints,
234 const unsigned int oper1_idx, const unsigned int oper2_idx,
235 const std::vector<Ref<DistArray4> >& Y_ints,
236 const std::vector<Ref<DistArray4> >& F_ints,
237 double* VX_array);
238 void accumulate_P_YxF(std::vector<std::string>& P_output,
239 std::vector<int>& b1b2_k1k2, std::vector<double>& P_prefactor,
240 const unsigned int oper1_idx, const unsigned int oper2_idx,
241 std::vector<Ref<DistArray4> >& Y_ints,
242 std::vector<Ref<DistArray4> >& F_ints,
243 double* P);
245 void compute_U(Ref<DistArray4> (&U)[NSpinCases2]);
246 void contract_VT1(const Ref<DistArray4>& V,
247 const int b1b2_k1k2, const bool swap_e12_V,
248 const double* const T1_array,
249 const int nv, const bool VT1_offset,
250 double* const VT1);
251
252#ifdef MPQC_NEW_FEATURES
254 double U2T2_ta();
256 double U1T1_plus_C1T1_ta();
258 double VTT_ta();
260 double PT2R12_ta();
261#endif
262
263 //
264 // compute the one electron density matrix for the diagonal ansatz
265 //RefSCMatrix D_ccsdf12_[NSpinCases1];
266 void compute_density_diag();
267
268 // functions needed for compute_density_diag() function:
269 void obtain_orbitals(const SpinCase2 spincase,
270 std::vector<Ref<OrbitalSpace> >& v_orbs1,
271 std::vector<Ref<OrbitalSpace> >& v_orbs2);
272
273 // activate three f12_ints for computing X
274 void activate_ints_X_f12(Ref<TwoBodyFourCenterMOIntsRuntime>& moints4_rtime, const std::string& index,
275 const std::vector< Ref<OrbitalSpace> >& v_orbs1,
276 const std::vector< Ref<OrbitalSpace> >& v_orbs2,
277 const std::string& descr_f12_key, std::vector<Ref<DistArray4> >& f12_ints);
278
279 // test function for D^i_i
280 void compute_Dii_test(const int nspincases1, const int nspincases2,
281 const int nocc_alpha, const int nocc_beta,
282 const double C_0, const double C_1);
283 // functions for compute_Dii_test:
284 // compute AlphaAlpha/BetaBeta: R^ij_ab R^ab_ij & R^ji_ab R^ab_ij
285 void compute_RRii_ii(std::vector<std::string>& output,
286 const std::vector< Ref<OrbitalSpace> >& v_orbs1,
287 const std::vector< Ref<OrbitalSpace> >& v_orbs2,
288 double* RRij_ij, double* RRji_ij);
289 // compute the openshell AlphaBeta:
290 // R^ij_ab R^ab_ij, R^ji_ab R^ab_ij, R^ij_ab R^ab_ji, & R^ji_ab R^ab_ji
291 void compute_Rii_ii(std::vector<std::string>& output,
292 const std::vector< Ref<OrbitalSpace> >& v_orbs1,
293 const std::vector< Ref<OrbitalSpace> >& v_orbs2,
294 double* RRi1i2_i1i2, double* RRi1i2_i2i1,
295 double* RRi2i1_i1i2, double* RRi2i1_i2i1);
296
297 // compute D^m_i
298 void compute_Dmi(const int nspincases1, const int nspincases2,
299 const double C_0, const double C_1,
300 const std::vector< Ref<OrbitalSpace> >& v_orbs1_ab,
301 const std::vector< Ref<OrbitalSpace> >& v_orbs2_ab,
302 double* const Dm_i_alpha, double* const Dm_i_beta);
303 // functions needed for computing D^m_i :
304 // compute R^ab_b1b2 R^k1k2 _ab (a, b represent complete virtual orbitals)
305 // sum over a, b, and index 3
306 // index 1,2: i, m; index 3: j e.g.: R^ab_ij R^jm_ab
307 void compute_RR_sum_abj2(const int RRb1b2_k1k2,
308 const int f12f12_idx, const int f12_idx,
309 const Ref<DistArray4>& f12f12_ints,
310 const std::vector< Ref<DistArray4> >& v_f12_ints1,
311 const std::vector< Ref<DistArray4> >& v_f12_ints2,
312 double* const RR_result);
313
314 // compute D^m_i through R^ij_a'b' R_ij^a'b' + 2 R^ij_ab' R_ij^ab'
315 void compute_Dmi_2(const int nspincases1, const int nspincases2,
316 const double C_0, const double C_1,
317 const std::vector< Ref<OrbitalSpace> >& v_orbs1_ab,
318 const std::vector< Ref<OrbitalSpace> >& v_orbs2_ab,
319 double* const Dm_i_alpha, double* const Dm_i_beta);
320
321 // test function fo D^c_b
322 void compute_Dcb(const int nspincases1, const int nspincases2,
323 const double C_0, const double C_1);
324
325 // test function for D^c'_b'
326 void compute_Dcpbp_a(const int nspincases1, const int nspincases2,
327 const double C_0, const double C_1);
328
329 // test function for D^c'_b' sum over a'
330 void compute_Dcpbp_ap(const int nspincases1, const int nspincases2,
331 const double C_0, const double C_1);
332
333 // test function for D^a'_a RR part
334 void compute_Dapa_RR(const int nspincases1, const int nspincases2,
335 const double C_0, const double C_1);
336
337 // \bar{\tilde{R}}^31_ij \bar{\tilde{R}}^ij_32 with spin orbitals
338 // which is for computing D^b'_c', D^a'_a
339 void compute_RR31_32_spin(const int orbitals_label,
340 const int nspincases1, const int nspincases2,
341 const double C_0, const double C_1,
342 const std::vector< Ref<OrbitalSpace> >& v_orbs1_ab,
343 const std::vector< Ref<OrbitalSpace> >& v_orbs2_ab,
344 double* const D_alhpha, double* const D_beta);
345
346 // function need for D^a'_a:
347 // compute R * T2 (CC amplitudes)
348 void compute_RT2_apa(const int nspincases1, const int nspincases2,
349 const double C_0, const double C_1,
350 const std::vector< Ref<OrbitalSpace> >& v_orbs1_ab,
351 const std::vector< Ref<OrbitalSpace> >& v_orbs2_ab,
352 double* RT2_alpha, double* RT2_beta);
353
354 // compute MP2 T2 amplitude (non-antisymmetrized)
355 void compute_T2_mp2(const std::vector< Ref<OrbitalSpace> >& v_orbs1,
356 const std::vector< Ref<OrbitalSpace> >& v_orbs2,
357 double* T2ab_ij);
358 // compute R * T2 (mp2 amplitudes)
359 void compute_RTmp2_apa(const int nspincases1, const int nspincases2,
360 const double C_0, const double C_1,
361 const std::vector< Ref<OrbitalSpace> >& v_orbs1_ab,
362 const std::vector< Ref<OrbitalSpace> >& v_orbs2_ab,
363 double* const RT2_alpha, double* const RT2_beta);
364
365 RefSCMatrix compute_D_CABS(SpinCase1 spin);
366 RefSCMatrix compute_D_CABS_test(SpinCase1 spin);
367 // compute MP2 one-electron density matrix
368 RefSCMatrix compute_1rdm_mp2(const SpinCase1 spin);
369 RefSCMatrix compute_1rdm_mp2_test(const SpinCase1 spin);
370 // MP2F12 one-electron density matrix
371 // D_MP2F12 = T(MP2)T(MP2) + 2 T(MP2)T(F12) + T(F12)T(F12)
372 void compute_1rdm_mp2f12(const int nspincases1, const int nspincases2,
373 const int C_0, const int C_1,
374 RefSCMatrix Dmp2f12[NSpinCases1]);
375 // compute contribution for MP2F12 one-electron density matrix: D_MP2F12
376 // compute T(MP2)T(MP2) or T(F12)T(F12)
377 RefSCMatrix compute_1rdm_mp2part(const SpinCase1 spin,
378 const int nocc1_act, const int nocc2_act,
379 const int nvir1, const int nvir2,
380 const double* const T2, const double* const T2_ab);
381 // compute T(MP2)T(F12) or T(MP2)T(F12)
382 RefSCMatrix compute_1rdm_mp2part(const SpinCase1 spin,
383 const int nocc1_act, const int nocc2_act,
384 const int nvir1, const int nvir2,
385 const double* const T2_left,const double* const T2_right,
386 const double* const T2_ab_left, const double* const T2_ab_right);
387 // compute MP2 T2 amplitude (antisymmetrized), stored in array (a,b,i,j)
388 void compute_T2_mp2(const SpinCase2 spincase,
389 double* const T2ab_ij);
390 // compute F12 corrected T2 amplitude (antisymmetrized), stored in array (a,b,i,j)
391 void compute_T2abij_f12corr(const SpinCase2 spincase,
392 const double C_0, const double C_1,
393 double* const T2ab_ij_f12corr);
394 // compute MP2F12 T2 amplitude = T2 (MP2) + T2 (F12 corrected)
395 void compute_T2abij_mp2f12(const int nocc1_act, const int nocc2_act,
396 const int nvir1, const int nvir2,
397 const double* const T2ab_ij_mp2,
398 const double* const T2ab_ij_f12corr,
399 double* const T2abij_mp2f12);
400 // compute MP2F12 T2 amplitude, stored in array (a,b,i,j)
401 void compute_T2abij_mp2f12(const SpinCase2 spincase,
402 const double C_0, const double C_1,
403 double* const T2ab_ij_mp2f12);
404 // transform the one-particle density matrix into the MPQC ordering
405 RefSCMatrix onepdm_transformed(const SpinCase1& spin, const bool frozen_core, const RefSCMatrix& D);
406 // test function for computing CCSD dipole moment
407 RefSCMatrix onepdm_transformed2(const SpinCase1& spin,const RefSCMatrix& D);
408 // compute D^a'_i = R * T1 (CC)
409 void compute_RT1_api(const int nspincases1, const int nspincases2,
410 const double C_0, const double C_1,
411 const std::vector< Ref<OrbitalSpace> >& v_orbs1_ab,
412 const std::vector< Ref<OrbitalSpace> >& v_orbs2_ab,
413 double* const D_alpha, double* const D_beta);
414
415 public:
418 bool include_obs_singles,
419 int debug);
421
423};
424
425Ref<MP2R12Energy> construct_MP2R12Energy(Ref<R12EnergyIntermediates> &r12intermediates,
426 bool include_obs_singles,
427 int debug,
428 bool diag);
429
430}
431
432#endif
433
434// Local Variables:
435// mode: c++
436// c-file-style: "CLJ"
437// End:
438
439
static std::ostream & out0()
Return an ostream that writes from node 0.
The class MP2R12Energy_Diag is an implementation of MP2R12Energy that supports Ten-no's diagonal orbi...
Definition mp2r12_energy.h:212
void save_data_state(StateOut &)
Save the base classes (with save_data_state) and the members in the same order that the StateIn CTOR ...
The class MP2R12Energy_SpinOrbital is the original implementation of MP2R12Energy It supports only th...
Definition mp2r12_energy.h:194
void save_data_state(StateOut &)
Save the base classes (with save_data_state) and the members in the same order that the StateIn CTOR ...
Class MP2R12Energy is the object that computes and maintains MP2-R12 energies.
Definition mp2r12_energy.h:117
void print(std::ostream &o=ExEnv::out0()) const
Print the object.
double energy()
total correlation energy (including OBS singles if include_obs_singles() == true)
void save_data_state(StateOut &)
Save the base classes (with save_data_state) and the members in the same order that the StateIn CTOR ...
void compute()
Computes the first-order R12 wave function and MP2-R12 energy.
RefSCMatrix C(SpinCase2 S)
Returns the matrix of first-order amplitudes of r12-multiplied occupied orbital pairs.
RefSCMatrix T2(SpinCase2 S)
Returns the matrix of first-order amplitudes of conventional orbital pairs.
The class R12EnergyIntermediates is the front-end to R12 intermediates.
Definition mp2r12_energy.h:46
void save_data_state(StateOut &so)
Save the base classes (with save_data_state) and the members in the same order that the StateIn CTOR ...
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
a 3-element version of SCVector
Definition vector3.h:44
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
Definition reftestx.h:32
Contains all MPQC code up to version 3.
Definition mpqcin.h:14
unsigned int nspincases2(bool spin_polarized)
Returns the number of unique combinations of 2 spin cases (1 or 3)
unsigned int nspincases1(bool spin_polarized)
Returns the number of unique spin cases (1 or 2)

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