1 /*
2    CheMPS2: a spin-adapted implementation of DMRG for ab initio quantum chemistry
3    Copyright (C) 2013-2018 Sebastian Wouters
4 
5    This program is free software; you can redistribute it and/or modify
6    it under the terms of the GNU General Public License as published by
7    the Free Software Foundation; either version 2 of the License, or
8    (at your option) any later version.
9 
10    This program is distributed in the hope that it will be useful,
11    but WITHOUT ANY WARRANTY; without even the implied warranty of
12    MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
13    GNU General Public License for more details.
14 
15    You should have received a copy of the GNU General Public License along
16    with this program; if not, write to the Free Software Foundation, Inc.,
17    51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
18 */
19 
20 #include <stdlib.h>
21 #include <iostream>
22 #include <math.h>
23 #include <assert.h>
24 #include <algorithm>
25 #include <sys/time.h>
26 
27 #include "CASPT2.h"
28 #include "Lapack.h"
29 #include "Options.h"
30 #include "ConjugateGradient.h"
31 #include "Davidson.h"
32 #include "Special.h"
33 
34 using std::cout;
35 using std::endl;
36 using std::min;
37 using std::max;
38 
CASPT2(DMRGSCFindices * idx,DMRGSCFintegrals * ints,DMRGSCFmatrix * oei,DMRGSCFmatrix * fock_in,double * one_dm,double * two_dm,double * three_dm,double * contract_4dm,const double IPEA)39 CheMPS2::CASPT2::CASPT2( DMRGSCFindices * idx, DMRGSCFintegrals * ints, DMRGSCFmatrix * oei, DMRGSCFmatrix * fock_in, double * one_dm, double * two_dm, double * three_dm, double * contract_4dm, const double IPEA ){
40 
41    indices    = idx;
42    fock       = fock_in;
43    one_rdm    = one_dm;
44    two_rdm    = two_dm;
45    three_rdm  = three_dm;
46    f_dot_4dm  = contract_4dm;
47    num_irreps = indices->getNirreps();
48 
49    struct timeval start, end;
50    gettimeofday( &start, NULL );
51 
52    create_f_dots();
53    vector_helper();
54 
55    make_AA_CC( true, 0.0 );
56    make_DD( true, 0.0 );
57    make_EE_GG( true, 0.0 );
58    make_BB_FF_singlet( true, 0.0 );
59    make_BB_FF_triplet( true, 0.0 );
60 
61    construct_rhs( oei, ints );
62 
63    make_AA_CC( false, IPEA );
64    make_DD( false, IPEA );
65    make_EE_GG( false, IPEA );
66    make_BB_FF_singlet( false, IPEA );
67    make_BB_FF_triplet( false, IPEA );
68    make_FAD_FCD();
69    make_FEH_FGH();
70    make_FAB_FCF_singlet();
71    make_FAB_FCF_triplet();
72    make_FBE_FFG_singlet();
73    make_FBE_FFG_triplet();
74    make_FDE_FDG();
75 
76    delete [] f_dot_3dm;
77    delete [] f_dot_2dm;
78 
79    gettimeofday( &end, NULL );
80    double elapsed = ( end.tv_sec - start.tv_sec ) + 1e-6 * ( end.tv_usec - start.tv_usec );
81    cout << "CASPT2 : Wall time tensors    = " << elapsed << " seconds" << endl;
82    gettimeofday( &start, NULL );
83 
84    recreate();
85 
86    gettimeofday( &end, NULL );
87    elapsed = ( end.tv_sec - start.tv_sec ) + 1e-6 * ( end.tv_usec - start.tv_usec );
88    cout << "CASPT2 : Wall time diag(ovlp) = " << elapsed << " seconds" << endl;
89 
90 }
91 
~CASPT2()92 CheMPS2::CASPT2::~CASPT2(){
93 
94    for ( int irrep = 0; irrep < num_irreps; irrep++ ){
95       delete [] FAA[ irrep ];
96       delete [] FCC[ irrep ];
97       delete [] FDD[ irrep ];
98       delete [] FEE[ irrep ];
99       delete [] FGG[ irrep ];
100       delete [] FBB_singlet[ irrep ];
101       delete [] FBB_triplet[ irrep ];
102       delete [] FFF_singlet[ irrep ];
103       delete [] FFF_triplet[ irrep ];
104    }
105    delete [] FAA;
106    delete [] FCC;
107    delete [] FDD;
108    delete [] FEE;
109    delete [] FGG;
110    delete [] FBB_singlet;
111    delete [] FBB_triplet;
112    delete [] FFF_singlet;
113    delete [] FFF_triplet;
114 
115    for ( int irrep_left = 0; irrep_left < num_irreps; irrep_left++ ){
116       for ( int irrep_right = 0; irrep_right < num_irreps; irrep_right++ ){
117          const int irrep_w = Irreps::directProd( irrep_left, irrep_right );
118          const int num_w   = indices->getNDMRG( irrep_w );
119          for ( int w = 0; w < num_w; w++ ){
120             delete [] FAD[ irrep_left ][ irrep_right ][ w ];
121             delete [] FCD[ irrep_left ][ irrep_right ][ w ];
122             delete [] FAB_singlet[ irrep_left ][ irrep_right ][ w ];
123             delete [] FAB_triplet[ irrep_left ][ irrep_right ][ w ];
124             delete [] FCF_singlet[ irrep_left ][ irrep_right ][ w ];
125             delete [] FCF_triplet[ irrep_left ][ irrep_right ][ w ];
126             delete [] FBE_singlet[ irrep_left ][ irrep_right ][ w ];
127             delete [] FBE_triplet[ irrep_left ][ irrep_right ][ w ];
128             delete [] FFG_singlet[ irrep_left ][ irrep_right ][ w ];
129             delete [] FFG_triplet[ irrep_left ][ irrep_right ][ w ];
130             delete [] FDE_singlet[ irrep_left ][ irrep_right ][ w ];
131             delete [] FDE_triplet[ irrep_left ][ irrep_right ][ w ];
132             delete [] FDG_singlet[ irrep_left ][ irrep_right ][ w ];
133             delete [] FDG_triplet[ irrep_left ][ irrep_right ][ w ];
134          }
135          delete [] FAD[ irrep_left ][ irrep_right ];
136          delete [] FCD[ irrep_left ][ irrep_right ];
137          delete [] FAB_singlet[ irrep_left ][ irrep_right ];
138          delete [] FAB_triplet[ irrep_left ][ irrep_right ];
139          delete [] FCF_singlet[ irrep_left ][ irrep_right ];
140          delete [] FCF_triplet[ irrep_left ][ irrep_right ];
141          delete [] FBE_singlet[ irrep_left ][ irrep_right ];
142          delete [] FBE_triplet[ irrep_left ][ irrep_right ];
143          delete [] FFG_singlet[ irrep_left ][ irrep_right ];
144          delete [] FFG_triplet[ irrep_left ][ irrep_right ];
145          delete [] FDE_singlet[ irrep_left ][ irrep_right ];
146          delete [] FDE_triplet[ irrep_left ][ irrep_right ];
147          delete [] FDG_singlet[ irrep_left ][ irrep_right ];
148          delete [] FDG_triplet[ irrep_left ][ irrep_right ];
149       }
150       delete [] FAD[ irrep_left ];
151       delete [] FCD[ irrep_left ];
152       delete [] FAB_singlet[ irrep_left ];
153       delete [] FAB_triplet[ irrep_left ];
154       delete [] FCF_singlet[ irrep_left ];
155       delete [] FCF_triplet[ irrep_left ];
156       delete [] FBE_singlet[ irrep_left ];
157       delete [] FBE_triplet[ irrep_left ];
158       delete [] FFG_singlet[ irrep_left ];
159       delete [] FFG_triplet[ irrep_left ];
160       delete [] FDE_singlet[ irrep_left ];
161       delete [] FDE_triplet[ irrep_left ];
162       delete [] FDG_singlet[ irrep_left ];
163       delete [] FDG_triplet[ irrep_left ];
164    }
165    delete [] FAD;
166    delete [] FCD;
167    delete [] FAB_singlet;
168    delete [] FAB_triplet;
169    delete [] FCF_singlet;
170    delete [] FCF_triplet;
171    delete [] FBE_singlet;
172    delete [] FBE_triplet;
173    delete [] FFG_singlet;
174    delete [] FFG_triplet;
175    delete [] FDE_singlet;
176    delete [] FDE_triplet;
177    delete [] FDG_singlet;
178    delete [] FDG_triplet;
179 
180    for ( int irrep = 0; irrep < num_irreps; irrep++ ){
181       const int num_w = indices->getNDMRG( irrep );
182       for ( int w = 0; w < num_w; w++ ){
183          delete [] FEH[ irrep ][ w ];
184          delete [] FGH[ irrep ][ w ];
185       }
186       delete [] FEH[ irrep ];
187       delete [] FGH[ irrep ];
188    }
189    delete [] FEH;
190    delete [] FGH;
191 
192    delete [] size_A;
193    delete [] size_C;
194    delete [] size_D;
195    delete [] size_E;
196    delete [] size_G;
197    delete [] size_B_singlet;
198    delete [] size_B_triplet;
199    delete [] size_F_singlet;
200    delete [] size_F_triplet;
201    delete [] jump;
202    delete [] vector_rhs;
203 
204 }
205 
solve(const double imag_shift,const bool CONJUGATE_GRADIENT) const206 double CheMPS2::CASPT2::solve( const double imag_shift, const bool CONJUGATE_GRADIENT ) const{
207 
208    struct timeval start, end;
209    gettimeofday( &start, NULL );
210 
211    // Normalizations of sectors   A  Bs Bt C  D  Es Et Fs Ft Gs Gt Hs Ht
212    const int normalizations[] = { 1, 2, 2, 1, 1, 2, 6, 2, 2, 2, 6, 4, 12 };
213    const bool apply_shift = (( fabs( imag_shift ) > 0.0 ) ? true : false );
214 
215    int total_size = jump[ CHEMPS2_CASPT2_NUM_CASES * num_irreps ];
216    double * diag_fock = new double[ total_size ];
217    diagonal( diag_fock );
218    double min_eig = diag_fock[ 0 ];
219    for ( int elem = 1; elem < total_size; elem++ ){ min_eig = min( min_eig, diag_fock[ elem ] ); }
220    cout << "CASPT2 : Solution algorithm   = " << (( CONJUGATE_GRADIENT ) ? "Conjugate Gradient" : "Davidson" ) << endl;
221    cout << "CASPT2 : Minimum(diagonal)    = " << min_eig << endl;
222 
223    ConjugateGradient * CG = (( CONJUGATE_GRADIENT ) ? new ConjugateGradient( total_size, CheMPS2::CONJ_GRADIENT_RTOL, CheMPS2::CONJ_GRADIENT_PRECOND_CUTOFF, false ) : NULL );
224    Davidson * DAVID = (( CONJUGATE_GRADIENT ) ? NULL : new Davidson( total_size,
225                                                                      CheMPS2::DAVIDSON_NUM_VEC,
226                                                                      CheMPS2::DAVIDSON_NUM_VEC_KEEP,
227                                                                      CheMPS2::CONJ_GRADIENT_RTOL,
228                                                                      CheMPS2::CONJ_GRADIENT_PRECOND_CUTOFF,
229                                                                      true, // debug_print
230                                                                      'L' )); // Linear problem
231    double ** pointers = new double*[ 3 ];
232    char instruction = (( CONJUGATE_GRADIENT ) ? CG->step( pointers ) : DAVID->FetchInstruction( pointers ));
233    assert( instruction == 'A' );
234    for ( int elem = 0; elem < total_size; elem++ ){ pointers[ 0 ][ elem ] = vector_rhs[ elem ] / diag_fock[ elem ]; } // Initial guess of F * x = V
235    for ( int elem = 0; elem < total_size; elem++ ){ pointers[ 1 ][ elem ] =  diag_fock[ elem ]; } // Diagonal of the operator F
236    for ( int elem = 0; elem < total_size; elem++ ){ pointers[ 2 ][ elem ] = vector_rhs[ elem ]; } // RHS of the linear problem F * x = V
237    int inc1 = 1;
238    const double E2_DIAGONAL = - ddot_( &total_size, pointers[ 0 ], &inc1, pointers[ 2 ], &inc1 );
239    instruction = (( CONJUGATE_GRADIENT ) ? CG->step( pointers ) : DAVID->FetchInstruction( pointers ));
240    assert( instruction == 'B' );
241    while ( instruction == 'B' ){
242       matvec( pointers[ 0 ], pointers[ 1 ], diag_fock );
243       if ( apply_shift ){ add_shift( pointers[ 0 ], pointers[ 1 ], diag_fock, imag_shift, normalizations ); }
244       instruction = (( CONJUGATE_GRADIENT ) ? CG->step( pointers ) : DAVID->FetchInstruction( pointers ));
245    }
246    assert( instruction == 'C' );
247    const double E2_NONVARIATIONAL = - ddot_( &total_size, pointers[ 0 ], &inc1, vector_rhs, &inc1 );
248    const double rnorm = pointers[ 1 ][ 0 ];
249    cout << "CASPT2 : Number of iterations = " << (( CONJUGATE_GRADIENT ) ? CG->get_num_matvec() : DAVID->GetNumMultiplications() ) << endl;
250    cout << "CASPT2 : Residual norm        = " << rnorm << endl;
251    matvec( pointers[ 0 ], pointers[ 1 ], diag_fock ); // pointers[ 1 ] is a WORK array when instruction == 'C'
252    const double E2_VARIATIONAL = 2 * E2_NONVARIATIONAL + ddot_( &total_size, pointers[ 0 ], &inc1, pointers[ 1 ], &inc1 );
253    delete [] diag_fock;
254 
255    const double inproduct = inproduct_vectors( pointers[ 0 ], pointers[ 0 ], normalizations );
256    const double reference_weight = 1.0 / ( 1.0 + inproduct );
257    cout << "CASPT2 : Reference weight     = " << reference_weight << endl;
258    //energy_per_sector( pointers[ 0 ] );
259    delete [] pointers;
260    if ( CG    != NULL ){ delete CG;    }
261    if ( DAVID != NULL ){ delete DAVID; }
262 
263    gettimeofday( &end, NULL );
264    double elapsed = ( end.tv_sec - start.tv_sec ) + 1e-6 * ( end.tv_usec - start.tv_usec );
265    cout << "CASPT2 : Wall time solution   = " << elapsed << " seconds" << endl;
266 
267    cout << "CASPT2 : E2 [DIAGONAL]        = " << E2_DIAGONAL << endl;
268    cout << "CASPT2 : E2 [NON-VARIATIONAL] = " << E2_NONVARIATIONAL << endl;
269    cout << "CASPT2 : E2 [VARIATIONAL]     = " << E2_VARIATIONAL << endl;
270    return E2_VARIATIONAL;
271 
272 }
273 
add_shift(double * vector,double * result,double * diag_fock,const double imag_shift,const int * normalizations) const274 void CheMPS2::CASPT2::add_shift( double * vector, double * result, double * diag_fock, const double imag_shift, const int * normalizations ) const{
275 
276    for ( int sector = 0; sector < CHEMPS2_CASPT2_NUM_CASES; sector++ ){
277       const int start = jump[ num_irreps * sector         ];
278       const int stop  = jump[ num_irreps * ( sector + 1 ) ];
279       const double factor = imag_shift * imag_shift * normalizations[ sector ] * normalizations[ sector ];
280       for ( int elem = start; elem < stop; elem ++ ){
281          result[ elem ] += factor * vector[ elem ] / diag_fock[ elem ];
282       }
283    }
284 
285 }
286 
inproduct_vectors(double * first,double * second,const int * normalizations) const287 double CheMPS2::CASPT2::inproduct_vectors( double * first, double * second, const int * normalizations ) const{
288 
289    int inc1 = 1;
290    double value = 0.0;
291    for ( int sector = 0; sector < CHEMPS2_CASPT2_NUM_CASES; sector++ ){
292       int pointer = jump[ num_irreps * sector         ];
293       int size    = jump[ num_irreps * ( sector + 1 ) ] - pointer;
294       value += normalizations[ sector ] * ddot_( &size, first + pointer, &inc1, second + pointer, &inc1 );
295    }
296    return value;
297 
298 }
299 
energy_per_sector(double * solution) const300 void CheMPS2::CASPT2::energy_per_sector( double * solution ) const{
301 
302    int inc1 = 1;
303    double energies[ CHEMPS2_CASPT2_NUM_CASES ];
304    for ( int sector = 0; sector < CHEMPS2_CASPT2_NUM_CASES; sector++ ){
305       int pointer = jump[ num_irreps * sector         ];
306       int size    = jump[ num_irreps * ( sector + 1 ) ] - pointer;
307       energies[ sector ] = - ddot_( &size, solution + pointer, &inc1, vector_rhs + pointer, &inc1 );
308    }
309    cout << "************************************************" << endl;
310    cout << "*   CASPT2 non-variational energy per sector   *" << endl;
311    cout << "************************************************" << endl;
312    cout << "   A or VJTU = " << energies[ CHEMPS2_CASPT2_A ] << endl;
313    cout << "   B or VJTI = " << energies[ CHEMPS2_CASPT2_B_SINGLET ] + energies[ CHEMPS2_CASPT2_B_TRIPLET ] << endl;
314    cout << "   C or ATVX = " << energies[ CHEMPS2_CASPT2_C ] << endl;
315    cout << "   D or AIVX = " << energies[ CHEMPS2_CASPT2_D ] << endl;
316    cout << "   E or VJAI = " << energies[ CHEMPS2_CASPT2_E_SINGLET ] + energies[ CHEMPS2_CASPT2_E_TRIPLET ] << endl;
317    cout << "   F or BVAT = " << energies[ CHEMPS2_CASPT2_F_SINGLET ] + energies[ CHEMPS2_CASPT2_F_TRIPLET ] << endl;
318    cout << "   G or BJAT = " << energies[ CHEMPS2_CASPT2_G_SINGLET ] + energies[ CHEMPS2_CASPT2_G_TRIPLET ] << endl;
319    cout << "   H or BJAI = " << energies[ CHEMPS2_CASPT2_H_SINGLET ] + energies[ CHEMPS2_CASPT2_H_TRIPLET ] << endl;
320    cout << "************************************************" << endl;
321 
322 }
323 
create_f_dots()324 void CheMPS2::CASPT2::create_f_dots(){
325 
326    const int LAS = indices->getDMRGcumulative( num_irreps );
327    f_dot_3dm = new double[ LAS * LAS * LAS * LAS ];
328    f_dot_2dm = new double[ LAS * LAS ];
329    f_dot_1dm = 0.0;
330 
331    for ( int irrep = 0; irrep < num_irreps; irrep++ ){
332       const int NDMRG = indices->getNDMRG( irrep );
333       const int NOCC  = indices->getNOCC( irrep );
334       const int jumpx = indices->getDMRGcumulative( irrep );
335 
336       double value = 0.0;
337       for ( int orb = 0; orb < NDMRG; orb++ ){
338          value += fock->get( irrep, NOCC + orb, NOCC + orb ) * one_rdm[ jumpx + orb + LAS * ( jumpx + orb ) ];
339       }
340       f_dot_1dm += value;
341    }
342 
343    for ( int irrep1 = 0; irrep1 < num_irreps; irrep1++ ){
344       const int lower1 = indices->getDMRGcumulative( irrep1 );
345       const int upper1 = lower1 + indices->getNDMRG( irrep1 );
346 
347       for ( int i1 = lower1; i1 < upper1; i1++ ){
348          for ( int i2 = lower1; i2 < upper1; i2++ ){
349 
350             double value = 0.0;
351             for ( int irrepx = 0; irrepx < num_irreps; irrepx++ ){
352                const int NOCCx  = indices->getNOCC( irrepx );
353                const int NDMRGx = indices->getNDMRG( irrepx );
354                const int jumpx  = indices->getDMRGcumulative( irrepx );
355 
356                for ( int orbx = 0; orbx < NDMRGx; orbx++ ){
357                   value += fock->get( irrepx, NOCCx + orbx, NOCCx + orbx ) * two_rdm[ i1 + LAS * ( jumpx + orbx + LAS * ( i2 + LAS * ( jumpx + orbx ))) ];
358                }
359             }
360             f_dot_2dm[ i1 + LAS * i2 ] = value;
361          }
362       }
363    }
364 
365    for ( int irrep1 = 0; irrep1 < num_irreps; irrep1++ ){
366       const int lower1 = indices->getDMRGcumulative( irrep1 );
367       const int upper1 = lower1 + indices->getNDMRG( irrep1 );
368       for ( int irrep2 = 0; irrep2 < num_irreps; irrep2++ ){
369          const int lower2 = indices->getDMRGcumulative( irrep2 );
370          const int upper2 = lower2 + indices->getNDMRG( irrep2 );
371          const int irr_12 = Irreps::directProd( irrep1, irrep2 );
372          for ( int irrep3 = 0; irrep3 < num_irreps; irrep3++ ){
373             const int irrep4 = Irreps::directProd( irr_12, irrep3 );
374             const int lower3 = indices->getDMRGcumulative( irrep3 );
375             const int lower4 = indices->getDMRGcumulative( irrep4 );
376             const int upper3 = lower3 + indices->getNDMRG( irrep3 );
377             const int upper4 = lower4 + indices->getNDMRG( irrep4 );
378 
379             for ( int i1 = lower1; i1 < upper1; i1++ ){
380                for ( int i2 = lower2; i2 < upper2; i2++ ){
381                   for ( int i3 = lower3; i3 < upper3; i3++ ){
382                      for ( int i4 = lower4; i4 < upper4; i4++ ){
383 
384                         double value = 0.0;
385                         for ( int irrepx = 0; irrepx < num_irreps; irrepx++ ){
386                            const int jumpx  = indices->getDMRGcumulative( irrepx );
387                            const int NOCCx  = indices->getNOCC( irrepx );
388                            const int NDMRGx = indices->getNDMRG( irrepx );
389 
390                            for ( int orbx = 0; orbx < NDMRGx; orbx++ ){
391                               value += ( fock->get( irrepx, NOCCx + orbx, NOCCx + orbx )
392                                        * three_rdm[ i1 + LAS*( i2 + LAS*( jumpx + orbx + LAS*( i3 + LAS*( i4 + LAS*( jumpx + orbx ))))) ] );
393                            }
394                         }
395                         f_dot_3dm[ i1 + LAS * ( i2 + LAS * ( i3 + LAS * i4 )) ] = value;
396                      }
397                   }
398                }
399             }
400          }
401       }
402    }
403 
404    /*
405    double sum_f_kk = 0.0;
406    for ( int irrep = 0; irrep < num_irreps; irrep++ ){
407       const int NOCC = indices->getNOCC( irrep );
408       for ( int orb = 0; orb < NOCC; orb++ ){
409          sum_f_kk += fock->get( irrep, orb, orb );
410       }
411    }
412 
413    const double E_FOCK = 2 * sum_f_kk + f_dot_1dm;
414    cout << "CASPT2 : < 0 | F | 0 >        = " << E_FOCK << endl;
415    */
416 
417 }
418 
vector_helper()419 int CheMPS2::CASPT2::vector_helper(){
420 
421    int * helper = new int[ CHEMPS2_CASPT2_NUM_CASES * num_irreps ];
422 
423    /*** Type A : c_tiuv E_ti E_uv | 0 >
424                  c_tiuv = vector[ jump[ irrep + num_irreps * CHEMPS2_CASPT2_A ] + count_tuv + size_A[ irrep ] * count_i ]
425         Type C : c_atuv E_at E_uv | 0 >
426                  c_atuv = vector[ jump[ irrep + num_irreps * CHEMPS2_CASPT2_C ] + count_tuv + size_C[ irrep ] * count_a ]
427 
428         1/ (TYPE A) count_i = 0 .. NOCC[ irrep ]
429            (TYPE C) count_a = 0 .. NVIRT[ irrep ]
430         2/ jump_tuv = 0
431            irrep_t = 0 .. num_irreps
432               irrep_u = 0 .. num_irreps
433                  irrep_v = irrep x irrep_t x irrep_u
434                     ---> count_tuv = jump_tuv + t + NDMRG[ irrep_t ] * ( u + NDMRG[ irrep_u ] * v )
435                     jump_tuv += NDMRG[ irrep_t ] * NDMRG[ irrep_u ] * NDMRG[ irrep_v ]
436    */
437 
438    size_A = new int[ num_irreps ];
439    size_C = new int[ num_irreps ];
440    for ( int irrep = 0; irrep < num_irreps; irrep++ ){
441       int linsize_AC = 0;
442       for ( int irrep_t = 0; irrep_t < num_irreps; irrep_t++ ){
443          for ( int irrep_u = 0; irrep_u < num_irreps; irrep_u++ ){
444             const int irrep_v = Irreps::directProd( Irreps::directProd( irrep, irrep_t ), irrep_u );
445             linsize_AC += indices->getNDMRG( irrep_t ) * indices->getNDMRG( irrep_u ) * indices->getNDMRG( irrep_v );
446          }
447       }
448       size_A[ irrep ] = linsize_AC;
449       size_C[ irrep ] = linsize_AC;
450       helper[ irrep + num_irreps * CHEMPS2_CASPT2_A ] = size_A[ irrep ] * indices->getNOCC( irrep );
451       helper[ irrep + num_irreps * CHEMPS2_CASPT2_C ] = size_C[ irrep ] * indices->getNVIRT( irrep );
452    }
453 
454    /*** Type D1 : c1_aitu E_ai E_tu | 0 >
455         Type D2 : c2_tiau E_ti E_au | 0 >
456 
457                   c1_aitu = vector[ jump[ irrep + num_irreps * CHEMPS2_CASPT2_D ] + count_tu +          size_D[ irrep ] * count_ai ]
458                   c2_tiau = vector[ jump[ irrep + num_irreps * CHEMPS2_CASPT2_D ] + count_tu + D2JUMP + size_D[ irrep ] * count_ai ]
459 
460                   D2JUMP = size_D[ irrep ] / 2
461 
462         1/ jump_ai = 0
463            irrep_i = 0 .. num_irreps
464               irrep_a = irrep_i x irrep
465               ---> count_ai = jump_ai + i + NOCC[ irrep_i ] * a
466               jump_ai += NOCC[ irrep_i ] * NVIRT[ irrep_a ]
467         2/ jump_tu = 0
468            irrep_t = 0 .. num_irreps
469               irrep_u = irrep x irrep_t
470               ---> count_tu = jump_tu + t + NDMRG[ irrep_t ] * u
471               jump_tu += NDMRG[ irrep_t ] * NDMRG[ irrep_u ]
472    */
473 
474    size_D = new int[ num_irreps ];
475    for ( int irrep = 0; irrep < num_irreps; irrep++ ){
476       int jump_tu = 0;
477       for ( int irrep_t = 0; irrep_t < num_irreps; irrep_t++ ){
478          const int irrep_u = Irreps::directProd( irrep, irrep_t );
479          const int nact_t = indices->getNDMRG( irrep_t );
480          const int nact_u = indices->getNDMRG( irrep_u );
481          jump_tu += nact_t * nact_u;
482       }
483       size_D[ irrep ] = 2 * jump_tu;
484       int jump_ai = 0;
485       for ( int irrep_i = 0; irrep_i < num_irreps; irrep_i++ ){
486          const int irrep_a = Irreps::directProd( irrep_i, irrep );
487          const int nocc_i  = indices->getNOCC( irrep_i );
488          const int nvirt_a = indices->getNVIRT( irrep_a );
489          jump_ai += nocc_i * nvirt_a;
490       }
491       helper[ irrep + num_irreps * CHEMPS2_CASPT2_D ] = jump_ai * size_D[ irrep ];
492    }
493 
494    /*** Type B singlet : c_tiuj ( E_ti E_uj + E_tj E_ui ) / sqrt( 1 + delta_ij ) | 0 > with i <= j and t <= u
495                          c_tiuj = vector[ jump[ irrep + num_irreps * CHEMPS2_CASPT2_B_SINGLET ] + count_tu + size_B_singlet[ irrep ] * count_ij ]
496         Type B triplet : c_tiuj ( E_ti E_uj - E_tj E_ui ) / sqrt( 1 + delta_ij ) | 0 > with i <  j and t <  u
497                          c_tiuj = vector[ jump[ irrep + num_irreps * CHEMPS2_CASPT2_B_TRIPLET ] + count_tu + size_B_triplet[ irrep ] * count_ij ]
498         Type F singlet : c_atbu ( E_at E_bu + E_bt E_au ) / sqrt( 1 + delta_ab ) | 0 > with a <= b and t <= u
499                          c_atbu = vector[ jump[ irrep + num_irreps * CHEMPS2_CASPT2_F_SINGLET ] + count_tu + size_F_singlet[ irrep ] * count_ab ]
500         Type F triplet : c_atbu ( E_at E_bu - E_bt E_au ) / sqrt( 1 + delta_ab ) | 0 > with a <  b and t <  u
501                          c_atbu = vector[ jump[ irrep + num_irreps * CHEMPS2_CASPT2_F_TRIPLET ] + count_tu + size_F_triplet[ irrep ] * count_ab ]
502 
503         1/ jump_ab = 0
504            if ( irrep == 0 ):
505               irrep_ab = 0 .. num_irreps
506                  (SINGLET) ---> count_ab = jump_ab + a + b(b+1)/2
507                  (SINGLET) jump_ab += NVIRT[ irrep_ab ] * ( NVIRT[ irrep_ab ] + 1 ) / 2
508                  (TRIPLET) ---> count_ab = jump_ab + a + b(b-1)/2
509                  (TRIPLET) jump_ab += NVIRT[ irrep_ab ] * ( NVIRT[ irrep_ab ] - 1 ) / 2
510            else:
511               irrep_a = 0 .. num_irreps
512                  irrep_b = irrep x irrep_a
513                  if ( irrep_a < irrep_b ):
514                     ---> count_ab = jump_ab + a + NVIRT[ irrep_a ] * b
515                     jump_ab += NVIRT[ irrep_i ] * NVIRT[ irrep_j ]
516 
517         2/ jump_ij = 0
518            if irrep == 0:
519               irrep_ij = 0 .. num_irreps
520                  (SINGLET) ---> count_ij = jump_ij + i + j(j+1)/2
521                  (SINGLET) jump_ij += NOCC[ irrep_ij ] * ( NOCC[ irrep_ij ] + 1 ) / 2
522                  (TRIPLET) ---> count_ij = jump_ij + i + j(j-1)/2
523                  (TRIPLET) jump_ij += NOCC[ irrep_ij ] * ( NOCC[ irrep_ij ] - 1 ) / 2
524            else:
525               irrep_i = 0 .. num_irreps
526                  irrep_j = irrep x irrep_i
527                  if ( irrep_i < irrep_j ):
528                     ---> count_ij = jump_ij + i + NOCC[ irrep_i ] * j
529                     jump_ij += NOCC[ irrep_i ] * NOCC[ irrep_j ]
530 
531         3/ jump_tu = 0
532            if irrep == 0:
533               irrep_tu = 0 .. num_irreps
534                  (SINGLET) ---> count_tu = jump_tu + t + u(u+1)/2
535                  (SINGLET) jump_tu += NDMRG[ irrep_tu ] * ( NDMRG[ irrep_tu ] + 1 ) / 2
536                  (TRIPLET) ---> count_tu = jump_tu + t + u(u-1)/2
537                  (TRIPLET) jump_tu += NDMRG[ irrep_tu ] * ( NDMRG[ irrep_tu ] - 1 ) / 2
538            else:
539               irrep_t = 0 .. num_irreps
540                  irrep_u = irrep x irrep_t
541                  if ( irrep_t < irrep_u ):
542                     ---> count_tu = jump_tu + t + NDMRG[ irrep_t ] * u
543                     jump_tu += NDMRG[ irrep_t ] * NDMRG[ irrep_u ]
544    */
545 
546    size_B_singlet = new int[ num_irreps ];
547    size_B_triplet = new int[ num_irreps ];
548    size_F_singlet = new int[ num_irreps ];
549    size_F_triplet = new int[ num_irreps ];
550    for ( int irrep = 0; irrep < num_irreps; irrep++ ){
551       int jump_tu_singlet = 0;
552       int jump_tu_triplet = 0;
553       if ( irrep == 0 ){
554          for ( int irrep_tu = 0; irrep_tu < num_irreps; irrep_tu++ ){ // irrep_u == irrep_t
555             const int nact_tu = indices->getNDMRG( irrep_tu );
556             jump_tu_singlet += ( nact_tu * ( nact_tu + 1 )) / 2;
557             jump_tu_triplet += ( nact_tu * ( nact_tu - 1 )) / 2;
558          }
559       } else {
560          for ( int irrep_t = 0; irrep_t < num_irreps; irrep_t++ ){
561             const int irrep_u = Irreps::directProd( irrep, irrep_t );
562             if ( irrep_t < irrep_u ){
563                const int nact_t = indices->getNDMRG( irrep_t );
564                const int nact_u = indices->getNDMRG( irrep_u );
565                jump_tu_singlet += nact_t * nact_u;
566                jump_tu_triplet += nact_t * nact_u;
567             }
568          }
569       }
570       size_B_singlet[ irrep ] = jump_tu_singlet;
571       size_B_triplet[ irrep ] = jump_tu_triplet;
572       size_F_singlet[ irrep ] = jump_tu_singlet;
573       size_F_triplet[ irrep ] = jump_tu_triplet;
574 
575       int linsize_B_singlet = 0;
576       int linsize_B_triplet = 0;
577       int linsize_F_singlet = 0;
578       int linsize_F_triplet = 0;
579       if ( irrep == 0 ){ // irrep_i == irrep_j    or    irrep_a == irrep_b
580          for ( int irrep_ijab = 0; irrep_ijab < num_irreps; irrep_ijab++ ){
581             const int nocc_ij = indices->getNOCC( irrep_ijab );
582             linsize_B_singlet += ( nocc_ij * ( nocc_ij + 1 ))/2;
583             linsize_B_triplet += ( nocc_ij * ( nocc_ij - 1 ))/2;
584             const int nvirt_ab = indices->getNVIRT( irrep_ijab );
585             linsize_F_singlet += ( nvirt_ab * ( nvirt_ab + 1 ))/2;
586             linsize_F_triplet += ( nvirt_ab * ( nvirt_ab - 1 ))/2;
587          }
588       } else { // irrep_i < irrep_j = irrep_i x irrep    or    irrep_a < irrep_b = irrep_a x irrep
589          for ( int irrep_ai = 0; irrep_ai < num_irreps; irrep_ai++ ){
590             const int irrep_bj = Irreps::directProd( irrep, irrep_ai );
591             if ( irrep_ai < irrep_bj ){
592                const int nocc_i = indices->getNOCC( irrep_ai );
593                const int nocc_j = indices->getNOCC( irrep_bj );
594                linsize_B_singlet += nocc_i * nocc_j;
595                linsize_B_triplet += nocc_i * nocc_j;
596                const int nvirt_a = indices->getNVIRT( irrep_ai );
597                const int nvirt_b = indices->getNVIRT( irrep_bj );
598                linsize_F_singlet += nvirt_a * nvirt_b;
599                linsize_F_triplet += nvirt_a * nvirt_b;
600             }
601          }
602       }
603       helper[ irrep + num_irreps * CHEMPS2_CASPT2_B_SINGLET ] = linsize_B_singlet * size_B_singlet[ irrep ];
604       helper[ irrep + num_irreps * CHEMPS2_CASPT2_B_TRIPLET ] = linsize_B_triplet * size_B_triplet[ irrep ];
605       helper[ irrep + num_irreps * CHEMPS2_CASPT2_F_SINGLET ] = linsize_F_singlet * size_F_singlet[ irrep ];
606       helper[ irrep + num_irreps * CHEMPS2_CASPT2_F_TRIPLET ] = linsize_F_triplet * size_F_triplet[ irrep ];
607    }
608 
609    /*** Type E singlet : c_tiaj ( E_ti E_aj + E_tj E_ai ) / sqrt( 1 + delta_ij ) | 0 > with i <= j
610                          c_tiaj = vector[ jump[ irrep + num_irreps * CHEMPS2_CASPT2_E_SINGLET ] + count_t + size_E[ irrep ] * count_aij ]
611         Type E triplet : c_tiaj ( E_ti E_aj - E_tj E_ai ) / sqrt( 1 + delta_ij ) | 0 > with i <  j
612                          c_tiaj = vector[ jump[ irrep + num_irreps * CHEMPS2_CASPT2_E_TRIPLET ] + count_t + size_E[ irrep ] * count_aij ]
613 
614         1/ jump_aij = 0
615            irrep_a = 0 .. num_irreps
616               irrep_occ = irrep_a x irrep
617               if ( irrep_occ == 0 ):
618                  irrep_ij = 0 .. num_irreps
619                     (SINGLET) ---> count_aij = jump_aij + a + NVIRT[ irrep_a ] * ( i + j(j+1)/2 )
620                     (SINGLET) jump_aij += NVIRT[ irrep_a ] * NOCC[ irrep_ij ] * ( NOCC[ irrep_ij ] + 1 ) / 2
621                     (TRIPLET) ---> count_aij = jump_aij + a + NVIRT[ irrep_a ] * ( i + j(j-1)/2 )
622                     (TRIPLET) jump_aij += NVIRT[ irrep_a ] * NOCC[ irrep_ij ] * ( NOCC[ irrep_ij ] - 1 ) / 2
623               else:
624                  irrep_i = 0 .. num_irreps
625                     irrep_j = irrep_occ x irrep_i
626                        if ( irrep_i < irrep_j ):
627                           ---> count_aij = jump_aij + a + NVIRT[ irrep_a ] * ( i + NOCC[ irrep_i ] * j )
628                           jump_aij += NVIRT[ irrep_a ] * NOCC[ irrep_i ] * NOCC[ irrep_j ]
629         2/ count_t = 0 .. NDMRG[ irrep ]
630    */
631 
632    size_E = new int[ num_irreps ];
633    for ( int irrep = 0; irrep < num_irreps; irrep++ ){
634       size_E[ irrep ] = indices->getNDMRG( irrep );
635       int linsize_E_singlet = 0;
636       int linsize_E_triplet = 0;
637       for ( int irrep_a = 0; irrep_a < num_irreps; irrep_a++ ){
638          const int nvirt_a = indices->getNVIRT( irrep_a );
639          const int irrep_occ = Irreps::directProd( irrep, irrep_a );
640          if ( irrep_occ == 0 ){ // irrep_i == irrep_j
641             for ( int irrep_ij = 0; irrep_ij < num_irreps; irrep_ij++ ){
642                const int nocc_ij = indices->getNOCC( irrep_ij );
643                linsize_E_singlet += ( nvirt_a * nocc_ij * ( nocc_ij + 1 )) / 2;
644                linsize_E_triplet += ( nvirt_a * nocc_ij * ( nocc_ij - 1 )) / 2;
645             }
646          } else { // irrep_i < irrep_j = irrep_i x irrep_occ
647             for ( int irrep_i = 0; irrep_i < num_irreps; irrep_i++ ){
648                const int irrep_j = Irreps::directProd( irrep_occ, irrep_i );
649                if ( irrep_i < irrep_j ){
650                   const int nocc_i = indices->getNOCC( irrep_i );
651                   const int nocc_j = indices->getNOCC( irrep_j );
652                   linsize_E_singlet += nvirt_a * nocc_i * nocc_j;
653                   linsize_E_triplet += nvirt_a * nocc_i * nocc_j;
654                }
655             }
656          }
657       }
658       helper[ irrep + num_irreps * CHEMPS2_CASPT2_E_SINGLET ] = linsize_E_singlet * size_E[ irrep ];
659       helper[ irrep + num_irreps * CHEMPS2_CASPT2_E_TRIPLET ] = linsize_E_triplet * size_E[ irrep ];
660    }
661 
662    /*** Type G singlet : c_aibt ( E_ai E_bt + E_bi E_at ) / sqrt( 1 + delta_ab ) | 0 > with a <= b
663                          c_aibt = vector[ jump[ irrep + num_irreps * CHEMPS2_CASPT2_G_SINGLET ] + count_t + size_G[ irrep ] * count_abi ]
664         Type G triplet : c_aibt ( E_ai E_bt - E_bi E_at ) / sqrt( 1 + delta_ab ) | 0 > with a <  b
665                          c_aibt = vector[ jump[ irrep + num_irreps * CHEMPS2_CASPT2_G_TRIPLET ] + count_t + size_G[ irrep ] * count_abi ]
666 
667         1/ jump_abi = 0
668            irrep_i = 0 .. num_irreps
669               irrep_virt = irrep_i x irrep
670               if irrep_virt == 0:
671                  irrep_ab = 0 .. num_irreps
672                     (SINGLET) ---> count_abi = jump_abi + i + NOCC[ irrep_i ] * ( a + b(b+1)/2 )
673                     (SINGLET) jump_abi += NOCC[ irrep_i ] * NVIRT[ irrep_ab ] * ( NVIRT[ irrep_ab ] + 1 ) / 2
674                     (TRIPLET) ---> count_abi = jump_abi + i + NOCC[ irrep_i ] * ( a + b(b-1)/2 )
675                     (TRIPLET) jump_abi += NOCC[ irrep_i ] * NVIRT[ irrep_ab ] * ( NVIRT[ irrep_ab ] - 1 ) / 2
676               else:
677                  irrep_a = 0 .. num_irreps
678                     irrep_b = irrep_virt x irrep_a
679                        if ( irrep_a < irrep_b ):
680                           ---> count_abi = jump_abi + i + NOCC[ irrep_i ] * ( a + NVIRT[ irrep_a ] * b )
681                           jump_abi += NOCC[ irrep_i ] * NVIRT[ irrep_a ] * NVIRT[ irrep_b ]
682         2/ count_t = 0 .. NDMRG[ irrep ]
683    */
684 
685    size_G = new int[ num_irreps ];
686    for ( int irrep = 0; irrep < num_irreps; irrep++ ){
687       size_G[ irrep ] = indices->getNDMRG( irrep );
688       int linsize_G_singlet = 0;
689       int linsize_G_triplet = 0;
690       for ( int irrep_i = 0; irrep_i < num_irreps; irrep_i++ ){
691          const int nocc_i = indices->getNOCC( irrep_i );
692          const int irrep_virt = Irreps::directProd( irrep, irrep_i );
693          if ( irrep_virt == 0 ){ // irrep_a == irrep_b
694             for ( int irrep_ab = 0; irrep_ab < num_irreps; irrep_ab++ ){
695                const int nvirt_ab = indices->getNVIRT( irrep_ab );
696                linsize_G_singlet += ( nocc_i * nvirt_ab * ( nvirt_ab + 1 )) / 2;
697                linsize_G_triplet += ( nocc_i * nvirt_ab * ( nvirt_ab - 1 )) / 2;
698             }
699          } else { // irrep_a < irrep_b = irrep_a x irrep_virt
700             for ( int irrep_a = 0; irrep_a < num_irreps; irrep_a++ ){
701                const int irrep_b = Irreps::directProd( irrep_virt, irrep_a );
702                if ( irrep_a < irrep_b ){
703                   const int nvirt_a = indices->getNVIRT( irrep_a );
704                   const int nvirt_b = indices->getNVIRT( irrep_b );
705                   linsize_G_singlet += nocc_i * nvirt_a * nvirt_b;
706                   linsize_G_triplet += nocc_i * nvirt_a * nvirt_b;
707                }
708             }
709          }
710       }
711       helper[ irrep + num_irreps * CHEMPS2_CASPT2_G_SINGLET ] = linsize_G_singlet * size_G[ irrep ];
712       helper[ irrep + num_irreps * CHEMPS2_CASPT2_G_TRIPLET ] = linsize_G_triplet * size_G[ irrep ];
713    }
714 
715    /*** Type H singlet : c_aibj ( E_ai E_bj + E_bi E_aj ) / sqrt( ( 1 + delta_ij ) * ( 1 + delta_ab ) ) | 0 > with a <= b and i <= j
716                          c_aibj = vector[ jump[ irrep + num_irreps * CHEMPS2_CASPT2_H_SINGLET ] + count_aibj ]
717         Type H triplet : c_aibj ( E_ai E_bj - E_bi E_aj ) / sqrt( ( 1 + delta_ij ) * ( 1 + delta_ab ) ) | 0 > with a <  b and i <  j
718                          c_aibj = vector[ jump[ irrep + num_irreps * CHEMPS2_CASPT2_H_TRIPLET ] + count_aibj ]
719 
720         1/ irrep = 0 .. num_irreps
721               if ( irrep == 0 ):
722                  jump_aibj = 0
723                  irrep_ij = 0 .. num_irreps
724                     (SINGLET) linsize_ij = NOCC[ irrep_ij ] * ( NOCC[ irrep_ij ] + 1 ) / 2
725                     (TRIPLET) linsize_ij = NOCC[ irrep_ij ] * ( NOCC[ irrep_ij ] - 1 ) / 2
726                     irrep_ab = 0 .. num_irreps
727                        (SINGLET) linsize_ab = NVIRT[ irrep_ab ] * ( NVIRT[ irrep_ab ] + 1 ) / 2
728                        (TRIPLET) linsize_ab = NVIRT[ irrep_ab ] * ( NVIRT[ irrep_ab ] - 1 ) / 2
729                        (SINGLET) ---> count_aibj = jump_aibj + i + j(j+1)/2 + linsize_ij * ( a + b(b+1)/2 )
730                        (TRIPLET) ---> count_aibj = jump_aibj + i + j(j-1)/2 + linsize_ij * ( a + b(b-1)/2 )
731                        jump_aibj += linsize_ij * linsize_ab
732               else:
733                  jump_aibj = 0
734                  irrep_i = 0 .. num_irreps
735                     irrep_j = irrep x irrep_i
736                     if ( irrep_i < irrep_j ):
737                        irrep_a = 0 .. num_irreps
738                           irrep_b = irrep x irrep_a
739                           if ( irrep_a < irrep_b ):
740                              ---> count_aibj = jump_aibj + i + NOCC[ irrep_i ] * ( j + NOCC[ irrep_j ] * ( a + NVIRT[ irrep_a ] * b ) )
741                              jump_aibj += NOCC[ irrep_i ] * NOCC[ irrep_j ] * NVIRT[ irrep_a ] * NVIRT[ irrep_b ]
742    */
743 
744    for ( int irrep = 0; irrep < num_irreps; irrep++ ){
745       int linsize_H_singlet = 0;
746       int linsize_H_triplet = 0;
747       if ( irrep == 0 ){ // irrep_i == irrep_j  and  irrep_a == irrep_b
748          int linsize_ij_singlet = 0;
749          int linsize_ij_triplet = 0;
750          for ( int irrep_ij = 0; irrep_ij < num_irreps; irrep_ij++ ){
751             const int nocc_ij = indices->getNOCC( irrep_ij );
752             linsize_ij_singlet += ( nocc_ij * ( nocc_ij + 1 )) / 2;
753             linsize_ij_triplet += ( nocc_ij * ( nocc_ij - 1 )) / 2;
754          }
755          int linsize_ab_singlet = 0;
756          int linsize_ab_triplet = 0;
757          for ( int irrep_ab = 0; irrep_ab < num_irreps; irrep_ab++ ){
758             const int nvirt_ab = indices->getNVIRT( irrep_ab );
759             linsize_ab_singlet += ( nvirt_ab * ( nvirt_ab + 1 )) / 2;
760             linsize_ab_triplet += ( nvirt_ab * ( nvirt_ab - 1 )) / 2;
761          }
762          linsize_H_singlet = linsize_ij_singlet * linsize_ab_singlet;
763          linsize_H_triplet = linsize_ij_triplet * linsize_ab_triplet;
764       } else { // irrep_i < irrep_j = irrep_i x irrep   and   irrep_a < irrep_b = irrep_a x irrep
765          int linsize_ij = 0;
766          for ( int irrep_i = 0; irrep_i < num_irreps; irrep_i++ ){
767             const int irrep_j = Irreps::directProd( irrep, irrep_i );
768             if ( irrep_i < irrep_j ){ linsize_ij += indices->getNOCC( irrep_i ) * indices->getNOCC( irrep_j ); }
769          }
770          int linsize_ab = 0;
771          for ( int irrep_a = 0; irrep_a < num_irreps; irrep_a++ ){
772             const int irrep_b = Irreps::directProd( irrep, irrep_a );
773             if ( irrep_a < irrep_b ){ linsize_ab += indices->getNVIRT( irrep_a ) * indices->getNVIRT( irrep_b ); }
774          }
775          linsize_H_singlet = linsize_ij * linsize_ab;
776          linsize_H_triplet = linsize_ij * linsize_ab;
777       }
778       helper[ irrep + num_irreps * CHEMPS2_CASPT2_H_SINGLET ] = linsize_H_singlet;
779       helper[ irrep + num_irreps * CHEMPS2_CASPT2_H_TRIPLET ] = linsize_H_triplet;
780    }
781 
782    jump = new int[ CHEMPS2_CASPT2_NUM_CASES * num_irreps + 1 ];
783    jump[ 0 ] = 0;
784    for ( int cnt = 0; cnt < CHEMPS2_CASPT2_NUM_CASES * num_irreps; cnt++ ){ jump[ cnt+1 ] = jump[ cnt ] + helper[ cnt ]; }
785    delete [] helper;
786    const int total_size = jump[ CHEMPS2_CASPT2_NUM_CASES * num_irreps ];
787    assert( total_size == vector_length( indices ) );
788    cout << "CASPT2 : Old size V_SD space  = " << total_size << endl;
789    return total_size;
790 
791 }
792 
vector_length(const DMRGSCFindices * idx)793 long long CheMPS2::CASPT2::vector_length( const DMRGSCFindices * idx ){
794 
795    const int NUM_IRREPS = idx->getNirreps();
796    long long length = 0;
797    for ( int i1 = 0; i1 < NUM_IRREPS; i1++ ){
798       const long long nocc1 = idx->getNOCC( i1 );
799       const long long nact1 = idx->getNDMRG( i1 );
800       const long long nvir1 = idx->getNVIRT( i1 );
801       for ( int i2 = 0; i2 < NUM_IRREPS; i2++ ){
802          const long long nocc2 = idx->getNOCC( i2 );
803          const long long nact2 = idx->getNDMRG( i2 );
804          for ( int i3 = 0; i3 < NUM_IRREPS; i3++ ){
805             const int i4 = Irreps::directProd( Irreps::directProd( i1, i2 ), i3 );
806             const long long nact3 = idx->getNDMRG( i3 );
807             const long long nvir3 = idx->getNVIRT( i3 );
808             const long long nocc4 = idx->getNOCC( i4 );
809             const long long nact4 = idx->getNDMRG( i4 );
810             const long long nvir4 = idx->getNVIRT( i4 );
811 
812             length +=     nocc1 * nact2 * nact3 * nact4; // A:  E_ti E_uv | 0 >
813             length +=     nact1 * nact2 * nact3 * nvir4; // C:  E_at E_uv | 0 >
814             length += 2 * nocc1 * nact2 * nact3 * nvir4; // D:  E_ai E_tu | 0 >  and  E_ti E_au | 0 >
815             length +=     nocc1 * nocc2 * nact3 * nvir4; // E:  E_ti E_aj | 0 >
816             length +=     nocc1 * nact2 * nvir3 * nvir4; // G:  E_ai E_bt | 0 >
817             if ( i2 < i4 ){ // 2 < 4 and irrep_2 < irrep_4
818                length += nact1 * nact3 * nocc2 * nocc4;  // B:  E_ti E_uj | 0 >
819                length += nvir1 * nvir3 * nact2 * nact4;  // F:  E_at E_bu | 0 >
820                length += nvir1 * nvir3 * nocc2 * nocc4;  // H:  E_ai E_bj | 0 >
821             }
822             if ( i2 == i4 ){ // i2 == i4 implies i1 == i3
823                // 2 < 4 and irrep_2 == irrep_4
824                length += ( nact1 * nact3 * nocc2 * ( nocc2 - 1 ) ) / 2; // B:  E_ti E_uj | 0 >
825                length += ( nvir1 * nvir3 * nact2 * ( nact2 - 1 ) ) / 2; // F:  E_at E_bu | 0 >
826                length += ( nvir1 * nvir3 * nocc2 * ( nocc2 - 1 ) ) / 2; // H:  E_ai E_bj | 0 >
827                // 2 == 4 and 1 <= 3
828                length += ( nact1 * ( nact3 + 1 ) * nocc2 ) / 2; // B:  E_ti E_uj | 0 >
829                length += ( nvir1 * ( nvir3 + 1 ) * nact2 ) / 2; // F:  E_at E_bu | 0 >
830                length += ( nvir1 * ( nvir3 + 1 ) * nocc2 ) / 2; // H:  E_ai E_bj | 0 >
831             }
832          }
833       }
834    }
835 
836    return length;
837 
838 }
839 
recreatehelper1(double * FOCK,double * OVLP,int SIZE,double * work,double * eigs,int lwork)840 int CheMPS2::CASPT2::recreatehelper1( double * FOCK, double * OVLP, int SIZE, double * work, double * eigs, int lwork ){
841 
842    if ( SIZE == 0 ){ return SIZE; }
843 
844    // S = U_S eigs_S U_S^T
845    char jobz = 'V';
846    char uplo = 'U';
847    int info;
848    dsyev_( &jobz, &uplo, &SIZE, OVLP, &SIZE, eigs, work, &lwork, &info ); // eigs in ascending order
849 
850    // Discard smallest eigenvalues
851    int skip = 0;
852    bool ctu = true;
853    while (( skip < SIZE ) && ( ctu )){
854       if ( eigs[ skip ] < CheMPS2::CASPT2_OVLP_CUTOFF ){
855          skip++;
856       } else {
857          ctu = false;
858       }
859    }
860    int NEWSIZE = SIZE - skip;
861 
862    if ( NEWSIZE == 0 ){ return NEWSIZE; }
863 
864    // OVLP  <---  U_S eigs_S^{-0.5}
865    for ( int col = skip; col < SIZE; col++ ){
866       const double prefactor = 1.0 / sqrt( eigs[ col ] );
867       for ( int row = 0; row < SIZE; row++ ){
868          OVLP[ row + SIZE * col ] *= prefactor;
869       }
870    }
871 
872    // FOCK  <---  FOCK_tilde = eigs_S^{-0.5} U_S^T FOCK U_S eigs_S^{-0.5}
873    char notrans = 'N';
874    char trans   = 'T';
875    double one   = 1.0;
876    double set   = 0.0;
877    dgemm_( &notrans, &notrans, &SIZE,    &NEWSIZE, &SIZE, &one, FOCK, &SIZE, OVLP + skip * SIZE, &SIZE, &set, work, &SIZE    ); // work = FOCK * V * eigs^{-0.5}
878    dgemm_( &trans,   &notrans, &NEWSIZE, &NEWSIZE, &SIZE, &one, OVLP + skip * SIZE, &SIZE, work, &SIZE, &set, FOCK, &NEWSIZE ); // FOCK = ( V * eigs^{-0.5} )^T * work
879 
880    // FOCK_tilde = U_F_tilde eigs_F_tilde U_F_tilde^T
881    dsyev_( &jobz, &uplo, &NEWSIZE, FOCK, &NEWSIZE, eigs, work, &lwork, &info ); // eigs in ascending order
882 
883    // OVLP  <---  U_S eigs_S^{-0.5} U_F_tilde
884    dgemm_( &notrans, &notrans, &SIZE, &NEWSIZE, &NEWSIZE, &one, OVLP + skip * SIZE, &SIZE, FOCK, &NEWSIZE, &set, work, &SIZE );
885    int size_copy = SIZE * NEWSIZE;
886    int inc1 = 1;
887    dcopy_( &size_copy, work, &inc1, OVLP, &inc1 );
888 
889    // FOCK  <---  eigs_F_tilde
890    dcopy_( &NEWSIZE, eigs, &inc1, FOCK, &inc1 );
891 
892    return NEWSIZE;
893 
894 }
895 
recreatehelper2(double * LEFT,double * RIGHT,double ** matrix,double * work,int OLD_LEFT,int NEW_LEFT,int OLD_RIGHT,int NEW_RIGHT,const int number)896 void CheMPS2::CASPT2::recreatehelper2( double * LEFT, double * RIGHT, double ** matrix, double * work, int OLD_LEFT, int NEW_LEFT, int OLD_RIGHT, int NEW_RIGHT, const int number ){
897 
898    double set = 0.0;
899    double one = 1.0;
900    char trans   = 'T';
901    char notrans = 'N';
902    for ( int count = 0; count < number; count++ ){
903       // work  <---  LEFT[ old_left, new_left ]^T * matrix[ count ][ old_left, old_right ]
904       dgemm_( &trans, &notrans, &NEW_LEFT, &OLD_RIGHT, &OLD_LEFT, &one, LEFT, &OLD_LEFT, matrix[ count ], &OLD_LEFT, &set, work, &NEW_LEFT );
905       // matrix[ count ]  <---  work[ new_left, old_right ] * RIGHT[ old_right, new_right ]
906       dgemm_( &notrans, &notrans, &NEW_LEFT, &NEW_RIGHT, &OLD_RIGHT, &one, work, &NEW_LEFT, RIGHT, &OLD_RIGHT, &set, matrix[ count ], &NEW_LEFT );
907    }
908 
909 }
910 
recreatehelper3(double * OVLP,int OLDSIZE,int NEWSIZE,double * rhs_old,double * rhs_new,const int num_rhs)911 void CheMPS2::CASPT2::recreatehelper3( double * OVLP, int OLDSIZE, int NEWSIZE, double * rhs_old, double * rhs_new, const int num_rhs ){
912 
913    // rhs <-- eigs^{-0.5} V^T rhs
914    int inc1 = 1;
915    double set = 0.0;
916    double one = 1.0;
917    char trans = 'T';
918    for ( int sector = 0; sector < num_rhs; sector++ ){
919       dgemv_( &trans, &OLDSIZE, &NEWSIZE, &one, OVLP, &OLDSIZE, rhs_old + OLDSIZE * sector, &inc1, &set, rhs_new + NEWSIZE * sector, &inc1 );
920    }
921 
922 }
923 
recreate()924 void CheMPS2::CASPT2::recreate(){
925 
926    const int maxsize = get_maxsize(); // Minimum 3
927    const int lwork = maxsize * maxsize;
928    double * work  = new double[ lwork ];
929    double * eigs  = new double[ maxsize ];
930 
931    int * newsize_A = new int[ num_irreps ];
932    int * newsize_C = new int[ num_irreps ];
933    int * newsize_D = new int[ num_irreps ];
934    int * newsize_E = new int[ num_irreps ];
935    int * newsize_G = new int[ num_irreps ];
936    int * newsize_B_singlet = new int[ num_irreps ];
937    int * newsize_B_triplet = new int[ num_irreps ];
938    int * newsize_F_singlet = new int[ num_irreps ];
939    int * newsize_F_triplet = new int[ num_irreps ];
940 
941    for ( int irrep = 0; irrep < num_irreps; irrep++ ){
942       newsize_A[ irrep ] = recreatehelper1( FAA[ irrep ], SAA[ irrep ], size_A[ irrep ], work, eigs, lwork );
943       newsize_C[ irrep ] = recreatehelper1( FCC[ irrep ], SCC[ irrep ], size_C[ irrep ], work, eigs, lwork );
944       newsize_D[ irrep ] = recreatehelper1( FDD[ irrep ], SDD[ irrep ], size_D[ irrep ], work, eigs, lwork );
945       newsize_E[ irrep ] = recreatehelper1( FEE[ irrep ], SEE[ irrep ], size_E[ irrep ], work, eigs, lwork );
946       newsize_G[ irrep ] = recreatehelper1( FGG[ irrep ], SGG[ irrep ], size_G[ irrep ], work, eigs, lwork );
947       newsize_B_singlet[ irrep ] = recreatehelper1( FBB_singlet[ irrep ], SBB_singlet[ irrep ], size_B_singlet[ irrep ], work, eigs, lwork );
948       newsize_B_triplet[ irrep ] = recreatehelper1( FBB_triplet[ irrep ], SBB_triplet[ irrep ], size_B_triplet[ irrep ], work, eigs, lwork );
949       newsize_F_singlet[ irrep ] = recreatehelper1( FFF_singlet[ irrep ], SFF_singlet[ irrep ], size_F_singlet[ irrep ], work, eigs, lwork );
950       newsize_F_triplet[ irrep ] = recreatehelper1( FFF_triplet[ irrep ], SFF_triplet[ irrep ], size_F_triplet[ irrep ], work, eigs, lwork );
951    }
952 
953    for ( int IL = 0; IL < num_irreps; IL++ ){
954       for ( int IR = 0; IR < num_irreps; IR++ ){
955 
956          const int irrep_w = Irreps::directProd( IL, IR );
957          const int num_w   = indices->getNDMRG( irrep_w );
958 
959          if ( newsize_A[ IL ] * newsize_D[ IR ] * num_w > 0 ){
960             recreatehelper2( SAA[ IL ], SDD[ IR ], FAD[ IL ][ IR ], work, size_A[ IL ], newsize_A[ IL ], size_D[ IR ], newsize_D[ IR ], num_w );
961          }
962 
963          if ( newsize_C[ IL ] * newsize_D[ IR ] * num_w > 0 ){
964             recreatehelper2( SCC[ IL ], SDD[ IR ], FCD[ IL ][ IR ], work, size_C[ IL ], newsize_C[ IL ], size_D[ IR ], newsize_D[ IR ], num_w );
965          }
966 
967          if ( newsize_A[ IL ] * newsize_B_singlet[ IR ] * num_w > 0 ){
968             recreatehelper2( SAA[ IL ], SBB_singlet[ IR ], FAB_singlet[ IL ][ IR ], work, size_A[ IL ], newsize_A[ IL ], size_B_singlet[ IR ], newsize_B_singlet[ IR ], num_w );
969          }
970 
971          if ( newsize_A[ IL ] * newsize_B_triplet[ IR ] * num_w > 0 ){
972             recreatehelper2( SAA[ IL ], SBB_triplet[ IR ], FAB_triplet[ IL ][ IR ], work, size_A[ IL ], newsize_A[ IL ], size_B_triplet[ IR ], newsize_B_triplet[ IR ], num_w );
973          }
974 
975          if ( newsize_C[ IL ] * newsize_F_singlet[ IR ] * num_w > 0 ){
976             recreatehelper2( SCC[ IL ], SFF_singlet[ IR ], FCF_singlet[ IL ][ IR ], work, size_C[ IL ], newsize_C[ IL ], size_F_singlet[ IR ], newsize_F_singlet[ IR ], num_w );
977          }
978 
979          if ( newsize_C[ IL ] * newsize_F_triplet[ IR ] * num_w > 0 ){
980             recreatehelper2( SCC[ IL ], SFF_triplet[ IR ], FCF_triplet[ IL ][ IR ], work, size_C[ IL ], newsize_C[ IL ], size_F_triplet[ IR ], newsize_F_triplet[ IR ], num_w );
981          }
982 
983          if ( newsize_B_singlet[ IL ] * newsize_E[ IR ] * num_w > 0 ){
984             recreatehelper2( SBB_singlet[ IL ], SEE[ IR ], FBE_singlet[ IL ][ IR ], work, size_B_singlet[ IL ], newsize_B_singlet[ IL ], size_E[ IR ], newsize_E[ IR ], num_w );
985          }
986 
987          if ( newsize_B_triplet[ IL ] * newsize_E[ IR ] * num_w > 0 ){
988             recreatehelper2( SBB_triplet[ IL ], SEE[ IR ], FBE_triplet[ IL ][ IR ], work, size_B_triplet[ IL ], newsize_B_triplet[ IL ], size_E[ IR ], newsize_E[ IR ], num_w );
989          }
990 
991          if ( newsize_F_singlet[ IL ] * newsize_G[ IR ] * num_w > 0 ){
992             recreatehelper2( SFF_singlet[ IL ], SGG[ IR ], FFG_singlet[ IL ][ IR ], work, size_F_singlet[ IL ], newsize_F_singlet[ IL ], size_G[ IR ], newsize_G[ IR ], num_w );
993          }
994 
995          if ( newsize_F_triplet[ IL ] * newsize_G[ IR ] * num_w > 0 ){
996             recreatehelper2( SFF_triplet[ IL ], SGG[ IR ], FFG_triplet[ IL ][ IR ], work, size_F_triplet[ IL ], newsize_F_triplet[ IL ], size_G[ IR ], newsize_G[ IR ], num_w );
997          }
998 
999          if ( newsize_D[ IL ] * newsize_E[ IR ] * num_w > 0 ){
1000             recreatehelper2( SDD[ IL ], SEE[ IR ], FDE_singlet[ IL ][ IR ], work, size_D[ IL ], newsize_D[ IL ], size_E[ IR ], newsize_E[ IR ], num_w );
1001             recreatehelper2( SDD[ IL ], SEE[ IR ], FDE_triplet[ IL ][ IR ], work, size_D[ IL ], newsize_D[ IL ], size_E[ IR ], newsize_E[ IR ], num_w );
1002          }
1003 
1004          if ( newsize_D[ IL ] * newsize_G[ IR ] * num_w > 0 ){
1005             recreatehelper2( SDD[ IL ], SGG[ IR ], FDG_singlet[ IL ][ IR ], work, size_D[ IL ], newsize_D[ IL ], size_G[ IR ], newsize_G[ IR ], num_w );
1006             recreatehelper2( SDD[ IL ], SGG[ IR ], FDG_triplet[ IL ][ IR ], work, size_D[ IL ], newsize_D[ IL ], size_G[ IR ], newsize_G[ IR ], num_w );
1007          }
1008 
1009          if ( IR == 0 ){ // IL == irrep_w
1010             if ( newsize_E[ IL ] * num_w > 0 ){
1011                double one = 1.0;
1012                recreatehelper2( SEE[ IL ], &one, FEH[ IL ], work, size_E[ IL ], newsize_E[ IL ], 1, 1, num_w );
1013             }
1014 
1015             if ( newsize_G[ IL ] * num_w > 0 ){
1016                double one = 1.0;
1017                recreatehelper2( SGG[ IL ], &one, FGH[ IL ], work, size_G[ IL ], newsize_G[ IL ], 1, 1, num_w );
1018             }
1019          }
1020       }
1021    }
1022 
1023    delete [] work;
1024    delete [] eigs;
1025 
1026    double * tempvector_rhs = new double[ jump[ num_irreps * CHEMPS2_CASPT2_NUM_CASES ] ];
1027    int * helper = new int[ num_irreps * CHEMPS2_CASPT2_NUM_CASES ];
1028    for ( int ptr = 0; ptr < num_irreps * CHEMPS2_CASPT2_NUM_CASES; ptr++ ){ helper[ ptr ] = 0; }
1029 
1030    for ( int irrep = 0; irrep < num_irreps; irrep++ ){
1031 
1032       if ( newsize_A[ irrep ] > 0 ){
1033          const int ptr = irrep + num_irreps * CHEMPS2_CASPT2_A;
1034          const int num_rhs = ( jump[ ptr + 1 ] - jump[ ptr ] ) / size_A[ irrep ];
1035          assert( num_rhs * size_A[ irrep ] == jump[ ptr + 1 ] - jump[ ptr ] );
1036          helper[ ptr ] = num_rhs * newsize_A[ irrep ];
1037          recreatehelper3( SAA[ irrep ], size_A[ irrep ], newsize_A[ irrep ], vector_rhs + jump[ ptr ], tempvector_rhs + jump[ ptr ], num_rhs );
1038       }
1039 
1040       if ( newsize_B_singlet[ irrep ] > 0 ){
1041          const int ptr = irrep + num_irreps * CHEMPS2_CASPT2_B_SINGLET;
1042          const int num_rhs = ( jump[ ptr + 1 ] - jump[ ptr ] ) / size_B_singlet[ irrep ];
1043          assert( num_rhs * size_B_singlet[ irrep ] == jump[ ptr + 1 ] - jump[ ptr ] );
1044          helper[ ptr ] = num_rhs * newsize_B_singlet[ irrep ];
1045          recreatehelper3( SBB_singlet[ irrep ], size_B_singlet[ irrep ], newsize_B_singlet[ irrep ], vector_rhs + jump[ ptr ], tempvector_rhs + jump[ ptr ], num_rhs );
1046       }
1047 
1048       if ( newsize_B_triplet[ irrep ] > 0 ){
1049          const int ptr = irrep + num_irreps * CHEMPS2_CASPT2_B_TRIPLET;
1050          const int num_rhs = ( jump[ ptr + 1 ] - jump[ ptr ] ) / size_B_triplet[ irrep ];
1051          assert( num_rhs * size_B_triplet[ irrep ] == jump[ ptr + 1 ] - jump[ ptr ] );
1052          helper[ ptr ] = num_rhs * newsize_B_triplet[ irrep ];
1053          recreatehelper3( SBB_triplet[ irrep ], size_B_triplet[ irrep ], newsize_B_triplet[ irrep ], vector_rhs + jump[ ptr ], tempvector_rhs + jump[ ptr ], num_rhs );
1054       }
1055 
1056       if ( newsize_C[ irrep ] > 0 ){
1057          const int ptr = irrep + num_irreps * CHEMPS2_CASPT2_C;
1058          const int num_rhs = ( jump[ ptr + 1 ] - jump[ ptr ] ) / size_C[ irrep ];
1059          assert( num_rhs * size_C[ irrep ] == jump[ ptr + 1 ] - jump[ ptr ] );
1060          helper[ ptr ] = num_rhs * newsize_C[ irrep ];
1061          recreatehelper3( SCC[ irrep ], size_C[ irrep ], newsize_C[ irrep ], vector_rhs + jump[ ptr ], tempvector_rhs + jump[ ptr ], num_rhs );
1062       }
1063 
1064       if ( newsize_D[ irrep ] > 0 ){
1065          const int ptr = irrep + num_irreps * CHEMPS2_CASPT2_D;
1066          const int num_rhs = ( jump[ ptr + 1 ] - jump[ ptr ] ) / size_D[ irrep ];
1067          assert( num_rhs * size_D[ irrep ] == jump[ ptr + 1 ] - jump[ ptr ] );
1068          helper[ ptr ] = num_rhs * newsize_D[ irrep ];
1069          recreatehelper3( SDD[ irrep ], size_D[ irrep ], newsize_D[ irrep ], vector_rhs + jump[ ptr ], tempvector_rhs + jump[ ptr ], num_rhs );
1070       }
1071 
1072       if ( newsize_E[ irrep ] > 0 ){
1073          const int ptr1 = irrep + num_irreps * CHEMPS2_CASPT2_E_SINGLET;
1074          const int num_rhs1 = ( jump[ ptr1 + 1 ] - jump[ ptr1 ] ) / size_E[ irrep ];
1075          assert( num_rhs1 * size_E[ irrep ] == jump[ ptr1 + 1 ] - jump[ ptr1 ] );
1076          helper[ ptr1 ] = num_rhs1 * newsize_E[ irrep ];
1077          recreatehelper3( SEE[ irrep ], size_E[ irrep ], newsize_E[ irrep ], vector_rhs + jump[ ptr1 ], tempvector_rhs + jump[ ptr1 ], num_rhs1 );
1078          const int ptr2 = irrep + num_irreps * CHEMPS2_CASPT2_E_TRIPLET;
1079          const int num_rhs2 = ( jump[ ptr2 + 1 ] - jump[ ptr2 ] ) / size_E[ irrep ];
1080          assert( num_rhs2 * size_E[ irrep ] == jump[ ptr2 + 1 ] - jump[ ptr2 ] );
1081          helper[ ptr2 ] = num_rhs2 * newsize_E[ irrep ];
1082          recreatehelper3( SEE[ irrep ], size_E[ irrep ], newsize_E[ irrep ], vector_rhs + jump[ ptr2 ], tempvector_rhs + jump[ ptr2 ], num_rhs2 );
1083       }
1084 
1085       if ( newsize_F_singlet[ irrep ] > 0 ){
1086          const int ptr = irrep + num_irreps * CHEMPS2_CASPT2_F_SINGLET;
1087          const int num_rhs = ( jump[ ptr + 1 ] - jump[ ptr ] ) / size_F_singlet[ irrep ];
1088          assert( num_rhs * size_F_singlet[ irrep ] == jump[ ptr + 1 ] - jump[ ptr ] );
1089          helper[ ptr ] = num_rhs * newsize_F_singlet[ irrep ];
1090          recreatehelper3( SFF_singlet[ irrep ], size_F_singlet[ irrep ], newsize_F_singlet[ irrep ], vector_rhs + jump[ ptr ], tempvector_rhs + jump[ ptr ], num_rhs );
1091       }
1092 
1093       if ( newsize_F_triplet[ irrep ] > 0 ){
1094          const int ptr = irrep + num_irreps * CHEMPS2_CASPT2_F_TRIPLET;
1095          const int num_rhs = ( jump[ ptr + 1 ] - jump[ ptr ] ) / size_F_triplet[ irrep ];
1096          assert( num_rhs * size_F_triplet[ irrep ] == jump[ ptr + 1 ] - jump[ ptr ] );
1097          helper[ ptr ] = num_rhs * newsize_F_triplet[ irrep ];
1098          recreatehelper3( SFF_triplet[ irrep ], size_F_triplet[ irrep ], newsize_F_triplet[ irrep ], vector_rhs + jump[ ptr ], tempvector_rhs + jump[ ptr ], num_rhs );
1099       }
1100 
1101       if ( newsize_G[ irrep ] > 0 ){
1102          const int ptr1 = irrep + num_irreps * CHEMPS2_CASPT2_G_SINGLET;
1103          const int num_rhs1 = ( jump[ ptr1 + 1 ] - jump[ ptr1 ] ) / size_G[ irrep ];
1104          assert( num_rhs1 * size_G[ irrep ] == jump[ ptr1 + 1 ] - jump[ ptr1 ] );
1105          helper[ ptr1 ] = num_rhs1 * newsize_G[ irrep ];
1106          recreatehelper3( SGG[ irrep ], size_G[ irrep ], newsize_G[ irrep ], vector_rhs + jump[ ptr1 ], tempvector_rhs + jump[ ptr1 ], num_rhs1 );
1107          const int ptr2 = irrep + num_irreps * CHEMPS2_CASPT2_G_TRIPLET;
1108          const int num_rhs2 = ( jump[ ptr2 + 1 ] - jump[ ptr2 ] ) / size_G[ irrep ];
1109          assert( num_rhs2 * size_G[ irrep ] == jump[ ptr2 + 1 ] - jump[ ptr2 ] );
1110          helper[ ptr2 ] = num_rhs2 * newsize_G[ irrep ];
1111          recreatehelper3( SGG[ irrep ], size_G[ irrep ], newsize_G[ irrep ], vector_rhs + jump[ ptr2 ], tempvector_rhs + jump[ ptr2 ], num_rhs2 );
1112       }
1113 
1114       const int ptr1 = irrep + num_irreps * CHEMPS2_CASPT2_H_SINGLET;
1115       helper[ ptr1 ] = jump[ ptr1 + 1 ] - jump[ ptr1 ];
1116       int inc1 = 1;
1117       dcopy_( helper + ptr1, vector_rhs + jump[ ptr1 ], &inc1, tempvector_rhs + jump[ ptr1 ], &inc1 );
1118       const int ptr2 = irrep + num_irreps * CHEMPS2_CASPT2_H_TRIPLET;
1119       helper[ ptr2 ] = jump[ ptr2 + 1 ] - jump[ ptr2 ];
1120       dcopy_( helper + ptr2, vector_rhs + jump[ ptr2 ], &inc1, tempvector_rhs + jump[ ptr2 ], &inc1 );
1121 
1122    }
1123 
1124    delete [] vector_rhs;
1125    delete [] size_A;         size_A = newsize_A;
1126    delete [] size_C;         size_C = newsize_C;
1127    delete [] size_D;         size_D = newsize_D;
1128    delete [] size_E;         size_E = newsize_E;
1129    delete [] size_G;         size_G = newsize_G;
1130    delete [] size_B_singlet; size_B_singlet = newsize_B_singlet;
1131    delete [] size_B_triplet; size_B_triplet = newsize_B_triplet;
1132    delete [] size_F_singlet; size_F_singlet = newsize_F_singlet;
1133    delete [] size_F_triplet; size_F_triplet = newsize_F_triplet;
1134 
1135    int * newjump = new int[ num_irreps * CHEMPS2_CASPT2_NUM_CASES + 1 ];
1136    newjump[ 0 ] = 0;
1137    for ( int cnt = 0; cnt < num_irreps * CHEMPS2_CASPT2_NUM_CASES; cnt++ ){
1138       newjump[ cnt + 1 ] = newjump[ cnt ] + helper[ cnt ];
1139    }
1140    vector_rhs = new double[ newjump[ num_irreps * CHEMPS2_CASPT2_NUM_CASES ] ];
1141    int inc1 = 1;
1142    for ( int cnt = 0; cnt < num_irreps * CHEMPS2_CASPT2_NUM_CASES; cnt++ ){
1143       dcopy_( helper + cnt, tempvector_rhs + jump[ cnt ], &inc1, vector_rhs + newjump[ cnt ], &inc1 );
1144    }
1145    delete [] tempvector_rhs;
1146    delete [] helper;
1147    delete [] jump;
1148    jump = newjump;
1149 
1150    for ( int irrep = 0; irrep < num_irreps; irrep++ ){
1151       delete [] SAA[ irrep ];
1152       delete [] SCC[ irrep ];
1153       delete [] SDD[ irrep ];
1154       delete [] SEE[ irrep ];
1155       delete [] SGG[ irrep ];
1156       delete [] SBB_singlet[ irrep ];
1157       delete [] SBB_triplet[ irrep ];
1158       delete [] SFF_singlet[ irrep ];
1159       delete [] SFF_triplet[ irrep ];
1160    }
1161    delete [] SAA;
1162    delete [] SCC;
1163    delete [] SDD;
1164    delete [] SEE;
1165    delete [] SGG;
1166    delete [] SBB_singlet;
1167    delete [] SBB_triplet;
1168    delete [] SFF_singlet;
1169    delete [] SFF_triplet;
1170 
1171    double * temp = NULL;
1172    for ( int irrep = 0; irrep < num_irreps; irrep++ ){
1173        temp = new double[ size_A[ irrep ] ]; dcopy_( size_A + irrep, FAA[ irrep ], &inc1, temp, &inc1 ); delete [] FAA[ irrep ]; FAA[ irrep ] = temp;
1174        temp = new double[ size_C[ irrep ] ]; dcopy_( size_C + irrep, FCC[ irrep ], &inc1, temp, &inc1 ); delete [] FCC[ irrep ]; FCC[ irrep ] = temp;
1175        temp = new double[ size_D[ irrep ] ]; dcopy_( size_D + irrep, FDD[ irrep ], &inc1, temp, &inc1 ); delete [] FDD[ irrep ]; FDD[ irrep ] = temp;
1176        temp = new double[ size_E[ irrep ] ]; dcopy_( size_E + irrep, FEE[ irrep ], &inc1, temp, &inc1 ); delete [] FEE[ irrep ]; FEE[ irrep ] = temp;
1177        temp = new double[ size_G[ irrep ] ]; dcopy_( size_G + irrep, FGG[ irrep ], &inc1, temp, &inc1 ); delete [] FGG[ irrep ]; FGG[ irrep ] = temp;
1178        temp = new double[ size_B_singlet[ irrep ] ]; dcopy_( size_B_singlet + irrep, FBB_singlet[ irrep ], &inc1, temp, &inc1 ); delete [] FBB_singlet[ irrep ]; FBB_singlet[ irrep ] = temp;
1179        temp = new double[ size_B_triplet[ irrep ] ]; dcopy_( size_B_triplet + irrep, FBB_triplet[ irrep ], &inc1, temp, &inc1 ); delete [] FBB_triplet[ irrep ]; FBB_triplet[ irrep ] = temp;
1180        temp = new double[ size_F_singlet[ irrep ] ]; dcopy_( size_F_singlet + irrep, FFF_singlet[ irrep ], &inc1, temp, &inc1 ); delete [] FFF_singlet[ irrep ]; FFF_singlet[ irrep ] = temp;
1181        temp = new double[ size_F_triplet[ irrep ] ]; dcopy_( size_F_triplet + irrep, FFF_triplet[ irrep ], &inc1, temp, &inc1 ); delete [] FFF_triplet[ irrep ]; FFF_triplet[ irrep ] = temp;
1182    }
1183 
1184    const int total_size = jump[ num_irreps * CHEMPS2_CASPT2_NUM_CASES ];
1185    cout << "CASPT2 : New size V_SD space  = " << total_size << endl;
1186 
1187 }
1188 
get_maxsize() const1189 int CheMPS2::CASPT2::get_maxsize() const{
1190 
1191    int maxsize = 0;
1192    for ( int irrep = 0; irrep < num_irreps; irrep++ ){
1193       maxsize = max( max( max( max( max( max( max( max( max(
1194                                         size_A[irrep],
1195                                         size_C[irrep] ),
1196                                         size_D[irrep] ),
1197                                         size_E[irrep] ),
1198                                         size_G[irrep] ),
1199                                         size_B_singlet[irrep] ),
1200                                         size_B_triplet[irrep] ),
1201                                         size_F_singlet[irrep] ),
1202                                         size_F_triplet[irrep] ),
1203                                         maxsize );
1204    }
1205    if ( maxsize <= 2 ){ maxsize = 3; }
1206    return maxsize;
1207 
1208 }
1209 
matmat(char totrans,int rowdim,int coldim,int sumdim,double alpha,double * matrix,int ldaM,double * origin,int ldaO,double * target,int ldaT)1210 void CheMPS2::CASPT2::matmat( char totrans, int rowdim, int coldim, int sumdim, double alpha, double * matrix, int ldaM, double * origin, int ldaO, double * target, int ldaT ){
1211 
1212    double add = 1.0;
1213    char notrans = 'N';
1214    dgemm_( &totrans, &notrans, &rowdim, &coldim, &sumdim, &alpha, matrix, &ldaM, origin, &ldaO, &add, target, &ldaT );
1215 
1216 }
1217 
matvec(double * vector,double * result,double * diag_fock) const1218 void CheMPS2::CASPT2::matvec( double * vector, double * result, double * diag_fock ) const{
1219 
1220    /*
1221          FOCK  | A  Bsinglet  Btriplet  C     D1     D2    Esinglet  Etriplet  Fsinglet  Ftriplet  Gsinglet  Gtriplet  Hsinglet  Htriplet
1222       ---------+-------------------------------------------------------------------------------------------------------------------------
1223       A        | OK    OK      OK       0     OK     OK    GRAD      GRAD      0         0         0         0         0         0
1224       Bsinglet | OK    OK      0        0     0      0     OK        0         0         0         0         0         0         0
1225       Btriplet | OK    0       OK       0     0      0     0         OK        0         0         0         0         0         0
1226       C        | 0     0       0        OK    OK     OK    0         0         OK        OK        GRAD      GRAD      0         0
1227       D1       | OK    0       0        OK    OK     OK    OK        OK        0         0         OK        OK        GRAD      GRAD
1228       D2       | OK    0       0        OK    OK     OK    OK        OK        0         0         OK        OK        GRAD      GRAD
1229       Esinglet | GRAD  OK      0        0     OK     OK    OK        0         0         0         0         0         OK        0
1230       Etriplet | GRAD  0       OK       0     OK     OK    0         OK        0         0         0         0         0         OK
1231       Fsinglet | 0     0       0        OK    0      0     0         0         OK        0         OK        0         0         0
1232       Ftriplet | 0     0       0        OK    0      0     0         0         0         OK        0         OK        0         0
1233       Gsinglet | 0     0       0        GRAD  OK     OK    0         0         OK        0         OK        0         OK        0
1234       Gtriplet | 0     0       0        GRAD  OK     OK    0         0         0         OK        0         OK        0         OK
1235       Hsinglet | 0     0       0        0     GRAD   GRAD  OK        0         0         0         OK        0         OK        0
1236       Htriplet | 0     0       0        0     GRAD   GRAD  0         OK        0         0         0         OK        0         OK
1237 
1238    */
1239 
1240    const int vectorlength = jump[ CHEMPS2_CASPT2_NUM_CASES * num_irreps ];
1241    #pragma omp simd
1242    for ( int elem = 0; elem < vectorlength; elem++ ){ result[ elem ] = diag_fock[ elem ] * vector[ elem ]; }
1243    const int maxlinsize = get_maxsize();
1244    double * workspace = new double[ maxlinsize * maxlinsize ];
1245    const double SQRT2 = sqrt( 2.0 );
1246 
1247    // FAD: < A(xjyz) E_wc D(aitu) > = delta_ac delta_ij FAD[ Ij ][ Ii x Ia ][ w ][ (xyz),(tu) ]
1248    for ( int IL = 0; IL < num_irreps; IL++ ){ // IL == Ii == Ij == Ix x Iy x Iz
1249       const int SIZE_L = size_A[ IL ];
1250       const int nocc_ij = indices->getNOCC( IL );
1251       for ( int IR = 0; IR < num_irreps; IR++ ){ // IR == It x Iu == Ii x Ia
1252          const int SIZE_R = size_D[ IR ];
1253          const int Iw = Irreps::directProd( IL, IR ); // Ia == Ic == Iw == IL x IR
1254          const int shift = shift_D_nonactive( indices, IL, Iw );
1255          const int nocc_w = indices->getNOCC( Iw );
1256          const int nact_w = indices->getNDMRG( Iw );
1257          const int n_oa_w = nocc_w + nact_w;
1258          const int nvir_w = indices->getNVIRT( Iw );
1259          int total_size = SIZE_L * SIZE_R;
1260          if ( total_size * nocc_ij * nact_w * nvir_w > 0 ){
1261             for ( int ac = 0; ac < nvir_w; ac++ ){
1262                for ( int cnt = 0; cnt < total_size; cnt++ ){ workspace[ cnt ] = 0.0; }
1263                for ( int w = 0; w < nact_w; w++ ){
1264                   double f_wc = fock->get( Iw, nocc_w + w, n_oa_w + ac );
1265                   int inc1 = 1;
1266                   daxpy_( &total_size, &f_wc, FAD[ IL ][ IR ][ w ], &inc1, workspace, &inc1 );
1267                }
1268                const int ptr_L = jump[ IL + num_irreps * CHEMPS2_CASPT2_A ];
1269                const int ptr_R = jump[ IR + num_irreps * CHEMPS2_CASPT2_D ] + SIZE_R * ( shift + nocc_ij * ac );
1270                matmat( 'N', SIZE_L, nocc_ij, SIZE_R, 1.0, workspace, SIZE_L, vector + ptr_R, SIZE_R, result + ptr_L, SIZE_L );
1271                matmat( 'T', SIZE_R, nocc_ij, SIZE_L, 1.0, workspace, SIZE_L, vector + ptr_L, SIZE_L, result + ptr_R, SIZE_R );
1272             }
1273          }
1274       }
1275    }
1276 
1277    // FCD: < C(bxyz) E_kw D(aitu) > = delta_ik delta_ab FCD[ Ib ][ Ii x Ia ][ w ][ (xyz),(tu) ]
1278    for ( int IL = 0; IL < num_irreps; IL++ ){ // IL == Ia == Ib
1279       const int SIZE_L = size_C[ IL ];
1280       const int nvir_ab = indices->getNVIRT( IL );
1281       for ( int IR = 0; IR < num_irreps; IR++ ){ // IR == It x Iu == Ii x Ia
1282          const int SIZE_R = size_D[ IR ];
1283          const int Iw = Irreps::directProd( IL, IR ); // Ii == Ik == Iw == IL x IR
1284          const int shift = shift_D_nonactive( indices, Iw, IL );
1285          const int nocc_w = indices->getNOCC( Iw );
1286          const int nact_w = indices->getNDMRG( Iw );
1287          int total_size = SIZE_L * SIZE_R;
1288          if ( total_size * nvir_ab * nact_w * nocc_w > 0 ){
1289             for ( int ik = 0; ik < nocc_w; ik++ ){
1290                for ( int cnt = 0; cnt < total_size; cnt++ ){ workspace[ cnt ] = 0.0; }
1291                for ( int w = 0; w < nact_w; w++ ){
1292                   double f_kw = fock->get( Iw, ik, nocc_w + w );
1293                   int inc1 = 1;
1294                   daxpy_( &total_size, &f_kw, FCD[ IL ][ IR ][ w ], &inc1, workspace, &inc1 );
1295                }
1296                const int ptr_L = jump[ IL + num_irreps * CHEMPS2_CASPT2_C ];
1297                const int ptr_R = jump[ IR + num_irreps * CHEMPS2_CASPT2_D ] + SIZE_R * ( shift + ik );
1298                const int LDA_R = SIZE_R * nocc_w;
1299                matmat( 'N', SIZE_L, nvir_ab, SIZE_R, 1.0, workspace, SIZE_L, vector + ptr_R, LDA_R,  result + ptr_L, SIZE_L );
1300                matmat( 'T', SIZE_R, nvir_ab, SIZE_L, 1.0, workspace, SIZE_L, vector + ptr_L, SIZE_L, result + ptr_R, LDA_R  );
1301             }
1302          }
1303       }
1304    }
1305 
1306    // FAB singlet: < A(xlyz) E_kw SB_tiuj > = ( delta_ik delta_jl + delta_jk delta_il ) / sqrt( 1 + delta_ij ) * FAB_singlet[ Il ][ Ii x Ij ][ w ][ (xyz),(tu) ]
1307    for ( int IL = 0; IL < num_irreps; IL++ ){ // IL == Il == Ix x Iy x Iz
1308       const int SIZE_L = size_A[ IL ];
1309       const int nocc_l = indices->getNOCC( IL );
1310       for ( int IR = 0; IR < num_irreps; IR++ ){ // IR == It x Iu == Ii x Ij
1311          const int SIZE_R = size_B_singlet[ IR ];
1312          const int Iw = Irreps::directProd( IL, IR ); // Iw == Ik
1313          const int shift = (( Iw < IL ) ? shift_B_nonactive( indices, Iw, IL, +1 ) : shift_B_nonactive( indices, IL, Iw, +1 ));
1314          const int nocc_w = indices->getNOCC( Iw );
1315          const int nact_w = indices->getNDMRG( Iw );
1316          int total_size = SIZE_L * SIZE_R;
1317          if ( total_size * nocc_l * nact_w * nocc_w > 0 ){
1318             for ( int k = 0; k < nocc_w; k++ ){
1319                for ( int cnt = 0; cnt < total_size; cnt++ ){ workspace[ cnt ] = 0.0; }
1320                for ( int w = 0; w < nact_w; w++ ){
1321                   double f_kw = fock->get( Iw, k, nocc_w + w );
1322                   int inc1 = 1;
1323                   daxpy_( &total_size, &f_kw, FAB_singlet[ IL ][ IR ][ w ], &inc1, workspace, &inc1 );
1324                }
1325                if ( IR == 0 ){ // Ii == Ij  and  Ik == Il
1326                   if ( k > 0 ){
1327                      const int ptr_L = jump[ IL + num_irreps * CHEMPS2_CASPT2_A         ];
1328                      const int ptr_R = jump[ IR + num_irreps * CHEMPS2_CASPT2_B_SINGLET ] + SIZE_R * ( shift + ( k * ( k + 1 ) ) / 2 );
1329                      matmat( 'N', SIZE_L, k, SIZE_R, 1.0, workspace, SIZE_L, vector + ptr_R, SIZE_R, result + ptr_L, SIZE_L );
1330                      matmat( 'T', SIZE_R, k, SIZE_L, 1.0, workspace, SIZE_L, vector + ptr_L, SIZE_L, result + ptr_R, SIZE_R );
1331                   }
1332                   #pragma omp parallel for schedule(static)
1333                   for ( int l = k; l < nocc_l; l++ ){
1334                      const int ptr_L = jump[ IL + num_irreps * CHEMPS2_CASPT2_A         ] + SIZE_L * l;
1335                      const int ptr_R = jump[ IR + num_irreps * CHEMPS2_CASPT2_B_SINGLET ] + SIZE_R * ( shift + k + ( l * ( l + 1 ) ) / 2 );
1336                      const double factor = (( k == l ) ? SQRT2 : 1.0 );
1337                      matmat( 'N', SIZE_L, 1, SIZE_R, factor, workspace, SIZE_L, vector + ptr_R, SIZE_R, result + ptr_L, SIZE_L );
1338                      matmat( 'T', SIZE_R, 1, SIZE_L, factor, workspace, SIZE_L, vector + ptr_L, SIZE_L, result + ptr_R, SIZE_R );
1339                   }
1340                } else {
1341                   const int ptr_L = jump[ IL + num_irreps * CHEMPS2_CASPT2_A         ];
1342                   const int ptr_R = jump[ IR + num_irreps * CHEMPS2_CASPT2_B_SINGLET ] + SIZE_R * ( shift + (( Iw < IL ) ? k : nocc_l * k ));
1343                   const int LDA_R = (( Iw < IL ) ? SIZE_R * nocc_w : SIZE_R );
1344                   matmat( 'N', SIZE_L, nocc_l, SIZE_R, 1.0, workspace, SIZE_L, vector + ptr_R, LDA_R,  result + ptr_L, SIZE_L );
1345                   matmat( 'T', SIZE_R, nocc_l, SIZE_L, 1.0, workspace, SIZE_L, vector + ptr_L, SIZE_L, result + ptr_R, LDA_R  );
1346                }
1347             }
1348          }
1349       }
1350    }
1351 
1352    // FAB triplet: < A(xlyz) E_kw TB_tiuj > = ( delta_ik delta_jl - delta_jk delta_il ) * FAB_triplet[ Il ][ Ii x Ij ][ w ][ (xyz),(tu) ]
1353    for ( int IL = 0; IL < num_irreps; IL++ ){ // IL == Il == Ix x Iy x Iz
1354       const int SIZE_L = size_A[ IL ];
1355       const int nocc_l = indices->getNOCC( IL );
1356       for ( int IR = 0; IR < num_irreps; IR++ ){ // IR == It x Iu == Ii x Ij
1357          const int SIZE_R = size_B_triplet[ IR ];
1358          const int Iw = Irreps::directProd( IL, IR ); // Iw == Ik
1359          const int shift = (( Iw < IL ) ? shift_B_nonactive( indices, Iw, IL, -1 ) : shift_B_nonactive( indices, IL, Iw, -1 ));
1360          const int nocc_w = indices->getNOCC( Iw );
1361          const int nact_w = indices->getNDMRG( Iw );
1362          int total_size = SIZE_L * SIZE_R;
1363          if ( total_size * nocc_l * nact_w * nocc_w > 0 ){
1364             for ( int k = 0; k < nocc_w; k++ ){
1365                for ( int cnt = 0; cnt < total_size; cnt++ ){ workspace[ cnt ] = 0.0; }
1366                for ( int w = 0; w < nact_w; w++ ){
1367                   double f_kw = fock->get( Iw, k, nocc_w + w );
1368                   int inc1 = 1;
1369                   daxpy_( &total_size, &f_kw, FAB_triplet[ IL ][ IR ][ w ], &inc1, workspace, &inc1 );
1370                }
1371                if ( IR == 0 ){ // Ii == Ij  and  Ik == Il
1372                   if ( k > 0 ){ // ( k > l  --->  - delta_jk delta_il )
1373                      const int ptr_L = jump[ IL + num_irreps * CHEMPS2_CASPT2_A         ];
1374                      const int ptr_R = jump[ IR + num_irreps * CHEMPS2_CASPT2_B_TRIPLET ] + SIZE_R * ( shift + ( k * ( k - 1 ) ) / 2 );
1375                      matmat( 'N', SIZE_L, k, SIZE_R, -1.0, workspace, SIZE_L, vector + ptr_R, SIZE_R, result + ptr_L, SIZE_L );
1376                      matmat( 'T', SIZE_R, k, SIZE_L, -1.0, workspace, SIZE_L, vector + ptr_L, SIZE_L, result + ptr_R, SIZE_R );
1377                   }
1378                   #pragma omp parallel for schedule(static)
1379                   for ( int l = k+1; l < nocc_l; l++ ){
1380                      const int ptr_L = jump[ IL + num_irreps * CHEMPS2_CASPT2_A         ] + SIZE_L * l;
1381                      const int ptr_R = jump[ IR + num_irreps * CHEMPS2_CASPT2_B_TRIPLET ] + SIZE_R * ( shift + k + ( l * ( l - 1 ) ) / 2 );
1382                      matmat( 'N', SIZE_L, 1, SIZE_R, 1.0, workspace, SIZE_L, vector + ptr_R, SIZE_R, result + ptr_L, SIZE_L ); // ( k < l  --->  + delta_ik delta_jl )
1383                      matmat( 'T', SIZE_R, 1, SIZE_L, 1.0, workspace, SIZE_L, vector + ptr_L, SIZE_L, result + ptr_R, SIZE_R );
1384                   }
1385                } else {
1386                   const int ptr_L = jump[ IL + num_irreps * CHEMPS2_CASPT2_A         ];
1387                   const int ptr_R = jump[ IR + num_irreps * CHEMPS2_CASPT2_B_TRIPLET ] + SIZE_R * ( shift + (( Iw < IL ) ? k : nocc_l * k ));
1388                   const int LDA_R = (( Iw < IL ) ? SIZE_R * nocc_w : SIZE_R );
1389                   const double factor = (( Iw < IL ) ? 1.0 : -1.0 ); // ( k < l  --->  + delta_ik delta_jl ) and ( k > l  --->  - delta_jk delta_il )
1390                   matmat( 'N', SIZE_L, nocc_l, SIZE_R, factor, workspace, SIZE_L, vector + ptr_R, LDA_R,  result + ptr_L, SIZE_L );
1391                   matmat( 'T', SIZE_R, nocc_l, SIZE_L, factor, workspace, SIZE_L, vector + ptr_L, SIZE_L, result + ptr_R, LDA_R  );
1392                }
1393             }
1394          }
1395       }
1396    }
1397 
1398    // FCF singlet: < C(dxyz) E_wc SF_atbu > = ( delta_ac delta_bd + delta_ad delta_bc ) / sqrt( 1 + delta_ab ) * FCF_singlet[ Id ][ Ia x Ib ][ w ][ (xyz),(tu) ]
1399    for ( int IL = 0; IL < num_irreps; IL++ ){ // IL == Id == Ix x Iy x Iz
1400       const int SIZE_L = size_C[ IL ];
1401       const int nvir_d = indices->getNVIRT( IL );
1402       for ( int IR = 0; IR < num_irreps; IR++ ){ // IR == It x Iu == Ia x Ib
1403          const int SIZE_R = size_F_singlet[ IR ];
1404          const int Iw = Irreps::directProd( IL, IR ); // Iw == Ic
1405          const int shift = (( Iw < IL ) ? shift_F_nonactive( indices, Iw, IL, +1 ) : shift_F_nonactive( indices, IL, Iw, +1 ));
1406          const int nocc_w = indices->getNOCC( Iw );
1407          const int nact_w = indices->getNDMRG( Iw );
1408          const int n_oa_w = nocc_w + nact_w;
1409          const int nvir_w = indices->getNVIRT( Iw );
1410          int total_size = SIZE_L * SIZE_R;
1411          if ( total_size * nvir_d * nact_w * nvir_w > 0 ){
1412             for ( int c = 0; c < nvir_w; c++ ){
1413                for ( int cnt = 0; cnt < total_size; cnt++ ){ workspace[ cnt ] = 0.0; }
1414                for ( int w = 0; w < nact_w; w++ ){
1415                   double f_wc = fock->get( Iw, nocc_w + w, n_oa_w + c );
1416                   int inc1 = 1;
1417                   daxpy_( &total_size, &f_wc, FCF_singlet[ IL ][ IR ][ w ], &inc1, workspace, &inc1 );
1418                }
1419                if ( IR == 0 ){ // Ia == Ib  and  Ic == Id
1420                   if ( c > 0 ){
1421                      const int ptr_L = jump[ IL + num_irreps * CHEMPS2_CASPT2_C         ];
1422                      const int ptr_R = jump[ IR + num_irreps * CHEMPS2_CASPT2_F_SINGLET ] + SIZE_R * ( shift + ( c * ( c + 1 ) ) / 2 );
1423                      matmat( 'N', SIZE_L, c, SIZE_R, 1.0, workspace, SIZE_L, vector + ptr_R, SIZE_R, result + ptr_L, SIZE_L );
1424                      matmat( 'T', SIZE_R, c, SIZE_L, 1.0, workspace, SIZE_L, vector + ptr_L, SIZE_L, result + ptr_R, SIZE_R );
1425                   }
1426                   #pragma omp parallel for schedule(static)
1427                   for ( int d = c; d < nvir_d; d++ ){
1428                      const int ptr_L = jump[ IL + num_irreps * CHEMPS2_CASPT2_C         ] + SIZE_L * d;
1429                      const int ptr_R = jump[ IR + num_irreps * CHEMPS2_CASPT2_F_SINGLET ] + SIZE_R * ( shift + c + ( d * ( d + 1 ) ) / 2 );
1430                      const double factor = (( c == d ) ? SQRT2 : 1.0 );
1431                      matmat( 'N', SIZE_L, 1, SIZE_R, factor, workspace, SIZE_L, vector + ptr_R, SIZE_R, result + ptr_L, SIZE_L );
1432                      matmat( 'T', SIZE_R, 1, SIZE_L, factor, workspace, SIZE_L, vector + ptr_L, SIZE_L, result + ptr_R, SIZE_R );
1433                   }
1434                } else {
1435                   const int ptr_L = jump[ IL + num_irreps * CHEMPS2_CASPT2_C         ];
1436                   const int ptr_R = jump[ IR + num_irreps * CHEMPS2_CASPT2_F_SINGLET ] + SIZE_R * ( shift + (( Iw < IL ) ? c : nvir_d * c ));
1437                   const int LDA_R = (( Iw < IL ) ? SIZE_R * nvir_w : SIZE_R );
1438                   matmat( 'N', SIZE_L, nvir_d, SIZE_R, 1.0, workspace, SIZE_L, vector + ptr_R, LDA_R,  result + ptr_L, SIZE_L );
1439                   matmat( 'T', SIZE_R, nvir_d, SIZE_L, 1.0, workspace, SIZE_L, vector + ptr_L, SIZE_L, result + ptr_R, LDA_R  );
1440                }
1441             }
1442          }
1443       }
1444    }
1445 
1446    // FCF triplet: < C(dxyz) E_wc TF_atbu > = ( delta_ac delta_bd - delta_ad delta_bc ) * FCF_triplet[ Id ][ Ia x Ib ][ w ][ (xyz),(tu) ]
1447    for ( int IL = 0; IL < num_irreps; IL++ ){ // IL == Id == Ix x Iy x Iz
1448       const int SIZE_L = size_C[ IL ];
1449       const int nvir_d = indices->getNVIRT( IL );
1450       for ( int IR = 0; IR < num_irreps; IR++ ){ // IR == It x Iu == Ia x Ib
1451          const int SIZE_R = size_F_triplet[ IR ];
1452          const int Iw = Irreps::directProd( IL, IR ); // Iw == Ic
1453          const int shift = (( Iw < IL ) ? shift_F_nonactive( indices, Iw, IL, -1 ) : shift_F_nonactive( indices, IL, Iw, -1 ));
1454          const int nocc_w = indices->getNOCC( Iw );
1455          const int nact_w = indices->getNDMRG( Iw );
1456          const int n_oa_w = nocc_w + nact_w;
1457          const int nvir_w = indices->getNVIRT( Iw );
1458          int total_size = SIZE_L * SIZE_R;
1459          if ( total_size * nvir_d * nact_w * nvir_w > 0 ){
1460             for ( int c = 0; c < nvir_w; c++ ){
1461                for ( int cnt = 0; cnt < total_size; cnt++ ){ workspace[ cnt ] = 0.0; }
1462                for ( int w = 0; w < nact_w; w++ ){
1463                   double f_wc = fock->get( Iw, nocc_w + w, n_oa_w + c );
1464                   int inc1 = 1;
1465                   daxpy_( &total_size, &f_wc, FCF_triplet[ IL ][ IR ][ w ], &inc1, workspace, &inc1 );
1466                }
1467                if ( IR == 0 ){ // Ia == Ib  and  Ic == Id
1468                   if ( c > 0 ){ // ( c > d  --->  - delta_ad delta_bc )
1469                      const int ptr_L = jump[ IL + num_irreps * CHEMPS2_CASPT2_C         ];
1470                      const int ptr_R = jump[ IR + num_irreps * CHEMPS2_CASPT2_F_TRIPLET ] + SIZE_R * ( shift + ( c * ( c - 1 ) ) / 2 );
1471                      matmat( 'N', SIZE_L, c, SIZE_R, -1.0, workspace, SIZE_L, vector + ptr_R, SIZE_R, result + ptr_L, SIZE_L );
1472                      matmat( 'T', SIZE_R, c, SIZE_L, -1.0, workspace, SIZE_L, vector + ptr_L, SIZE_L, result + ptr_R, SIZE_R );
1473                   }
1474                   #pragma omp parallel for schedule(static)
1475                   for ( int d = c+1; d < nvir_d; d++ ){
1476                      const int ptr_L = jump[ IL + num_irreps * CHEMPS2_CASPT2_C         ] + SIZE_L * d;
1477                      const int ptr_R = jump[ IR + num_irreps * CHEMPS2_CASPT2_F_TRIPLET ] + SIZE_R * ( shift + c + ( d * ( d - 1 ) ) / 2 );
1478                      matmat( 'N', SIZE_L, 1, SIZE_R, 1.0, workspace, SIZE_L, vector + ptr_R, SIZE_R, result + ptr_L, SIZE_L ); // ( c < d  --->  + delta_ac delta_bd )
1479                      matmat( 'T', SIZE_R, 1, SIZE_L, 1.0, workspace, SIZE_L, vector + ptr_L, SIZE_L, result + ptr_R, SIZE_R );
1480                   }
1481                } else {
1482                   const double factor = (( Iw < IL ) ? 1.0 : -1.0 ); // ( c < d  --->  + delta_ac delta_bd ) and ( c > d  --->  - delta_ad delta_bc )
1483                   const int ptr_L = jump[ IL + num_irreps * CHEMPS2_CASPT2_C         ];
1484                   const int ptr_R = jump[ IR + num_irreps * CHEMPS2_CASPT2_F_TRIPLET ] + SIZE_R * ( shift + (( Iw < IL ) ? c : nvir_d * c ));
1485                   const int LDA_R = (( Iw < IL ) ? SIZE_R * nvir_w : SIZE_R );
1486                   matmat( 'N', SIZE_L, nvir_d, SIZE_R, factor, workspace, SIZE_L, vector + ptr_R, LDA_R,  result + ptr_L, SIZE_L );
1487                   matmat( 'T', SIZE_R, nvir_d, SIZE_L, factor, workspace, SIZE_L, vector + ptr_L, SIZE_L, result + ptr_R, LDA_R  );
1488                }
1489             }
1490          }
1491       }
1492    }
1493 
1494    // FBE singlet: < SB_xkyl E_wc SE_tiaj > = 2 delta_ac delta_ik delta_jl FBE_singlet[ Ik x Il ][ It ][ w ][ xy, t ]
1495    for ( int IL = 0; IL < num_irreps; IL++ ){ // IL == Iik x Ijl == Ix x Iy
1496       const int SIZE_L = size_B_singlet[ IL ];
1497       for ( int IR = 0; IR < num_irreps; IR++ ){ // IR == It == Iik x Ijl x Iac
1498          const int SIZE_R = size_E[ IR ];
1499          const int Iw = Irreps::directProd( IL, IR ); // Iw == Iac
1500          const int nocc_w = indices->getNOCC( Iw );
1501          const int nact_w = indices->getNDMRG( Iw );
1502          const int n_oa_w = nocc_w + nact_w;
1503          const int nvir_w = indices->getNVIRT( Iw );
1504          int linsize = 0;
1505          for ( int Iik = 0; Iik < num_irreps; Iik++ ){
1506             const int Ijl = Irreps::directProd( Iik, IL );
1507             if ( Iik <= Ijl ){
1508                const int nocc_ik = indices->getNOCC( Iik );
1509                linsize += (( IL == 0 ) ? ( nocc_ik * ( nocc_ik + 1 ) ) / 2 : nocc_ik * indices->getNOCC( Ijl ) );
1510             }
1511          }
1512          const int size_ij = linsize;
1513          int total_size = SIZE_L * SIZE_R;
1514          if ( total_size * nact_w * nvir_w * size_ij > 0 ){
1515             const int shift_E = shift_E_nonactive( indices, Iw, 0, IL, +1 );
1516             for ( int ac = 0; ac < nvir_w; ac++ ){
1517                for ( int cnt = 0; cnt < total_size; cnt++ ){ workspace[ cnt ] = 0.0; }
1518                for ( int w = 0; w < nact_w; w++ ){
1519                   double f_wc = fock->get( Iw, nocc_w + w, n_oa_w + ac );
1520                   int inc1 = 1;
1521                   daxpy_( &total_size, &f_wc, FBE_singlet[ IL ][ IR ][ w ], &inc1, workspace, &inc1 );
1522                }
1523                const int ptr_L = jump[ IL + num_irreps * CHEMPS2_CASPT2_B_SINGLET ];
1524                const int ptr_R = jump[ IR + num_irreps * CHEMPS2_CASPT2_E_SINGLET ] + SIZE_R * ( shift_E + ac );
1525                const int LDA_R = SIZE_R * nvir_w;
1526                matmat( 'N', SIZE_L, size_ij, SIZE_R, 2.0, workspace, SIZE_L, vector + ptr_R, LDA_R,  result + ptr_L, SIZE_L );
1527                matmat( 'T', SIZE_R, size_ij, SIZE_L, 2.0, workspace, SIZE_L, vector + ptr_L, SIZE_L, result + ptr_R, LDA_R  );
1528             }
1529          }
1530       }
1531    }
1532 
1533    // FBE triplet: < TB_xkyl E_wc TE_tiaj > = 2 delta_ac delta_ik delta_jl FBE_triplet[ Ik x Il ][ It ][ w ][ xy, t ]
1534    for ( int IL = 0; IL < num_irreps; IL++ ){ // IL == Iik x Ijl == Ix x Iy
1535       const int SIZE_L = size_B_triplet[ IL ];
1536       for ( int IR = 0; IR < num_irreps; IR++ ){ // IR == It == Iik x Ijl x Iac
1537          const int SIZE_R = size_E[ IR ];
1538          const int Iw = Irreps::directProd( IL, IR ); // Iw == Iac
1539          const int nocc_w = indices->getNOCC( Iw );
1540          const int nact_w = indices->getNDMRG( Iw );
1541          const int n_oa_w = nocc_w + nact_w;
1542          const int nvir_w = indices->getNVIRT( Iw );
1543          int linsize = 0;
1544          for ( int Iik = 0; Iik < num_irreps; Iik++ ){
1545             const int Ijl = Irreps::directProd( Iik, IL );
1546             if ( Iik <= Ijl ){
1547                const int nocc_ik = indices->getNOCC( Iik );
1548                linsize += (( IL == 0 ) ? ( nocc_ik * ( nocc_ik - 1 ) ) / 2 : nocc_ik * indices->getNOCC( Ijl ) );
1549             }
1550          }
1551          const int size_ij = linsize;
1552          int total_size = SIZE_L * SIZE_R;
1553          if ( total_size * nact_w * nvir_w * size_ij > 0 ){
1554             const int shift_E = shift_E_nonactive( indices, Iw, 0, IL, -1 );
1555             for ( int ac = 0; ac < nvir_w; ac++ ){
1556                for ( int cnt = 0; cnt < total_size; cnt++ ){ workspace[ cnt ] = 0.0; }
1557                for ( int w = 0; w < nact_w; w++ ){
1558                   double f_wc = fock->get( Iw, nocc_w + w, n_oa_w + ac );
1559                   int inc1 = 1;
1560                   daxpy_( &total_size, &f_wc, FBE_triplet[ IL ][ IR ][ w ], &inc1, workspace, &inc1 );
1561                }
1562                const int ptr_L = jump[ IL + num_irreps * CHEMPS2_CASPT2_B_TRIPLET ];
1563                const int ptr_R = jump[ IR + num_irreps * CHEMPS2_CASPT2_E_TRIPLET ] + SIZE_R * ( shift_E + ac );
1564                const int LDA_R = SIZE_R * nvir_w;
1565                matmat( 'N', SIZE_L, size_ij, SIZE_R, 2.0, workspace, SIZE_L, vector + ptr_R, LDA_R,  result + ptr_L, SIZE_L );
1566                matmat( 'T', SIZE_R, size_ij, SIZE_L, 2.0, workspace, SIZE_L, vector + ptr_L, SIZE_L, result + ptr_R, LDA_R  );
1567             }
1568          }
1569       }
1570    }
1571 
1572    // FFG singlet: < SF_cxdy E_kw SG_aibt > = 2 delta_ac delta_bd delta_ik FFG_singlet[ Ic x Id ][ It ][ w ][ xy, t ]
1573    for ( int IL = 0; IL < num_irreps; IL++ ){ // IL == Iac x Ibd == Ix x Iy
1574       const int SIZE_L = size_F_singlet[ IL ];
1575       for ( int IR = 0; IR < num_irreps; IR++ ){ // IR == It == Iac x Ibd x Iik
1576          const int SIZE_R = size_G[ IR ];
1577          const int Iw = Irreps::directProd( IL, IR ); // Iw == Iik
1578          const int nocc_w = indices->getNOCC( Iw );
1579          const int nact_w = indices->getNDMRG( Iw );
1580          int linsize = 0;
1581          for ( int Iac = 0; Iac < num_irreps; Iac++ ){
1582             const int Ibd = Irreps::directProd( Iac, IL );
1583             if ( Iac <= Ibd ){
1584                const int nvir_ac = indices->getNVIRT( Iac );
1585                linsize += (( IL == 0 ) ? ( nvir_ac * ( nvir_ac + 1 ) ) / 2 : nvir_ac * indices->getNVIRT( Ibd ) );
1586             }
1587          }
1588          const int size_ab = linsize;
1589          int total_size = SIZE_L * SIZE_R;
1590          if ( total_size * nact_w * nocc_w * size_ab > 0 ){
1591             const int shift_G = shift_G_nonactive( indices, Iw, 0, IL, +1 );
1592             for ( int ik = 0; ik < nocc_w; ik++ ){
1593                for ( int cnt = 0; cnt < total_size; cnt++ ){ workspace[ cnt ] = 0.0; }
1594                for ( int w = 0; w < nact_w; w++ ){
1595                   double f_kw = fock->get( Iw, ik, nocc_w + w );
1596                   int inc1 = 1;
1597                   daxpy_( &total_size, &f_kw, FFG_singlet[ IL ][ IR ][ w ], &inc1, workspace, &inc1 );
1598                }
1599                const int ptr_L = jump[ IL + num_irreps * CHEMPS2_CASPT2_F_SINGLET ];
1600                const int ptr_R = jump[ IR + num_irreps * CHEMPS2_CASPT2_G_SINGLET ] + SIZE_R * ( shift_G + ik );
1601                const int LDA_R = SIZE_R * nocc_w;
1602                matmat( 'N', SIZE_L, size_ab, SIZE_R, 2.0, workspace, SIZE_L, vector + ptr_R, LDA_R,  result + ptr_L, SIZE_L );
1603                matmat( 'T', SIZE_R, size_ab, SIZE_L, 2.0, workspace, SIZE_L, vector + ptr_L, SIZE_L, result + ptr_R, LDA_R  );
1604             }
1605          }
1606       }
1607    }
1608 
1609    // FFG triplet: < TF_cxdy E_kw TG_aibt > = 2 delta_ac delta_bd delta_ik FFG_triplet[ Ic x Id ][ It ][ w ][ xy, t ]
1610    for ( int IL = 0; IL < num_irreps; IL++ ){ // IL == Iac x Ibd == Ix x Iy
1611       const int SIZE_L = size_F_triplet[ IL ];
1612       for ( int IR = 0; IR < num_irreps; IR++ ){ // IR == It == Iac x Ibd x Iik
1613          const int SIZE_R = size_G[ IR ];
1614          const int Iw = Irreps::directProd( IL, IR ); // Iw == Iik
1615          const int nocc_w = indices->getNOCC( Iw );
1616          const int nact_w = indices->getNDMRG( Iw );
1617          int linsize = 0;
1618          for ( int Iac = 0; Iac < num_irreps; Iac++ ){
1619             const int Ibd = Irreps::directProd( Iac, IL );
1620             if ( Iac <= Ibd ){
1621                const int nvir_ac = indices->getNVIRT( Iac );
1622                linsize += (( IL == 0 ) ? ( nvir_ac * ( nvir_ac - 1 ) ) / 2 : nvir_ac * indices->getNVIRT( Ibd ) );
1623             }
1624          }
1625          const int size_ab = linsize;
1626          int total_size = SIZE_L * SIZE_R;
1627          if ( total_size * nact_w * nocc_w * size_ab > 0 ){
1628             const int shift_G = shift_G_nonactive( indices, Iw, 0, IL, -1 );
1629             for ( int ik = 0; ik < nocc_w; ik++ ){
1630                for ( int cnt = 0; cnt < total_size; cnt++ ){ workspace[ cnt ] = 0.0; }
1631                for ( int w = 0; w < nact_w; w++ ){
1632                   double f_kw = fock->get( Iw, ik, nocc_w + w );
1633                   int inc1 = 1;
1634                   daxpy_( &total_size, &f_kw, FFG_triplet[ IL ][ IR ][ w ], &inc1, workspace, &inc1 );
1635                }
1636                const int ptr_L = jump[ IL + num_irreps * CHEMPS2_CASPT2_F_TRIPLET ];
1637                const int ptr_R = jump[ IR + num_irreps * CHEMPS2_CASPT2_G_TRIPLET ] + SIZE_R * ( shift_G + ik );
1638                const int LDA_R = SIZE_R * nocc_w;
1639                matmat( 'N', SIZE_L, size_ab, SIZE_R, 2.0, workspace, SIZE_L, vector + ptr_R, LDA_R,  result + ptr_L, SIZE_L );
1640                matmat( 'T', SIZE_R, size_ab, SIZE_L, 2.0, workspace, SIZE_L, vector + ptr_L, SIZE_L, result + ptr_R, LDA_R  );
1641             }
1642          }
1643       }
1644    }
1645 
1646    // FEH singlet: < SE_xkdl E_wc SH_aibj > = 2 delta_ik delta_jl ( delta_ac delta_bd + delta_ad delta_bc ) / sqrt( 1 + delta_ab ) FEH[ Ix ][ w ][ x ]
1647    for ( int IL = 0; IL < num_irreps; IL++ ){ // IL == Ixw == Ic == Iik x Ijl x Id
1648       int SIZE = size_E[ IL ];
1649       const int nocc_w = indices->getNOCC( IL );
1650       const int nact_w = indices->getNDMRG( IL );
1651       const int n_oa_w = nocc_w + nact_w;
1652       const int nvir_w = indices->getNVIRT( IL );
1653       if ( SIZE * nact_w * nvir_w > 0 ){
1654          for ( int c = 0; c < nvir_w; c++ ){
1655 
1656             // workspace_c[ x ] = sum_w f_wc FEH[ Ixw == IL == Ic ][ w ][ x ]
1657             for ( int cnt = 0; cnt < SIZE; cnt++ ){ workspace[ cnt ] = 0.0; }
1658             for ( int w = 0; w < nact_w; w++ ){
1659                double f_wc = fock->get( IL, nocc_w + w, n_oa_w + c );
1660                int inc1 = 1;
1661                daxpy_( &SIZE, &f_wc, FEH[ IL ][ w ], &inc1, workspace, &inc1 );
1662             }
1663 
1664             for ( int Id = 0; Id < num_irreps; Id++ ){
1665 
1666                const int Icenter = Irreps::directProd( IL, Id );
1667                const int nvir_d = indices->getNVIRT( Id );
1668                const int LDA_E = SIZE * nvir_d;
1669 
1670                if ( IL < Id ){ // Ic < Id
1671                   for ( int Ii = 0; Ii < num_irreps; Ii++ ){
1672                      const int Ij = Irreps::directProd( Ii, Icenter );
1673                      if ( Ii < Ij ){
1674                         const int jump_E = jump[ IL + num_irreps * CHEMPS2_CASPT2_E_SINGLET ] + SIZE * shift_E_nonactive( indices, Id, Ii, Ij, +1 );
1675                         const int jump_H = jump[ Icenter + num_irreps * CHEMPS2_CASPT2_H_SINGLET ] + shift_H_nonactive( indices, Ii, Ij, IL, Id, +1 );
1676                         const int size_ij = indices->getNOCC( Ii ) * indices->getNOCC( Ij );
1677                         #pragma omp parallel for schedule(static)
1678                         for ( int d = 0; d < nvir_d; d++ ){
1679                            const int ptr_H = jump_H + size_ij * ( c + nvir_w * d );
1680                            const int ptr_E = jump_E + SIZE * d;
1681                            matmat( 'N', SIZE, size_ij, 1,    2.0, workspace, SIZE, vector + ptr_H, 1,     result + ptr_E, LDA_E );
1682                            matmat( 'T', 1,    size_ij, SIZE, 2.0, workspace, SIZE, vector + ptr_E, LDA_E, result + ptr_H, 1     );
1683                         }
1684                      }
1685                   }
1686                }
1687 
1688                if ( IL == Id ){ // Ic == Id == Ia == Ib --> Iik == Ijl
1689                   for ( int Iij = 0; Iij < num_irreps; Iij++ ){
1690                      const int jump_E = jump[ IL + num_irreps * CHEMPS2_CASPT2_E_SINGLET ] + SIZE * shift_E_nonactive( indices, Id,  Iij, Iij, +1 );
1691                      const int jump_H = jump[ Icenter + num_irreps * CHEMPS2_CASPT2_H_SINGLET ] + shift_H_nonactive( indices, Iij, Iij, IL,  Id, +1 );
1692                      const int size_ij = ( indices->getNOCC( Iij ) * ( indices->getNOCC( Iij ) + 1 ) ) / 2;
1693                      #pragma omp parallel for schedule(static)
1694                      for ( int d = 0; d < c; d++ ){
1695                         const int ptr_H = jump_H + size_ij * ( d + ( c * ( c + 1 ) ) / 2 );
1696                         const int ptr_E = jump_E + SIZE * d;
1697                         matmat( 'N', SIZE, size_ij, 1,    2.0, workspace, SIZE, vector + ptr_H, 1,     result + ptr_E, LDA_E );
1698                         matmat( 'T', 1,    size_ij, SIZE, 2.0, workspace, SIZE, vector + ptr_E, LDA_E, result + ptr_H, 1     );
1699                      }
1700                      #pragma omp parallel for schedule(static)
1701                      for ( int d = c; d < nvir_d; d++ ){
1702                         const int ptr_H = jump_H + size_ij * ( c + ( d * ( d + 1 ) ) / 2 );
1703                         const int ptr_E = jump_E + SIZE * d;
1704                         const double factor = 2 * (( c == d ) ? SQRT2 : 1.0 );
1705                         matmat( 'N', SIZE, size_ij, 1,    factor, workspace, SIZE, vector + ptr_H, 1,     result + ptr_E, LDA_E );
1706                         matmat( 'T', 1,    size_ij, SIZE, factor, workspace, SIZE, vector + ptr_E, LDA_E, result + ptr_H, 1     );
1707                      }
1708                   }
1709                }
1710 
1711                if ( IL > Id ){ // Ic > Id
1712                   for ( int Ii = 0; Ii < num_irreps; Ii++ ){
1713                      const int Ij = Irreps::directProd( Ii, Icenter );
1714                      if ( Ii < Ij ){
1715                         const int jump_E = jump[ IL + num_irreps * CHEMPS2_CASPT2_E_SINGLET ] + SIZE * shift_E_nonactive( indices, Id, Ii, Ij, +1 );
1716                         const int jump_H = jump[ Icenter + num_irreps * CHEMPS2_CASPT2_H_SINGLET ] + shift_H_nonactive( indices, Ii, Ij, Id, IL, +1 );
1717                         const int size_ij = indices->getNOCC( Ii ) * indices->getNOCC( Ij );
1718                         #pragma omp parallel for schedule(static)
1719                         for ( int d = 0; d < nvir_d; d++ ){
1720                            const int ptr_H = jump_H + size_ij * ( d + nvir_d * c );
1721                            const int ptr_E = jump_E + SIZE * d;
1722                            matmat( 'N', SIZE, size_ij, 1,    2.0, workspace, SIZE, vector + ptr_H, 1,     result + ptr_E, LDA_E );
1723                            matmat( 'T', 1,    size_ij, SIZE, 2.0, workspace, SIZE, vector + ptr_E, LDA_E, result + ptr_H, 1     );
1724                         }
1725                      }
1726                   }
1727                }
1728             }
1729          }
1730       }
1731    }
1732 
1733    // FEH triplet: < TE_xkdl E_wc TH_aibj > = 6 delta_ik delta_jl ( delta_ac delta_bd - delta_ad delta_bc ) / sqrt( 1 + delta_ab ) FEH[ Ix ][ w ][ x ]
1734    for ( int IL = 0; IL < num_irreps; IL++ ){ // IL == Ixw == Ic == Iik x Ijl x Id
1735       int SIZE = size_E[ IL ];
1736       const int nocc_w = indices->getNOCC( IL );
1737       const int nact_w = indices->getNDMRG( IL );
1738       const int n_oa_w = nocc_w + nact_w;
1739       const int nvir_w = indices->getNVIRT( IL );
1740       if ( SIZE * nact_w * nvir_w > 0 ){
1741          for ( int c = 0; c < nvir_w; c++ ){
1742 
1743             // workspace_c[ x ] = sum_w f_wc FEH[ Ixw == IL == Ic ][ w ][ x ]
1744             for ( int cnt = 0; cnt < SIZE; cnt++ ){ workspace[ cnt ] = 0.0; }
1745             for ( int w = 0; w < nact_w; w++ ){
1746                double f_wc = fock->get( IL, nocc_w + w, n_oa_w + c );
1747                int inc1 = 1;
1748                daxpy_( &SIZE, &f_wc, FEH[ IL ][ w ], &inc1, workspace, &inc1 );
1749             }
1750 
1751             for ( int Id = 0; Id < num_irreps; Id++ ){
1752 
1753                const int Icenter = Irreps::directProd( IL, Id );
1754                const int nvir_d = indices->getNVIRT( Id );
1755                const int LDA_E = SIZE * nvir_d;
1756 
1757                if ( IL < Id ){ // Ic < Id
1758                   for ( int Ii = 0; Ii < num_irreps; Ii++ ){
1759                      const int Ij = Irreps::directProd( Ii, Icenter );
1760                      if ( Ii < Ij ){
1761                         const int jump_E = jump[ IL + num_irreps * CHEMPS2_CASPT2_E_TRIPLET ] + SIZE * shift_E_nonactive( indices, Id, Ii, Ij, -1 );
1762                         const int jump_H = jump[ Icenter + num_irreps * CHEMPS2_CASPT2_H_TRIPLET ] + shift_H_nonactive( indices, Ii, Ij, IL, Id, -1 );
1763                         const int size_ij = indices->getNOCC( Ii ) * indices->getNOCC( Ij );
1764                         #pragma omp parallel for schedule(static)
1765                         for ( int d = 0; d < nvir_d; d++ ){
1766                            const int ptr_H = jump_H + size_ij * ( c + nvir_w * d );
1767                            const int ptr_E = jump_E + SIZE * d;
1768                            matmat( 'N', SIZE, size_ij, 1,    6.0, workspace, SIZE, vector + ptr_H, 1,     result + ptr_E, LDA_E );
1769                            matmat( 'T', 1,    size_ij, SIZE, 6.0, workspace, SIZE, vector + ptr_E, LDA_E, result + ptr_H, 1     );
1770                         }
1771                      }
1772                   }
1773                }
1774 
1775                if ( IL == Id ){ // Ic == Id == Ia == Ib --> Iik == Ijl
1776                   for ( int Iij = 0; Iij < num_irreps; Iij++ ){
1777                      const int jump_E = jump[ IL + num_irreps * CHEMPS2_CASPT2_E_TRIPLET ] + SIZE * shift_E_nonactive( indices, Id,  Iij, Iij, -1 );
1778                      const int jump_H = jump[ Icenter + num_irreps * CHEMPS2_CASPT2_H_TRIPLET ] + shift_H_nonactive( indices, Iij, Iij, IL,  Id, -1 );
1779                      const int size_ij = ( indices->getNOCC( Iij ) * ( indices->getNOCC( Iij ) - 1 ) ) / 2;
1780                      #pragma omp parallel for schedule(static)
1781                      for ( int d = 0; d < c; d++ ){
1782                         const int ptr_H = jump_H + size_ij * ( d + ( c * ( c - 1 ) ) / 2 );
1783                         const int ptr_E = jump_E + SIZE * d;
1784                         matmat( 'N', SIZE, size_ij, 1,    -6.0, workspace, SIZE, vector + ptr_H, 1,     result + ptr_E, LDA_E );
1785                         matmat( 'T', 1,    size_ij, SIZE, -6.0, workspace, SIZE, vector + ptr_E, LDA_E, result + ptr_H, 1     );
1786                      }
1787                      #pragma omp parallel for schedule(static)
1788                      for ( int d = c+1; d < nvir_d; d++ ){
1789                         const int ptr_H = jump_H + size_ij * ( c + ( d * ( d - 1 ) ) / 2 );
1790                         const int ptr_E = jump_E + SIZE * d;
1791                         matmat( 'N', SIZE, size_ij, 1,    6.0, workspace, SIZE, vector + ptr_H, 1,     result + ptr_E, LDA_E );
1792                         matmat( 'T', 1,    size_ij, SIZE, 6.0, workspace, SIZE, vector + ptr_E, LDA_E, result + ptr_H, 1     );
1793                      }
1794                   }
1795                }
1796 
1797                if ( IL > Id ){ // Ic > Id
1798                   for ( int Ii = 0; Ii < num_irreps; Ii++ ){
1799                      const int Ij = Irreps::directProd( Ii, Icenter );
1800                      if ( Ii < Ij ){
1801                         const int jump_E = jump[ IL + num_irreps * CHEMPS2_CASPT2_E_TRIPLET ] + SIZE * shift_E_nonactive( indices, Id, Ii, Ij, -1 );
1802                         const int jump_H = jump[ Icenter + num_irreps * CHEMPS2_CASPT2_H_TRIPLET ] + shift_H_nonactive( indices, Ii, Ij, Id, IL, -1 );
1803                         const int size_ij = indices->getNOCC( Ii ) * indices->getNOCC( Ij );
1804                         #pragma omp parallel for schedule(static)
1805                         for ( int d = 0; d < nvir_d; d++ ){
1806                            const int ptr_H = jump_H + size_ij * ( d + nvir_d * c );
1807                            const int ptr_E = jump_E + SIZE * d;
1808                            matmat( 'N', SIZE, size_ij, 1,    -6.0, workspace, SIZE, vector + ptr_H, 1,     result + ptr_E, LDA_E );
1809                            matmat( 'T', 1,    size_ij, SIZE, -6.0, workspace, SIZE, vector + ptr_E, LDA_E, result + ptr_H, 1     );
1810                         }
1811                      }
1812                   }
1813                }
1814             }
1815          }
1816       }
1817    }
1818 
1819    // FGH singlet: < SG_cldx E_kw SH_aibj > = 2 delta_ac delta_bd ( delta_il delta_jk + delta_ik delta_jl ) / sqrt( 1 + delta_ij ) FGH[ Ix ][ w ][ x ]
1820    for ( int IL = 0; IL < num_irreps; IL++ ){ // IL == Ixw == Ik == Iac x Ibd x Il
1821       int SIZE = size_G[ IL ];
1822       const int nocc_w = indices->getNOCC( IL );
1823       const int nact_w = indices->getNDMRG( IL );
1824       if ( SIZE * nocc_w * nact_w > 0 ){
1825          for ( int k = 0; k < nocc_w; k++ ){
1826 
1827             // workspace_k[ x ] = sum_w f_kw FGH[ Ixw == IL == Ik ][ w ][ x ]
1828             for ( int cnt = 0; cnt < SIZE; cnt++ ){ workspace[ cnt ] = 0.0; }
1829             for ( int w = 0; w < nact_w; w++ ){
1830                double f_kw = fock->get( IL, k, nocc_w + w );
1831                int inc1 = 1;
1832                daxpy_( &SIZE, &f_kw, FGH[ IL ][ w ], &inc1, workspace, &inc1 );
1833             }
1834 
1835             for ( int Il = 0; Il < num_irreps; Il++ ){
1836 
1837                const int Icenter = Irreps::directProd( IL, Il );
1838                const int nocc_l = indices->getNOCC( Il );
1839 
1840                if ( IL < Il ){ // Ik < Il
1841                   for ( int Ia = 0; Ia < num_irreps; Ia++ ){
1842                      const int Ib = Irreps::directProd( Ia, Icenter );
1843                      if ( Ia < Ib ){
1844                         const int colsize = nocc_l * indices->getNVIRT( Ia ) * indices->getNVIRT( Ib );
1845                         if ( colsize > 0 ){
1846                            const int jump_G = jump[ IL + num_irreps * CHEMPS2_CASPT2_G_SINGLET ] + SIZE * shift_G_nonactive( indices, Il, Ia, Ib, +1 );
1847                            const int jump_H = jump[ Icenter + num_irreps * CHEMPS2_CASPT2_H_SINGLET ] + shift_H_nonactive( indices, IL, Il, Ia, Ib, +1 ) + k;
1848                            matmat( 'N', SIZE, colsize, 1,    2.0, workspace, SIZE, vector + jump_H, nocc_w, result + jump_G, SIZE   );
1849                            matmat( 'T', 1,    colsize, SIZE, 2.0, workspace, SIZE, vector + jump_G, SIZE,   result + jump_H, nocc_w );
1850                         }
1851                      }
1852                   }
1853                }
1854 
1855                if ( IL == Il ){ // Ik == Il == Ii == Ij --> Iac == Ibd   and   nocc_l == nocc_w
1856                   for ( int Iab = 0; Iab < num_irreps; Iab++ ){
1857                      const int jump_G = jump[ IL + num_irreps * CHEMPS2_CASPT2_G_SINGLET ] + SIZE * shift_G_nonactive( indices, Il, Iab, Iab, +1 );
1858                      const int jump_H = jump[ Icenter + num_irreps * CHEMPS2_CASPT2_H_SINGLET ] + shift_H_nonactive( indices, IL, Il, Iab, Iab, +1 );
1859                      const int size_ij = ( nocc_w * ( nocc_w + 1 ) ) / 2;
1860                      const int size_ab = ( indices->getNVIRT( Iab ) * ( indices->getNVIRT( Iab ) + 1 ) ) / 2;
1861                      const int LDA_G = SIZE * nocc_l;
1862                      #pragma omp parallel for schedule(static)
1863                      for ( int l = 0; l < k; l++ ){
1864                         const int ptr_H = jump_H + ( l + ( k * ( k + 1 ) ) / 2 );
1865                         const int ptr_G = jump_G + SIZE * l;
1866                         matmat( 'N', SIZE, size_ab, 1,    2.0, workspace, SIZE, vector + ptr_H, size_ij, result + ptr_G, LDA_G   );
1867                         matmat( 'T', 1,    size_ab, SIZE, 2.0, workspace, SIZE, vector + ptr_G, LDA_G,   result + ptr_H, size_ij );
1868                      }
1869                      #pragma omp parallel for schedule(static)
1870                      for ( int l = k; l < nocc_l; l++ ){
1871                         const int ptr_H = jump_H + ( k + ( l * ( l + 1 ) ) / 2 );
1872                         const int ptr_G = jump_G + SIZE * l;
1873                         const double factor = 2 * (( k == l ) ? SQRT2 : 1.0 );
1874                         matmat( 'N', SIZE, size_ab, 1,    factor, workspace, SIZE, vector + ptr_H, size_ij, result + ptr_G, LDA_G   );
1875                         matmat( 'T', 1,    size_ab, SIZE, factor, workspace, SIZE, vector + ptr_G, LDA_G,   result + ptr_H, size_ij );
1876                      }
1877                   }
1878                }
1879 
1880                if ( IL > Il ){ // Ik > Il
1881                   for ( int Ia = 0; Ia < num_irreps; Ia++ ){
1882                      const int Ib = Irreps::directProd( Ia, Icenter );
1883                      if ( Ia < Ib ){
1884                         const int jump_G = jump[ IL + num_irreps * CHEMPS2_CASPT2_G_SINGLET ] + SIZE * shift_G_nonactive( indices, Il, Ia, Ib, +1 );
1885                         const int jump_H = jump[ Icenter + num_irreps * CHEMPS2_CASPT2_H_SINGLET ] + shift_H_nonactive( indices, Il, IL, Ia, Ib, +1 );
1886                         const int size_ab = indices->getNVIRT( Ia ) * indices->getNVIRT( Ib );
1887                         #pragma omp parallel for schedule(static)
1888                         for ( int ab = 0; ab < size_ab; ab++ ){
1889                            const int ptr_H = jump_H + nocc_l * ( k + nocc_w * ab );
1890                            const int ptr_G = jump_G + SIZE * nocc_l * ab;
1891                            matmat( 'N', SIZE, nocc_l, 1,    2.0, workspace, SIZE, vector + ptr_H, 1,    result + ptr_G, SIZE );
1892                            matmat( 'T', 1,    nocc_l, SIZE, 2.0, workspace, SIZE, vector + ptr_G, SIZE, result + ptr_H, 1    );
1893                         }
1894                      }
1895                   }
1896                }
1897             }
1898          }
1899       }
1900    }
1901 
1902    // FGH triplet: < TG_cldx E_kw TH_aibj > = 6 delta_ac delta_bd ( delta_il delta_jk - delta_ik delta_jl ) / sqrt( 1 + delta_ij ) FGH[ Ix ][ w ][ x ]
1903    for ( int IL = 0; IL < num_irreps; IL++ ){ // IL == Ixw == Ik == Iac x Ibd x Il
1904       int SIZE = size_G[ IL ];
1905       const int nocc_w = indices->getNOCC( IL );
1906       const int nact_w = indices->getNDMRG( IL );
1907       if ( SIZE * nocc_w * nact_w > 0 ){
1908          for ( int k = 0; k < nocc_w; k++ ){
1909 
1910             // workspace_k[ x ] = sum_w f_kw FGH[ Ixw == IL == Ik ][ w ][ x ]
1911             for ( int cnt = 0; cnt < SIZE; cnt++ ){ workspace[ cnt ] = 0.0; }
1912             for ( int w = 0; w < nact_w; w++ ){
1913                double f_kw = fock->get( IL, k, nocc_w + w );
1914                int inc1 = 1;
1915                daxpy_( &SIZE, &f_kw, FGH[ IL ][ w ], &inc1, workspace, &inc1 );
1916             }
1917 
1918             for ( int Il = 0; Il < num_irreps; Il++ ){
1919 
1920                const int Icenter = Irreps::directProd( IL, Il );
1921                const int nocc_l = indices->getNOCC( Il );
1922 
1923                if ( IL < Il ){ // Ik < Il
1924                   for ( int Ia = 0; Ia < num_irreps; Ia++ ){
1925                      const int Ib = Irreps::directProd( Ia, Icenter );
1926                      if ( Ia < Ib ){
1927                         const int colsize = nocc_l * indices->getNVIRT( Ia ) * indices->getNVIRT( Ib );
1928                         if ( colsize > 0 ){
1929                            const int jump_G = jump[ IL + num_irreps * CHEMPS2_CASPT2_G_TRIPLET ] + SIZE * shift_G_nonactive( indices, Il, Ia, Ib, -1 );
1930                            const int jump_H = jump[ Icenter + num_irreps * CHEMPS2_CASPT2_H_TRIPLET ] + shift_H_nonactive( indices, IL, Il, Ia, Ib, -1 ) + k;
1931                            matmat( 'N', SIZE, colsize, 1,    -6.0, workspace, SIZE, vector + jump_H, nocc_w, result + jump_G, SIZE   );
1932                            matmat( 'T', 1,    colsize, SIZE, -6.0, workspace, SIZE, vector + jump_G, SIZE,   result + jump_H, nocc_w );
1933                         }
1934                      }
1935                   }
1936                }
1937 
1938                if ( IL == Il ){ // Ik == Il == Ii == Ij --> Iac == Ibd   and   nocc_l == nocc_w
1939                   for ( int Iab = 0; Iab < num_irreps; Iab++ ){
1940                      const int jump_G = jump[ IL + num_irreps * CHEMPS2_CASPT2_G_TRIPLET ] + SIZE * shift_G_nonactive( indices, Il, Iab, Iab, -1 );
1941                      const int jump_H = jump[ Icenter + num_irreps * CHEMPS2_CASPT2_H_TRIPLET ] + shift_H_nonactive( indices, IL, Il, Iab, Iab, -1 );
1942                      const int size_ij = ( nocc_w * ( nocc_w - 1 ) ) / 2;
1943                      const int size_ab = ( indices->getNVIRT( Iab ) * ( indices->getNVIRT( Iab ) - 1 ) ) / 2;
1944                      const int LDA_G = SIZE * nocc_l;
1945                      #pragma omp parallel for schedule(static)
1946                      for ( int l = 0; l < k; l++ ){
1947                         const int ptr_H = jump_H + ( l + ( k * ( k - 1 ) ) / 2 );
1948                         const int ptr_G = jump_G + SIZE * l;
1949                         matmat( 'N', SIZE, size_ab, 1,    6.0, workspace, SIZE, vector + ptr_H, size_ij, result + ptr_G, LDA_G   );
1950                         matmat( 'T', 1,    size_ab, SIZE, 6.0, workspace, SIZE, vector + ptr_G, LDA_G,   result + ptr_H, size_ij );
1951                      }
1952                      #pragma omp parallel for schedule(static)
1953                      for ( int l = k+1; l < nocc_l; l++ ){
1954                         const int ptr_H = jump_H + ( k + ( l * ( l - 1 ) ) / 2 );
1955                         const int ptr_G = jump_G + SIZE * l;
1956                         matmat( 'N', SIZE, size_ab, 1,    -6.0, workspace, SIZE, vector + ptr_H, size_ij, result + ptr_G, LDA_G   );
1957                         matmat( 'T', 1,    size_ab, SIZE, -6.0, workspace, SIZE, vector + ptr_G, LDA_G,   result + ptr_H, size_ij );
1958                      }
1959                   }
1960                }
1961 
1962                if ( IL > Il ){ // Ik > Il
1963                   for ( int Ia = 0; Ia < num_irreps; Ia++ ){
1964                      const int Ib = Irreps::directProd( Ia, Icenter );
1965                      if ( Ia < Ib ){
1966                         const int jump_G = jump[ IL + num_irreps * CHEMPS2_CASPT2_G_TRIPLET ] + SIZE * shift_G_nonactive( indices, Il, Ia, Ib, -1 );
1967                         const int jump_H = jump[ Icenter + num_irreps * CHEMPS2_CASPT2_H_TRIPLET ] + shift_H_nonactive( indices, Il, IL, Ia, Ib, -1 );
1968                         const int size_ab = indices->getNVIRT( Ia ) * indices->getNVIRT( Ib );
1969                         #pragma omp parallel for schedule(static)
1970                         for ( int ab = 0; ab < size_ab; ab++ ){
1971                            const int ptr_H = jump_H + nocc_l * ( k + nocc_w * ab );
1972                            const int ptr_G = jump_G + SIZE * nocc_l * ab;
1973                            matmat( 'N', SIZE, nocc_l, 1,    6.0, workspace, SIZE, vector + ptr_H, 1,    result + ptr_G, SIZE );
1974                            matmat( 'T', 1,    nocc_l, SIZE, 6.0, workspace, SIZE, vector + ptr_G, SIZE, result + ptr_H, 1    );
1975                         }
1976                      }
1977                   }
1978                }
1979             }
1980          }
1981       }
1982    }
1983 
1984    // FDE singlet: < D(blxy) E_kw SE_tiaj > = 1 delta_ab ( delta_ik delta_jl + delta_il delta_jk ) / sqrt( 1 + delta_ij ) FDE_singlet[ Ib x Il ][ It ][ w ][ xy, t ]
1985    for ( int IL = 0; IL < num_irreps; IL++ ){ // IL == Ix x Iy == Ib x Il
1986       const int SIZE_L = size_D[ IL ];
1987       for ( int IR = 0; IR < num_irreps; IR++ ){ // IR = It = Ia x Ii x Ij
1988          const int SIZE_R = size_E[ IR ];
1989          const int Ikw = Irreps::directProd( IL, IR ); // Ikw == Ik == Iw
1990          const int nocc_kw = indices->getNOCC( Ikw );
1991          const int nact_kw = indices->getNDMRG( Ikw );
1992          int total_size = SIZE_L * SIZE_R;
1993          if ( total_size * nocc_kw * nact_kw > 0 ){
1994             for ( int k = 0; k < nocc_kw; k++ ){
1995                for ( int cnt = 0; cnt < total_size; cnt++ ){ workspace[ cnt ] = 0.0; }
1996                for ( int w = 0; w < nact_kw; w++ ){
1997                   double f_kw = fock->get( Ikw, k, nocc_kw + w );
1998                   int inc1 = 1;
1999                   daxpy_( &total_size, &f_kw, FDE_singlet[ IL ][ IR ][ w ], &inc1, workspace, &inc1 );
2000                }
2001                for ( int Iab = 0; Iab < num_irreps; Iab++ ){
2002                   const int nvir_ab = indices->getNVIRT( Iab );
2003                   const int Il = Irreps::directProd( Iab, IL );
2004                   const int nocc_l = indices->getNOCC( Il );
2005                   if ( nvir_ab * nocc_l > 0 ){
2006                      const int jump_D = jump[ IL + num_irreps * CHEMPS2_CASPT2_D         ] + SIZE_L * shift_D_nonactive( indices, Il, Iab );
2007                      const int jump_E = jump[ IR + num_irreps * CHEMPS2_CASPT2_E_SINGLET ] + SIZE_R * (( Ikw <= Il ) ? shift_E_nonactive( indices, Iab, Ikw, Il,  +1 )
2008                                                                                                                      : shift_E_nonactive( indices, Iab, Il,  Ikw, +1 ));
2009                      if ( Ikw == Il ){ // irrep_k == irrep_l
2010                         #pragma omp parallel for schedule(static)
2011                         for ( int l = 0; l < nocc_l; l++ ){
2012                            const double factor = (( k == l ) ? SQRT2 : 1.0 );
2013                            const int ptr_L = jump_D + SIZE_L * l;
2014                            const int ptr_R = jump_E + SIZE_R * nvir_ab * (( k < l ) ? ( k + ( l * ( l + 1 ) ) / 2 ) : ( l + ( k * ( k + 1 ) ) / 2 ));
2015                            const int LDA_L = SIZE_L * nocc_l;
2016                            matmat( 'N', SIZE_L, nvir_ab, SIZE_R, factor, workspace, SIZE_L, vector + ptr_R, SIZE_R, result + ptr_L, LDA_L  );
2017                            matmat( 'T', SIZE_R, nvir_ab, SIZE_L, factor, workspace, SIZE_L, vector + ptr_L, LDA_L,  result + ptr_R, SIZE_R );
2018                         }
2019                      } else { // irrep_k != irrep_l
2020                         #pragma omp parallel for schedule(static)
2021                         for ( int l = 0; l < nocc_l; l++ ){
2022                            const int ptr_L = jump_D + SIZE_L * l;
2023                            const int ptr_R = jump_E + SIZE_R * nvir_ab * (( Ikw < Il ) ? ( k + nocc_kw * l ) : ( l + nocc_l * k ));
2024                            const int LDA_L = SIZE_L * nocc_l;
2025                            matmat( 'N', SIZE_L, nvir_ab, SIZE_R, 1.0, workspace, SIZE_L, vector + ptr_R, SIZE_R, result + ptr_L, LDA_L  );
2026                            matmat( 'T', SIZE_R, nvir_ab, SIZE_L, 1.0, workspace, SIZE_L, vector + ptr_L, LDA_L,  result + ptr_R, SIZE_R );
2027                         }
2028                      }
2029                   }
2030                }
2031             }
2032          }
2033       }
2034    }
2035 
2036    // FDE triplet: < D(blxy) E_kw TE_tiaj > = 3 delta_ab ( delta_ik delta_jl - delta_il delta_jk ) / sqrt( 1 + delta_ij ) FDE_triplet[ Ib x Il ][ It ][ w ][ xy, t ]
2037    for ( int IL = 0; IL < num_irreps; IL++ ){ // IL == Ix x Iy == Ib x Il
2038       const int SIZE_L = size_D[ IL ];
2039       for ( int IR = 0; IR < num_irreps; IR++ ){ // IR = It = Ia x Ii x Ij
2040          const int SIZE_R = size_E[ IR ];
2041          const int Ikw = Irreps::directProd( IL, IR ); // Ikw == Ik == Iw
2042          const int nocc_kw = indices->getNOCC( Ikw );
2043          const int nact_kw = indices->getNDMRG( Ikw );
2044          int total_size = SIZE_L * SIZE_R;
2045          if ( total_size * nocc_kw * nact_kw > 0 ){
2046             for ( int k = 0; k < nocc_kw; k++ ){
2047                for ( int cnt = 0; cnt < total_size; cnt++ ){ workspace[ cnt ] = 0.0; }
2048                for ( int w = 0; w < nact_kw; w++ ){
2049                   double f_kw = fock->get( Ikw, k, nocc_kw + w );
2050                   int inc1 = 1;
2051                   daxpy_( &total_size, &f_kw, FDE_triplet[ IL ][ IR ][ w ], &inc1, workspace, &inc1 );
2052                }
2053                for ( int Iab = 0; Iab < num_irreps; Iab++ ){
2054                   const int nvir_ab = indices->getNVIRT( Iab );
2055                   const int Il = Irreps::directProd( Iab, IL );
2056                   const int nocc_l = indices->getNOCC( Il );
2057                   if ( nvir_ab * nocc_l > 0 ){
2058                      const int jump_D = jump[ IL + num_irreps * CHEMPS2_CASPT2_D         ] + SIZE_L * shift_D_nonactive( indices, Il, Iab );
2059                      const int jump_E = jump[ IR + num_irreps * CHEMPS2_CASPT2_E_TRIPLET ] + SIZE_R * (( Ikw <= Il ) ? shift_E_nonactive( indices, Iab, Ikw, Il,  -1 )
2060                                                                                                                      : shift_E_nonactive( indices, Iab, Il,  Ikw, -1 ));
2061                      if ( Ikw == Il ){ // irrep_k == irrep_l
2062                         #pragma omp parallel for schedule(static)
2063                         for ( int l = 0; l < k; l++ ){
2064                            const int ptr_L = jump_D + SIZE_L * l;
2065                            const int ptr_R = jump_E + SIZE_R * nvir_ab * ( l + ( k * ( k - 1 ) ) / 2 );
2066                            const int LDA_L = SIZE_L * nocc_l;
2067                            matmat( 'N', SIZE_L, nvir_ab, SIZE_R, -3.0, workspace, SIZE_L, vector + ptr_R, SIZE_R, result + ptr_L, LDA_L  );
2068                            matmat( 'T', SIZE_R, nvir_ab, SIZE_L, -3.0, workspace, SIZE_L, vector + ptr_L, LDA_L,  result + ptr_R, SIZE_R );
2069                         }
2070                         #pragma omp parallel for schedule(static)
2071                         for ( int l = k+1; l < nocc_l; l++ ){
2072                            const int ptr_L = jump_D + SIZE_L * l;
2073                            const int ptr_R = jump_E + SIZE_R * nvir_ab * ( k + ( l * ( l - 1 ) ) / 2 );
2074                            const int LDA_L = SIZE_L * nocc_l;
2075                            matmat( 'N', SIZE_L, nvir_ab, SIZE_R, 3.0, workspace, SIZE_L, vector + ptr_R, SIZE_R, result + ptr_L, LDA_L  );
2076                            matmat( 'T', SIZE_R, nvir_ab, SIZE_L, 3.0, workspace, SIZE_L, vector + ptr_L, LDA_L,  result + ptr_R, SIZE_R );
2077                         }
2078                      } else { // irrep_k != irrep_l
2079                         #pragma omp parallel for schedule(static)
2080                         for ( int l = 0; l < nocc_l; l++ ){
2081                            const double factor = (( Ikw < Il ) ? 3.0 : -3.0 );
2082                            const int ptr_L = jump_D + SIZE_L * l;
2083                            const int ptr_R = jump_E + SIZE_R * nvir_ab * (( Ikw < Il ) ? ( k + nocc_kw * l ) : ( l + nocc_l * k ));
2084                            const int LDA_L = SIZE_L * nocc_l;
2085                            matmat( 'N', SIZE_L, nvir_ab, SIZE_R, factor, workspace, SIZE_L, vector + ptr_R, SIZE_R, result + ptr_L, LDA_L  );
2086                            matmat( 'T', SIZE_R, nvir_ab, SIZE_L, factor, workspace, SIZE_L, vector + ptr_L, LDA_L,  result + ptr_R, SIZE_R );
2087                         }
2088                      }
2089                   }
2090                }
2091             }
2092          }
2093       }
2094    }
2095 
2096    // FDG singlet: < D(djxy) E_wc SG_aibt > = 1 delta_ij ( delta_ac delta_bd + delta_ad delta_bc ) / sqrt( 1 + delta_ab ) FDG_singlet[ Ij x Id ][ It ][ w ][ xy, t ]
2097    for ( int IL = 0; IL < num_irreps; IL++ ){ // IL == Ix x Iy == Id x Ij
2098       const int SIZE_L = size_D[ IL ];
2099       for ( int IR = 0; IR < num_irreps; IR++ ){ // IR = It = Ia x Ib x Ii
2100          const int SIZE_R = size_E[ IR ];
2101          const int Iwc = Irreps::directProd( IL, IR ); // Iwc == Ic == Iw
2102          const int nocc_wc = indices->getNOCC( Iwc );
2103          const int nact_wc = indices->getNDMRG( Iwc );
2104          const int n_oa_wc = nocc_wc + nact_wc;
2105          const int nvir_wc = indices->getNVIRT( Iwc );
2106          int total_size = SIZE_L * SIZE_R;
2107          if ( total_size * nvir_wc * nact_wc > 0 ){
2108             for ( int c = 0; c < nvir_wc; c++ ){
2109                for ( int cnt = 0; cnt < total_size; cnt++ ){ workspace[ cnt ] = 0.0; }
2110                for ( int w = 0; w < nact_wc; w++ ){
2111                   double f_wc = fock->get( Iwc, nocc_wc + w, n_oa_wc + c );
2112                   int inc1 = 1;
2113                   daxpy_( &total_size, &f_wc, FDG_singlet[ IL ][ IR ][ w ], &inc1, workspace, &inc1 );
2114                }
2115                for ( int Iij = 0; Iij < num_irreps; Iij++ ){
2116                   const int nocc_ij = indices->getNOCC( Iij );
2117                   const int Id = Irreps::directProd( Iij, IL );
2118                   const int nvir_d = indices->getNVIRT( Id );
2119                   if ( nvir_d * nocc_ij > 0 ){
2120                      const int jump_D = jump[ IL + num_irreps * CHEMPS2_CASPT2_D ] + SIZE_L * shift_D_nonactive( indices, Iij, Id );
2121                      const int jump_G = jump[ IR + num_irreps * CHEMPS2_CASPT2_G_SINGLET ] + SIZE_R * (( Iwc <= Id ) ? shift_G_nonactive( indices, Iij, Iwc, Id,  +1 )
2122                                                                                                                      : shift_G_nonactive( indices, Iij, Id,  Iwc, +1 ));
2123                      if ( Iwc == Id ){ // irrep_c == irrep_d
2124                         if ( c > 0 ){
2125                            const int ptr_L = jump_D;
2126                            const int ptr_R = jump_G + SIZE_R * nocc_ij * ( c * ( c + 1 ) ) / 2;
2127                            matmat( 'N', SIZE_L, c * nocc_ij, SIZE_R, 1.0, workspace, SIZE_L, vector + ptr_R, SIZE_R, result + ptr_L, SIZE_L );
2128                            matmat( 'T', SIZE_R, c * nocc_ij, SIZE_L, 1.0, workspace, SIZE_L, vector + ptr_L, SIZE_L, result + ptr_R, SIZE_R );
2129                         }
2130                         #pragma omp parallel for schedule(static)
2131                         for ( int d = c; d < nvir_d; d++ ){
2132                            const double factor = (( c == d ) ? SQRT2 : 1.0 );
2133                            const int ptr_L = jump_D + SIZE_L * nocc_ij * d;
2134                            const int ptr_R = jump_G + SIZE_R * nocc_ij * ( c + ( d * ( d + 1 ) ) / 2 );
2135                            matmat( 'N', SIZE_L, nocc_ij, SIZE_R, factor, workspace, SIZE_L, vector + ptr_R, SIZE_R, result + ptr_L, SIZE_L );
2136                            matmat( 'T', SIZE_R, nocc_ij, SIZE_L, factor, workspace, SIZE_L, vector + ptr_L, SIZE_L, result + ptr_R, SIZE_R );
2137                         }
2138                      } else { // irrep_c != irrep_d
2139                         if ( Iwc < Id ){
2140                            #pragma omp parallel for schedule(static)
2141                            for ( int d = 0; d < nvir_d; d++ ){
2142                               const int ptr_L = jump_D + SIZE_L * nocc_ij * d;
2143                               const int ptr_R = jump_G + SIZE_R * nocc_ij * ( c + nvir_wc * d );
2144                               matmat( 'N', SIZE_L, nocc_ij, SIZE_R, 1.0, workspace, SIZE_L, vector + ptr_R, SIZE_R, result + ptr_L, SIZE_L );
2145                               matmat( 'T', SIZE_R, nocc_ij, SIZE_L, 1.0, workspace, SIZE_L, vector + ptr_L, SIZE_L, result + ptr_R, SIZE_R );
2146                            }
2147                         } else {
2148                            const int ptr_L = jump_D;
2149                            const int ptr_R = jump_G + SIZE_R * nocc_ij * nvir_d * c;
2150                            matmat( 'N', SIZE_L, nvir_d * nocc_ij, SIZE_R, 1.0, workspace, SIZE_L, vector + ptr_R, SIZE_R, result + ptr_L, SIZE_L );
2151                            matmat( 'T', SIZE_R, nvir_d * nocc_ij, SIZE_L, 1.0, workspace, SIZE_L, vector + ptr_L, SIZE_L, result + ptr_R, SIZE_R );
2152                         }
2153                      }
2154                   }
2155                }
2156             }
2157          }
2158       }
2159    }
2160 
2161    // FDG triplet: < D(djxy) E_wc TG_aibt > = 3 delta_ij ( delta_ac delta_bd - delta_ad delta_bc ) / sqrt( 1 + delta_ab ) FDG_triplet[ Ij x Id ][ It ][ w ][ xy, t ]
2162    for ( int IL = 0; IL < num_irreps; IL++ ){ // IL == Ix x Iy == Id x Ij
2163       const int SIZE_L = size_D[ IL ];
2164       for ( int IR = 0; IR < num_irreps; IR++ ){ // IR = It = Ia x Ib x Ii
2165          const int SIZE_R = size_E[ IR ];
2166          const int Iwc = Irreps::directProd( IL, IR ); // Iwc == Ic == Iw
2167          const int nocc_wc = indices->getNOCC( Iwc );
2168          const int nact_wc = indices->getNDMRG( Iwc );
2169          const int n_oa_wc = nocc_wc + nact_wc;
2170          const int nvir_wc = indices->getNVIRT( Iwc );
2171          int total_size = SIZE_L * SIZE_R;
2172          if ( total_size * nvir_wc * nact_wc > 0 ){
2173             for ( int c = 0; c < nvir_wc; c++ ){
2174                for ( int cnt = 0; cnt < total_size; cnt++ ){ workspace[ cnt ] = 0.0; }
2175                for ( int w = 0; w < nact_wc; w++ ){
2176                   double f_wc = fock->get( Iwc, nocc_wc + w, n_oa_wc + c );
2177                   int inc1 = 1;
2178                   daxpy_( &total_size, &f_wc, FDG_triplet[ IL ][ IR ][ w ], &inc1, workspace, &inc1 );
2179                }
2180                for ( int Iij = 0; Iij < num_irreps; Iij++ ){
2181                   const int nocc_ij = indices->getNOCC( Iij );
2182                   const int Id = Irreps::directProd( Iij, IL );
2183                   const int nvir_d = indices->getNVIRT( Id );
2184                   if ( nvir_d * nocc_ij > 0 ){
2185                      const int jump_D = jump[ IL + num_irreps * CHEMPS2_CASPT2_D ] + SIZE_L * shift_D_nonactive( indices, Iij, Id );
2186                      const int jump_G = jump[ IR + num_irreps * CHEMPS2_CASPT2_G_TRIPLET ] + SIZE_R * (( Iwc <= Id ) ? shift_G_nonactive( indices, Iij, Iwc, Id,  -1 )
2187                                                                                                                      : shift_G_nonactive( indices, Iij, Id,  Iwc, -1 ));
2188                      if ( Iwc == Id ){ // irrep_c == irrep_d
2189                         if ( c > 0 ){
2190                            const int ptr_L = jump_D;
2191                            const int ptr_R = jump_G + SIZE_R * nocc_ij * ( c * ( c - 1 ) ) / 2;
2192                            matmat( 'N', SIZE_L, c * nocc_ij, SIZE_R, -3.0, workspace, SIZE_L, vector + ptr_R, SIZE_R, result + ptr_L, SIZE_L );
2193                            matmat( 'T', SIZE_R, c * nocc_ij, SIZE_L, -3.0, workspace, SIZE_L, vector + ptr_L, SIZE_L, result + ptr_R, SIZE_R );
2194                         }
2195                         #pragma omp parallel for schedule(static)
2196                         for ( int d = c+1; d < nvir_d; d++ ){
2197                            const int ptr_L = jump_D + SIZE_L * nocc_ij * d;
2198                            const int ptr_R = jump_G + SIZE_R * nocc_ij * ( c + ( d * ( d - 1 ) ) / 2 );
2199                            matmat( 'N', SIZE_L, nocc_ij, SIZE_R, 3.0, workspace, SIZE_L, vector + ptr_R, SIZE_R, result + ptr_L, SIZE_L );
2200                            matmat( 'T', SIZE_R, nocc_ij, SIZE_L, 3.0, workspace, SIZE_L, vector + ptr_L, SIZE_L, result + ptr_R, SIZE_R );
2201                         }
2202                      } else { // irrep_c != irrep_d
2203                         if ( Iwc < Id ){
2204                            #pragma omp parallel for schedule(static)
2205                            for ( int d = 0; d < nvir_d; d++ ){
2206                               const int ptr_L = jump_D + SIZE_L * nocc_ij * d;
2207                               const int ptr_R = jump_G + SIZE_R * nocc_ij * ( c + nvir_wc * d );
2208                               matmat( 'N', SIZE_L, nocc_ij, SIZE_R, 3.0, workspace, SIZE_L, vector + ptr_R, SIZE_R, result + ptr_L, SIZE_L );
2209                               matmat( 'T', SIZE_R, nocc_ij, SIZE_L, 3.0, workspace, SIZE_L, vector + ptr_L, SIZE_L, result + ptr_R, SIZE_R );
2210                            }
2211                         } else {
2212                            const int ptr_L = jump_D;
2213                            const int ptr_R = jump_G + SIZE_R * nocc_ij * nvir_d * c;
2214                            matmat( 'N', SIZE_L, nvir_d * nocc_ij, SIZE_R, -3.0, workspace, SIZE_L, vector + ptr_R, SIZE_R, result + ptr_L, SIZE_L );
2215                            matmat( 'T', SIZE_R, nvir_d * nocc_ij, SIZE_L, -3.0, workspace, SIZE_L, vector + ptr_L, SIZE_L, result + ptr_R, SIZE_R );
2216                         }
2217                      }
2218                   }
2219                }
2220             }
2221          }
2222       }
2223    }
2224 
2225    delete [] workspace;
2226 
2227 }
2228 
shift_H_nonactive(const DMRGSCFindices * idx,const int irrep_i,const int irrep_j,const int irrep_a,const int irrep_b,const int ST)2229 int CheMPS2::CASPT2::shift_H_nonactive( const DMRGSCFindices * idx, const int irrep_i, const int irrep_j, const int irrep_a, const int irrep_b, const int ST ){
2230 
2231    assert( irrep_i <= irrep_j );
2232    assert( irrep_a <= irrep_b );
2233    const int irrep_prod = Irreps::directProd( irrep_i, irrep_j );
2234    const int irrep_virt = Irreps::directProd( irrep_a, irrep_b );
2235    assert( irrep_prod == irrep_virt );
2236    const int n_irreps = idx->getNirreps();
2237 
2238    int shift = 0;
2239    if ( irrep_prod == 0 ){
2240       for ( int Iij = 0; Iij < n_irreps; Iij++ ){
2241          for ( int Iab = 0; Iab < n_irreps; Iab++ ){
2242             if (( irrep_i == Iij ) && ( irrep_j == Iij ) && ( irrep_a == Iab ) && ( irrep_b == Iab )){
2243                Iij = n_irreps;
2244                Iab = n_irreps;
2245             } else {
2246                shift += ( idx->getNOCC( Iij ) * ( idx->getNOCC( Iij ) + ST ) * idx->getNVIRT( Iab ) * ( idx->getNVIRT( Iab ) + ST ) ) / 4;
2247             }
2248          }
2249       }
2250    } else {
2251       for ( int Ii = 0; Ii < n_irreps; Ii++ ){
2252          const int Ij = Irreps::directProd( irrep_prod, Ii );
2253          if ( Ii < Ij ){
2254             for ( int Ia = 0; Ia < n_irreps; Ia++ ){
2255                const int Ib = Irreps::directProd( irrep_prod, Ia );
2256                if ( Ia < Ib ){
2257                   if (( irrep_i == Ii ) && ( irrep_j == Ij ) && ( irrep_a == Ia ) && ( irrep_b == Ib )){
2258                      Ii = n_irreps;
2259                      Ia = n_irreps;
2260                   } else {
2261                      shift += idx->getNOCC( Ii ) * idx->getNOCC( Ij ) * idx->getNVIRT( Ia ) * idx->getNVIRT( Ib );
2262                   }
2263                }
2264             }
2265          }
2266       }
2267    }
2268    return shift;
2269 
2270 }
2271 
shift_G_nonactive(const DMRGSCFindices * idx,const int irrep_i,const int irrep_a,const int irrep_b,const int ST)2272 int CheMPS2::CASPT2::shift_G_nonactive( const DMRGSCFindices * idx, const int irrep_i, const int irrep_a, const int irrep_b, const int ST ){
2273 
2274    assert( irrep_a <= irrep_b );
2275    const int irrep_virt = Irreps::directProd( irrep_a,    irrep_b );
2276    const int irrep_prod = Irreps::directProd( irrep_virt, irrep_i );
2277    const int n_irreps = idx->getNirreps();
2278 
2279    int shift = 0;
2280    for ( int Ii = 0; Ii < n_irreps; Ii++ ){
2281       const int Ivirt = Irreps::directProd( Ii, irrep_prod );
2282       if ( Ivirt == 0 ){
2283          for ( int Iab = 0; Iab < n_irreps; Iab++ ){
2284             if (( irrep_i == Ii ) && ( irrep_a == Iab ) && ( irrep_b == Iab )){
2285                Ii  = n_irreps;
2286                Iab = n_irreps;
2287             } else {
2288                shift += ( idx->getNOCC( Ii ) * idx->getNVIRT( Iab ) * ( idx->getNVIRT( Iab ) + ST ) ) / 2;
2289             }
2290          }
2291       } else {
2292          for ( int Ia = 0; Ia < n_irreps; Ia++ ){
2293             const int Ib = Irreps::directProd( Ivirt, Ia );
2294             if ( Ia < Ib ){
2295                if (( irrep_i == Ii ) && ( irrep_a == Ia ) && ( irrep_b == Ib )){
2296                   Ii = n_irreps;
2297                   Ia = n_irreps;
2298                } else {
2299                   shift += idx->getNOCC( Ii ) * idx->getNVIRT( Ia ) * idx->getNVIRT( Ib );
2300                }
2301             }
2302          }
2303       }
2304    }
2305    return shift;
2306 
2307 }
2308 
shift_E_nonactive(const DMRGSCFindices * idx,const int irrep_a,const int irrep_i,const int irrep_j,const int ST)2309 int CheMPS2::CASPT2::shift_E_nonactive( const DMRGSCFindices * idx, const int irrep_a, const int irrep_i, const int irrep_j, const int ST ){
2310 
2311    assert( irrep_i <= irrep_j );
2312    const int irrep_occ  = Irreps::directProd( irrep_i,   irrep_j );
2313    const int irrep_prod = Irreps::directProd( irrep_occ, irrep_a );
2314    const int n_irreps = idx->getNirreps();
2315 
2316    int shift = 0;
2317    for ( int Ia = 0; Ia < n_irreps; Ia++ ){
2318       const int Iocc = Irreps::directProd( Ia, irrep_prod );
2319       if ( Iocc == 0 ){
2320          for ( int Iij = 0; Iij < n_irreps; Iij++ ){
2321             if (( irrep_a == Ia ) && ( irrep_i == Iij ) && ( irrep_j == Iij )){
2322                Ia  = n_irreps;
2323                Iij = n_irreps;
2324             } else {
2325                shift += ( idx->getNVIRT( Ia ) * idx->getNOCC( Iij ) * ( idx->getNOCC( Iij ) + ST ) ) / 2;
2326             }
2327          }
2328       } else {
2329          for ( int Ii = 0; Ii < n_irreps; Ii++ ){
2330             const int Ij = Irreps::directProd( Iocc, Ii );
2331             if ( Ii < Ij ){
2332                if (( irrep_a == Ia ) && ( irrep_i == Ii ) && ( irrep_j == Ij )){
2333                   Ia = n_irreps;
2334                   Ii = n_irreps;
2335                } else {
2336                   shift += idx->getNVIRT( Ia ) * idx->getNOCC( Ii ) * idx->getNOCC( Ij );
2337                }
2338             }
2339          }
2340       }
2341    }
2342    return shift;
2343 
2344 }
2345 
shift_F_nonactive(const DMRGSCFindices * idx,const int irrep_a,const int irrep_b,const int ST)2346 int CheMPS2::CASPT2::shift_F_nonactive( const DMRGSCFindices * idx, const int irrep_a, const int irrep_b, const int ST ){
2347 
2348    assert( irrep_a <= irrep_b );
2349    const int irr_prod = Irreps::directProd( irrep_a, irrep_b );
2350    const int n_irreps = idx->getNirreps();
2351 
2352    int shift = 0;
2353    if ( irr_prod == 0 ){
2354       for ( int Iab = 0; Iab < n_irreps; Iab++ ){
2355          if (( irrep_a == Iab ) && ( irrep_b == Iab )){
2356             Iab = n_irreps;
2357          } else {
2358             shift += ( idx->getNVIRT( Iab ) * ( idx->getNVIRT( Iab ) + ST ) ) / 2;
2359          }
2360       }
2361    } else {
2362       for ( int Ia = 0; Ia < n_irreps; Ia++ ){
2363          const int Ib = Irreps::directProd( Ia, irr_prod );
2364          if ( Ia < Ib ){
2365             if (( irrep_a == Ia ) && ( irrep_b == Ib )){
2366                Ia = n_irreps;
2367             } else {
2368                shift += idx->getNVIRT( Ia ) * idx->getNVIRT( Ib );
2369             }
2370          }
2371       }
2372    }
2373    return shift;
2374 
2375 }
2376 
shift_B_nonactive(const DMRGSCFindices * idx,const int irrep_i,const int irrep_j,const int ST)2377 int CheMPS2::CASPT2::shift_B_nonactive( const DMRGSCFindices * idx, const int irrep_i, const int irrep_j, const int ST ){
2378 
2379    assert( irrep_i <= irrep_j );
2380    const int irr_prod = Irreps::directProd( irrep_i, irrep_j );
2381    const int n_irreps = idx->getNirreps();
2382 
2383    int shift = 0;
2384    if ( irr_prod == 0 ){
2385       for ( int Iij = 0; Iij < n_irreps; Iij++ ){
2386          if (( Iij == irrep_i ) && ( Iij == irrep_j )){
2387             Iij = n_irreps;
2388          } else {
2389             shift += ( idx->getNOCC( Iij ) * ( idx->getNOCC( Iij ) + ST ) ) / 2;
2390          }
2391       }
2392    } else {
2393       for ( int Ii = 0; Ii < n_irreps; Ii++ ){
2394          const int Ij = Irreps::directProd( irr_prod, Ii );
2395          if ( Ii < Ij ){
2396             if (( Ii == irrep_i ) && ( Ij == irrep_j )){
2397                Ii = n_irreps;
2398             } else {
2399                shift += idx->getNOCC( Ii ) * idx->getNOCC( Ij );
2400             }
2401          }
2402       }
2403    }
2404    return shift;
2405 
2406 }
2407 
shift_D_nonactive(const DMRGSCFindices * idx,const int irrep_i,const int irrep_a)2408 int CheMPS2::CASPT2::shift_D_nonactive( const DMRGSCFindices * idx, const int irrep_i, const int irrep_a ){
2409 
2410    const int irrep_ia = Irreps::directProd( irrep_i, irrep_a );
2411    const int n_irreps = idx->getNirreps();
2412 
2413    int shift = 0;
2414    for ( int Ii = 0; Ii < n_irreps; Ii++ ){
2415       const int Ia = Irreps::directProd( irrep_ia, Ii );
2416       if (( Ii == irrep_i ) && ( Ia == irrep_a )){
2417          Ii = n_irreps;
2418       } else {
2419          shift += idx->getNOCC( Ii ) * idx->getNVIRT( Ia );
2420       }
2421    }
2422    return shift;
2423 
2424 }
2425 
diagonal(double * result) const2426 void CheMPS2::CASPT2::diagonal( double * result ) const{
2427 
2428    #pragma omp parallel
2429    {
2430 
2431       // FAA: < E_zy E_jx | F | E_ti E_uv > = delta_ij * ( FAA[ Ii ][ xyztuv ] + ( 2 sum_k f_kk - f_ii ) SAA[ Ii ][ xyztuv ] )
2432       for ( int irrep = 0; irrep < num_irreps; irrep++ ){
2433          const int SIZE = size_A[ irrep ];
2434          if ( SIZE > 0 ){
2435             const int NOCC = indices->getNOCC( irrep );
2436             #pragma omp for schedule(static)
2437             for ( int count = 0; count < NOCC; count++ ){
2438                const double beta = - f_dot_1dm - fock->get( irrep, count, count );
2439                double * target = result + jump[ irrep + num_irreps * CHEMPS2_CASPT2_A ] + SIZE * count;
2440                #pragma omp simd
2441                for ( int elem = 0; elem < SIZE; elem++ ){ target[ elem ] = FAA[ irrep ][ elem ] + beta; }
2442             }
2443          }
2444       }
2445 
2446       // FCC: < E_zy E_xb | F | E_at E_uv > = delta_ab * ( FCC[ Ia ][ xyztuv ] + ( 2 sum_k f_kk + f_aa ) SCC[ Ia ][ xyztuv ] )
2447       for ( int irrep = 0; irrep < num_irreps; irrep++ ){
2448          const int SIZE = size_C[ irrep ];
2449          if ( SIZE > 0 ){
2450             const int NVIR = indices->getNVIRT( irrep );
2451             const int N_OA = indices->getNOCC( irrep ) + indices->getNDMRG( irrep );
2452             #pragma omp for schedule(static)
2453             for ( int count = 0; count < NVIR; count++ ){
2454                const double beta = - f_dot_1dm + fock->get( irrep, N_OA + count, N_OA + count );
2455                double * target = result + jump[ irrep + num_irreps * CHEMPS2_CASPT2_C ] + SIZE * count;
2456                #pragma omp simd
2457                for ( int elem = 0; elem < SIZE; elem++ ){ target[ elem ] = FCC[ irrep ][ elem ] + beta; }
2458             }
2459          }
2460       }
2461 
2462       // FDD: < E_yx E_jb | F | E_ai E_tu > = delta_ab delta_ij ( FDD[ xytu] + ( 2 sum_k f_kk + f_aa - f_ii ) SDD[ xytu ] )
2463       for ( int irrep = 0; irrep < num_irreps; irrep++ ){
2464          const int SIZE = size_D[ irrep ];
2465          if ( SIZE > 0 ){
2466             int shift = 0;
2467             for ( int irrep_i = 0; irrep_i < num_irreps; irrep_i++ ){
2468                const int irrep_a = Irreps::directProd( irrep_i, irrep );
2469                const int NOCC_i  = indices->getNOCC( irrep_i );
2470                const int N_OA_a  = indices->getNOCC( irrep_a ) + indices->getNDMRG( irrep_a );
2471                const int SIZE_ia = NOCC_i * indices->getNVIRT( irrep_a );
2472                #pragma omp for schedule(static)
2473                for ( int combined = 0; combined < SIZE_ia; combined++ ){
2474                   const int i = combined % NOCC_i;
2475                   const int a = combined / NOCC_i;
2476                   const double f_ii = fock->get( irrep_i, i, i );
2477                   const double f_aa = fock->get( irrep_a, N_OA_a + a, N_OA_a + a );
2478                   const double beta = - f_dot_1dm + f_aa - f_ii;
2479                   double * target = result + jump[ irrep + num_irreps * CHEMPS2_CASPT2_D ] + SIZE * ( shift + combined );
2480                   #pragma omp simd
2481                   for ( int elem = 0; elem < SIZE; elem++ ){ target[ elem ] = FDD[ irrep ][ elem ] + beta; }
2482                }
2483                shift += SIZE_ia;
2484             }
2485          }
2486       }
2487 
2488       int triangle_idx[] = { 0, 0 };
2489 
2490       // FBB singlet: < SB_xkyl | F | SB_tiuj > = 2 delta_ik delta_jl ( FBB_singlet[ Iij ][ xytu ] + ( 2 sum_n f_nn - f_ii - f_jj ) * SBB_singlet[ Iij ][ xytu ] )
2491       {
2492          const int SIZE = size_B_singlet[ 0 ];
2493          if ( SIZE > 0 ){
2494             int shift = 0;
2495             for ( int irrep_ij = 0; irrep_ij < num_irreps; irrep_ij++ ){
2496                const int nocc_ij = indices->getNOCC( irrep_ij );
2497                const int size_ij = ( nocc_ij * ( nocc_ij + 1 ) ) / 2;
2498                #pragma omp for schedule(static)
2499                for ( int combined = 0; combined < size_ij; combined++ ){
2500                   Special::invert_triangle_two( combined, triangle_idx );
2501                   const int i = triangle_idx[ 0 ];
2502                   const int j = triangle_idx[ 1 ];
2503                   const double f_ii = fock->get( irrep_ij, i, i );
2504                   const double f_jj = fock->get( irrep_ij, j, j );
2505                   const double beta = - f_dot_1dm - f_ii - f_jj;
2506                   double * target = result + jump[ num_irreps * CHEMPS2_CASPT2_B_SINGLET ] + SIZE * ( shift + combined );
2507                   #pragma omp simd
2508                   for ( int elem = 0; elem < SIZE; elem++ ){ target[ elem ] = 2 * ( FBB_singlet[ 0 ][ elem ] + beta ); }
2509                }
2510                shift += size_ij;
2511             }
2512          }
2513       }
2514       for ( int irrep = 1; irrep < num_irreps; irrep++ ){
2515          const int SIZE = size_B_singlet[ irrep ];
2516          if ( SIZE > 0 ){
2517             int shift = 0;
2518             for ( int irrep_i = 0; irrep_i < num_irreps; irrep_i++ ){
2519                const int irrep_j = Irreps::directProd( irrep, irrep_i );
2520                if ( irrep_i < irrep_j ){
2521                   const int nocc_i  = indices->getNOCC( irrep_i );
2522                   const int size_ij = nocc_i * indices->getNOCC( irrep_j );
2523                   #pragma omp for schedule(static)
2524                   for ( int combined = 0; combined < size_ij; combined++ ){
2525                      const int i = combined % nocc_i;
2526                      const int j = combined / nocc_i;
2527                      const double f_ii = fock->get( irrep_i, i, i );
2528                      const double f_jj = fock->get( irrep_j, j, j );
2529                      const double beta = - f_dot_1dm - f_ii - f_jj;
2530                      double * target = result + jump[ irrep + num_irreps * CHEMPS2_CASPT2_B_SINGLET ] + SIZE * ( shift + combined );
2531                      #pragma omp simd
2532                      for ( int elem = 0; elem < SIZE; elem++ ){ target[ elem ] = 2 * ( FBB_singlet[ irrep ][ elem ] + beta ); }
2533                   }
2534                   shift += size_ij;
2535                }
2536             }
2537          }
2538       }
2539 
2540       // FBB triplet: < TB_xkyl | F | TB_tiuj > = 2 delta_ik delta_jl ( FBB_triplet[ Iij ][ xytu ] + ( 2 sum_n f_nn - f_ii - f_jj ) * SBB_triplet[ Iij ][ xytu ] )
2541       {
2542          const int SIZE = size_B_triplet[ 0 ];
2543          if ( SIZE > 0 ){
2544             int shift = 0;
2545             for ( int irrep_ij = 0; irrep_ij < num_irreps; irrep_ij++ ){
2546                const int nocc_ij = indices->getNOCC( irrep_ij );
2547                const int size_ij = ( nocc_ij * ( nocc_ij - 1 ) ) / 2;
2548                #pragma omp for schedule(static)
2549                for ( int combined = 0; combined < size_ij; combined++ ){
2550                   Special::invert_lower_triangle_two( combined, triangle_idx );
2551                   const int i = triangle_idx[ 0 ];
2552                   const int j = triangle_idx[ 1 ];
2553                   const double f_ii = fock->get( irrep_ij, i, i );
2554                   const double f_jj = fock->get( irrep_ij, j, j );
2555                   const double beta = - f_dot_1dm - f_ii - f_jj;
2556                   double * target = result + jump[ num_irreps * CHEMPS2_CASPT2_B_TRIPLET ] + SIZE * ( shift + combined );
2557                   #pragma omp simd
2558                   for ( int elem = 0; elem < SIZE; elem++ ){ target[ elem ] = 2 * ( FBB_triplet[ 0 ][ elem ] + beta ); }
2559                }
2560                shift += size_ij;
2561             }
2562          }
2563       }
2564       for ( int irrep = 1; irrep < num_irreps; irrep++ ){
2565          const int SIZE = size_B_triplet[ irrep ];
2566          if ( SIZE > 0 ){
2567             int shift = 0;
2568             for ( int irrep_i = 0; irrep_i < num_irreps; irrep_i++ ){
2569                const int irrep_j = Irreps::directProd( irrep, irrep_i );
2570                if ( irrep_i < irrep_j ){
2571                   const int nocc_i  = indices->getNOCC( irrep_i );
2572                   const int size_ij = nocc_i * indices->getNOCC( irrep_j );
2573                   #pragma omp for schedule(static)
2574                   for ( int combined = 0; combined < size_ij; combined++ ){
2575                      const int i = combined % nocc_i;
2576                      const int j = combined / nocc_i;
2577                      const double f_ii = fock->get( irrep_i, i, i );
2578                      const double f_jj = fock->get( irrep_j, j, j );
2579                      const double beta = - f_dot_1dm - f_ii - f_jj;
2580                      double * target = result + jump[ irrep + num_irreps * CHEMPS2_CASPT2_B_TRIPLET ] + SIZE * ( shift + combined );
2581                      #pragma omp simd
2582                      for ( int elem = 0; elem < SIZE; elem++ ){ target[ elem ] = 2 * ( FBB_triplet[ irrep ][ elem ] + beta ); }
2583                   }
2584                   shift += size_ij;
2585                }
2586             }
2587          }
2588       }
2589 
2590       // FFF singlet: < SF_cxdy | F | SF_atbu > = 2 delta_ac delta_bd ( FFF_singlet[ Iab ][ xytu ] + ( 2 sum_n f_nn + f_aa + f_bb ) * SFF_singlet[ Iab ][ xytu ] )
2591       {
2592          const int SIZE = size_F_singlet[ 0 ];
2593          if ( SIZE > 0 ){
2594             int shift = 0;
2595             for ( int irrep_ab = 0; irrep_ab < num_irreps; irrep_ab++ ){
2596                const int NVIR_ab = indices->getNVIRT( irrep_ab );
2597                const int N_OA_ab = indices->getNOCC( irrep_ab ) + indices->getNDMRG( irrep_ab );
2598                const int SIZE_ab = ( NVIR_ab * ( NVIR_ab + 1 ) ) / 2;
2599                #pragma omp for schedule(static)
2600                for ( int combined = 0; combined < SIZE_ab; combined++ ){
2601                   Special::invert_triangle_two( combined, triangle_idx );
2602                   const int a = triangle_idx[ 0 ];
2603                   const int b = triangle_idx[ 1 ];
2604                   const double f_aa = fock->get( irrep_ab, N_OA_ab + a, N_OA_ab + a );
2605                   const double f_bb = fock->get( irrep_ab, N_OA_ab + b, N_OA_ab + b );
2606                   const double beta = - f_dot_1dm + f_aa + f_bb;
2607                   double * target = result + jump[ num_irreps * CHEMPS2_CASPT2_F_SINGLET ] + SIZE * ( shift + combined );
2608                   #pragma omp simd
2609                   for ( int elem = 0; elem < SIZE; elem++ ){ target[ elem ] = 2 * ( FFF_singlet[ 0 ][ elem ] + beta ); }
2610                }
2611                shift += SIZE_ab;
2612             }
2613          }
2614       }
2615       for ( int irrep = 1; irrep < num_irreps; irrep++ ){
2616          const int SIZE = size_F_singlet[ irrep ];
2617          if ( SIZE > 0 ){
2618             int shift = 0;
2619             for ( int irrep_a = 0; irrep_a < num_irreps; irrep_a++ ){
2620                const int irrep_b = Irreps::directProd( irrep, irrep_a );
2621                if ( irrep_a < irrep_b ){
2622                   const int NVIR_a  = indices->getNVIRT( irrep_a );
2623                   const int SIZE_ab = NVIR_a * indices->getNVIRT( irrep_b );
2624                   const int N_OA_a  = indices->getNOCC( irrep_a ) + indices->getNDMRG( irrep_a );
2625                   const int N_OA_b  = indices->getNOCC( irrep_b ) + indices->getNDMRG( irrep_b );
2626                   #pragma omp for schedule(static)
2627                   for ( int combined = 0; combined < SIZE_ab; combined++ ){
2628                      const int a = combined % NVIR_a;
2629                      const int b = combined / NVIR_a;
2630                      const double f_aa = fock->get( irrep_a, N_OA_a + a, N_OA_a + a );
2631                      const double f_bb = fock->get( irrep_b, N_OA_b + b, N_OA_b + b );
2632                      const double beta = - f_dot_1dm + f_aa + f_bb;
2633                      double * target = result + jump[ irrep + num_irreps * CHEMPS2_CASPT2_F_SINGLET ] + SIZE * ( shift + combined );
2634                      #pragma omp simd
2635                      for ( int elem = 0; elem < SIZE; elem++ ){ target[ elem ] = 2 * ( FFF_singlet[ irrep ][ elem ] + beta ); }
2636                   }
2637                   shift += SIZE_ab;
2638                }
2639             }
2640          }
2641       }
2642 
2643       // FFF triplet: < TF_cxdy | F | TF_atbu > = 2 delta_ac delta_bd ( FFF_triplet[ Iab ][ xytu ] + ( 2 sum_n f_nn + f_aa + f_bb ) * SFF_triplet[ Iab ][ xytu ] )
2644       {
2645          const int SIZE = size_F_triplet[ 0 ];
2646          if ( SIZE > 0 ){
2647             int shift = 0;
2648             for ( int irrep_ab = 0; irrep_ab < num_irreps; irrep_ab++ ){
2649                const int NVIR_ab = indices->getNVIRT( irrep_ab );
2650                const int N_OA_ab = indices->getNOCC( irrep_ab ) + indices->getNDMRG( irrep_ab );
2651                const int SIZE_ab = ( NVIR_ab * ( NVIR_ab - 1 ) ) / 2;
2652                #pragma omp for schedule(static)
2653                for ( int combined = 0; combined < SIZE_ab; combined++ ){
2654                   Special::invert_lower_triangle_two( combined, triangle_idx );
2655                   const int a = triangle_idx[ 0 ];
2656                   const int b = triangle_idx[ 1 ];
2657                   const double f_aa = fock->get( irrep_ab, N_OA_ab + a, N_OA_ab + a );
2658                   const double f_bb = fock->get( irrep_ab, N_OA_ab + b, N_OA_ab + b );
2659                   const double beta = - f_dot_1dm + f_aa + f_bb;
2660                   double * target = result + jump[ num_irreps * CHEMPS2_CASPT2_F_TRIPLET ] + SIZE * ( shift + combined );
2661                   #pragma omp simd
2662                   for ( int elem = 0; elem < SIZE; elem++ ){ target[ elem ] = 2 * ( FFF_triplet[ 0 ][ elem ] + beta ); }
2663                }
2664                shift += SIZE_ab;
2665             }
2666          }
2667       }
2668       for ( int irrep = 1; irrep < num_irreps; irrep++ ){
2669          const int SIZE = size_F_triplet[ irrep ];
2670          if ( SIZE > 0 ){
2671             int shift = 0;
2672             for ( int irrep_a = 0; irrep_a < num_irreps; irrep_a++ ){
2673                const int irrep_b = Irreps::directProd( irrep, irrep_a );
2674                if ( irrep_a < irrep_b ){
2675                   const int NVIR_a  = indices->getNVIRT( irrep_a );
2676                   const int SIZE_ab = NVIR_a * indices->getNVIRT( irrep_b );
2677                   const int N_OA_a  = indices->getNOCC( irrep_a ) + indices->getNDMRG( irrep_a );
2678                   const int N_OA_b  = indices->getNOCC( irrep_b ) + indices->getNDMRG( irrep_b );
2679                   #pragma omp for schedule(static)
2680                   for ( int combined = 0; combined < SIZE_ab; combined++ ){
2681                      const int a = combined % NVIR_a;
2682                      const int b = combined / NVIR_a;
2683                      const double f_aa = fock->get( irrep_a, N_OA_a + a, N_OA_a + a );
2684                      const double f_bb = fock->get( irrep_b, N_OA_b + b, N_OA_b + b );
2685                      const double beta = - f_dot_1dm + f_aa + f_bb;
2686                      double * target = result + jump[ irrep + num_irreps * CHEMPS2_CASPT2_F_TRIPLET ] + SIZE * ( shift + combined );
2687                      #pragma omp simd
2688                      for ( int elem = 0; elem < SIZE; elem++ ){ target[ elem ] = 2 * ( FFF_triplet[ irrep ][ elem ] + beta ); }
2689                   }
2690                   shift += SIZE_ab;
2691                }
2692             }
2693          }
2694       }
2695 
2696       // FEE singlet: < SE_ukbl | F | SE_tiaj > = 2 delta_ab delta_ik delta_jl ( FEE[ It ][ ut ] + ( 2 sum_k f_kk + f_aa - f_ii - f_jj ) SEE[ It ][ ut ] )
2697       for ( int irrep = 0; irrep < num_irreps; irrep++ ){
2698          const int SIZE = size_E[ irrep ];
2699          if ( SIZE > 0 ){
2700             int shift = 0;
2701             for ( int irrep_a = 0; irrep_a < num_irreps; irrep_a++ ){
2702                const int NVIR_a = indices->getNVIRT( irrep_a );
2703                const int N_OA_a = indices->getNOCC( irrep_a ) + indices->getNDMRG( irrep_a );
2704                const int irrep_occ = Irreps::directProd( irrep_a, irrep );
2705                for ( int irrep_i = 0; irrep_i < num_irreps; irrep_i++ ){
2706                   const int NOCC_i = indices->getNOCC( irrep_i );
2707                   const int irrep_j = Irreps::directProd( irrep_occ, irrep_i );
2708                   if ( irrep_i == irrep_j ){
2709                      const int SIZE_aij = ( NVIR_a * NOCC_i * ( NOCC_i + 1 ) ) / 2;
2710                      #pragma omp for schedule(static)
2711                      for ( int combined = 0; combined < SIZE_aij; combined++ ){
2712                         const int a = combined % NVIR_a;
2713                         Special::invert_triangle_two( combined / NVIR_a, triangle_idx );
2714                         const int i = triangle_idx[ 0 ];
2715                         const int j = triangle_idx[ 1 ];
2716                         const double f_ii = fock->get( irrep_i, i, i );
2717                         const double f_jj = fock->get( irrep_j, j, j );
2718                         const double f_aa = fock->get( irrep_a, N_OA_a + a, N_OA_a + a );
2719                         const double beta = - f_dot_1dm + f_aa - f_ii - f_jj;
2720                         double * target = result + jump[ irrep + num_irreps * CHEMPS2_CASPT2_E_SINGLET ] + SIZE * ( shift + combined );
2721                         #pragma omp simd
2722                         for ( int elem = 0; elem < SIZE; elem++ ){ target[ elem ] = 2 * ( FEE[ irrep ][ elem ] + beta ); }
2723                      }
2724                      shift += SIZE_aij;
2725                   }
2726                   if ( irrep_i < irrep_j ){
2727                      const int SIZE_aij = NVIR_a * NOCC_i * indices->getNOCC( irrep_j );
2728                      #pragma omp for schedule(static)
2729                      for ( int combined = 0; combined < SIZE_aij; combined++ ){
2730                         const int a = combined % NVIR_a;
2731                         const int temp = combined / NVIR_a;
2732                         const int i = temp % NOCC_i;
2733                         const int j = temp / NOCC_i;
2734                         const double f_ii = fock->get( irrep_i, i, i );
2735                         const double f_jj = fock->get( irrep_j, j, j );
2736                         const double f_aa = fock->get( irrep_a, N_OA_a + a, N_OA_a + a );
2737                         const double beta = - f_dot_1dm + f_aa - f_ii - f_jj;
2738                         double * target = result + jump[ irrep + num_irreps * CHEMPS2_CASPT2_E_SINGLET ] + SIZE * ( shift + combined );
2739                         #pragma omp simd
2740                         for ( int elem = 0; elem < SIZE; elem++ ){ target[ elem ] = 2 * ( FEE[ irrep ][ elem ] + beta ); }
2741                      }
2742                      shift += SIZE_aij;
2743                   }
2744                }
2745             }
2746          }
2747       }
2748 
2749       // FEE triplet: < TE_ukbl | F | TE_tiaj > = 6 delta_ab delta_ik delta_jl ( FEE[ It ][ ut ] + ( 2 sum_k f_kk + f_aa - f_ii - f_jj ) SEE[ It ][ ut ] )
2750       for ( int irrep = 0; irrep < num_irreps; irrep++ ){
2751          const int SIZE = size_E[ irrep ];
2752          if ( SIZE > 0 ){
2753             int shift = 0;
2754             for ( int irrep_a = 0; irrep_a < num_irreps; irrep_a++ ){
2755                const int NVIR_a = indices->getNVIRT( irrep_a );
2756                const int N_OA_a = indices->getNOCC( irrep_a ) + indices->getNDMRG( irrep_a );
2757                const int irrep_occ = Irreps::directProd( irrep_a, irrep );
2758                for ( int irrep_i = 0; irrep_i < num_irreps; irrep_i++ ){
2759                   const int NOCC_i = indices->getNOCC( irrep_i );
2760                   const int irrep_j = Irreps::directProd( irrep_occ, irrep_i );
2761                   if ( irrep_i == irrep_j ){
2762                      const int SIZE_aij = ( NVIR_a * NOCC_i * ( NOCC_i - 1 ) ) / 2;
2763                      #pragma omp for schedule(static)
2764                      for ( int combined = 0; combined < SIZE_aij; combined++ ){
2765                         const int a = combined % NVIR_a;
2766                         Special::invert_lower_triangle_two( combined / NVIR_a, triangle_idx );
2767                         const int i = triangle_idx[ 0 ];
2768                         const int j = triangle_idx[ 1 ];
2769                         const double f_ii = fock->get( irrep_i, i, i );
2770                         const double f_jj = fock->get( irrep_j, j, j );
2771                         const double f_aa = fock->get( irrep_a, N_OA_a + a, N_OA_a + a );
2772                         const double beta = - f_dot_1dm + f_aa - f_ii - f_jj;
2773                         double * target = result + jump[ irrep + num_irreps * CHEMPS2_CASPT2_E_TRIPLET ] + SIZE * ( shift + combined );
2774                         #pragma omp simd
2775                         for ( int elem = 0; elem < SIZE; elem++ ){ target[ elem ] = 6 * ( FEE[ irrep ][ elem ] + beta ); }
2776                      }
2777                      shift += SIZE_aij;
2778                   }
2779                   if ( irrep_i < irrep_j ){
2780                      const int SIZE_aij = NVIR_a * NOCC_i * indices->getNOCC( irrep_j );
2781                      #pragma omp for schedule(static)
2782                      for ( int combined = 0; combined < SIZE_aij; combined++ ){
2783                         const int a = combined % NVIR_a;
2784                         const int temp = combined / NVIR_a;
2785                         const int i = temp % NOCC_i;
2786                         const int j = temp / NOCC_i;
2787                         const double f_ii = fock->get( irrep_i, i, i );
2788                         const double f_jj = fock->get( irrep_j, j, j );
2789                         const double f_aa = fock->get( irrep_a, N_OA_a + a, N_OA_a + a );
2790                         const double beta = - f_dot_1dm + f_aa - f_ii - f_jj;
2791                         double * target = result + jump[ irrep + num_irreps * CHEMPS2_CASPT2_E_TRIPLET ] + SIZE * ( shift + combined );
2792                         #pragma omp simd
2793                         for ( int elem = 0; elem < SIZE; elem++ ){ target[ elem ] = 6 * ( FEE[ irrep ][ elem ] + beta ); }
2794                      }
2795                      shift += SIZE_aij;
2796                   }
2797                }
2798             }
2799          }
2800       }
2801 
2802       // FGG singlet: < SG_cjdu | F | SG_aibt > = 2 delta_ij delta_ac delta_bd ( FGG[ It ][ ut ] + ( 2 sum_k f_kk + f_aa + f_bb - f_ii ) SGG[ It ][ ut ] )
2803       for ( int irrep = 0; irrep < num_irreps; irrep++ ){
2804          const int SIZE = size_G[ irrep ];
2805          if ( SIZE > 0 ){
2806             int shift = 0;
2807             for ( int irrep_i = 0; irrep_i < num_irreps; irrep_i++ ){
2808                const int NOCC_i = indices->getNOCC( irrep_i );
2809                const int irrep_vir = Irreps::directProd( irrep_i, irrep );
2810                for ( int irrep_a = 0; irrep_a < num_irreps; irrep_a++ ){
2811                   const int NVIR_a = indices->getNVIRT( irrep_a );
2812                   const int N_OA_a = indices->getNOCC( irrep_a ) + indices->getNDMRG( irrep_a );
2813                   const int irrep_b = Irreps::directProd( irrep_vir, irrep_a );
2814                   if ( irrep_a == irrep_b ){
2815                      const int SIZE_iab = ( NOCC_i * NVIR_a * ( NVIR_a + 1 ) ) / 2;
2816                      #pragma omp for schedule(static)
2817                      for ( int combined = 0; combined < SIZE_iab; combined++ ){
2818                         const int i = combined % NOCC_i;
2819                         Special::invert_triangle_two( combined / NOCC_i, triangle_idx );
2820                         const int a = triangle_idx[ 0 ];
2821                         const int b = triangle_idx[ 1 ];
2822                         const double f_ii = fock->get( irrep_i, i, i );
2823                         const double f_aa = fock->get( irrep_a, N_OA_a + a, N_OA_a + a );
2824                         const double f_bb = fock->get( irrep_b, N_OA_a + b, N_OA_a + b );
2825                         const double beta = - f_dot_1dm + f_aa + f_bb - f_ii;
2826                         double * target = result + jump[ irrep + num_irreps * CHEMPS2_CASPT2_G_SINGLET ] + SIZE * ( shift + combined );
2827                         #pragma omp simd
2828                         for ( int elem = 0; elem < SIZE; elem++ ){ target[ elem ] = 2 * ( FGG[ irrep ][ elem ] + beta ); }
2829                      }
2830                      shift += SIZE_iab;
2831                   }
2832                   if ( irrep_a < irrep_b ){
2833                      const int SIZE_iab = NOCC_i * NVIR_a * indices->getNVIRT( irrep_b );
2834                      const int N_OA_b = indices->getNOCC( irrep_b ) + indices->getNDMRG( irrep_b );
2835                      #pragma omp for schedule(static)
2836                      for ( int combined = 0; combined < SIZE_iab; combined++ ){
2837                         const int i = combined % NOCC_i;
2838                         const int temp = combined / NOCC_i;
2839                         const int a = temp % NVIR_a;
2840                         const int b = temp / NVIR_a;
2841                         const double f_ii = fock->get( irrep_i, i, i );
2842                         const double f_aa = fock->get( irrep_a, N_OA_a + a, N_OA_a + a );
2843                         const double f_bb = fock->get( irrep_b, N_OA_b + b, N_OA_b + b );
2844                         const double beta = - f_dot_1dm + f_aa + f_bb - f_ii;
2845                         double * target = result + jump[ irrep + num_irreps * CHEMPS2_CASPT2_G_SINGLET ] + SIZE * ( shift + combined );
2846                         #pragma omp simd
2847                         for ( int elem = 0; elem < SIZE; elem++ ){ target[ elem ] = 2 * ( FGG[ irrep ][ elem ] + beta ); }
2848                      }
2849                      shift += SIZE_iab;
2850                   }
2851                }
2852             }
2853          }
2854       }
2855 
2856       // FGG triplet: < TG_cjdu | F | TG_aibt > = 6 delta_ij delta_ac delta_bd ( FGG[ It ][ ut ] + ( 2 sum_k f_kk + f_aa + f_bb - f_ii ) SGG[ It ][ ut ] )
2857       for ( int irrep = 0; irrep < num_irreps; irrep++ ){
2858          const int SIZE = size_G[ irrep ];
2859          if ( SIZE > 0 ){
2860             int shift = 0;
2861             for ( int irrep_i = 0; irrep_i < num_irreps; irrep_i++ ){
2862                const int NOCC_i = indices->getNOCC( irrep_i );
2863                const int irrep_vir = Irreps::directProd( irrep_i, irrep );
2864                for ( int irrep_a = 0; irrep_a < num_irreps; irrep_a++ ){
2865                   const int NVIR_a = indices->getNVIRT( irrep_a );
2866                   const int N_OA_a = indices->getNOCC( irrep_a ) + indices->getNDMRG( irrep_a );
2867                   const int irrep_b = Irreps::directProd( irrep_vir, irrep_a );
2868                   if ( irrep_a == irrep_b ){
2869                      const int SIZE_iab = ( NOCC_i * NVIR_a * ( NVIR_a - 1 ) ) / 2;
2870                      #pragma omp for schedule(static)
2871                      for ( int combined = 0; combined < SIZE_iab; combined++ ){
2872                         const int i = combined % NOCC_i;
2873                         Special::invert_lower_triangle_two( combined / NOCC_i, triangle_idx );
2874                         const int a = triangle_idx[ 0 ];
2875                         const int b = triangle_idx[ 1 ];
2876                         const double f_ii = fock->get( irrep_i, i, i );
2877                         const double f_aa = fock->get( irrep_a, N_OA_a + a, N_OA_a + a );
2878                         const double f_bb = fock->get( irrep_b, N_OA_a + b, N_OA_a + b );
2879                         const double beta = - f_dot_1dm + f_aa + f_bb - f_ii;
2880                         double * target = result + jump[ irrep + num_irreps * CHEMPS2_CASPT2_G_TRIPLET ] + SIZE * ( shift + combined );
2881                         #pragma omp simd
2882                         for ( int elem = 0; elem < SIZE; elem++ ){ target[ elem ] = 6 * ( FGG[ irrep ][ elem ] + beta ); }
2883                      }
2884                      shift += SIZE_iab;
2885                   }
2886                   if ( irrep_a < irrep_b ){
2887                      const int SIZE_iab = NOCC_i * NVIR_a * indices->getNVIRT( irrep_b );
2888                      const int N_OA_b = indices->getNOCC( irrep_b ) + indices->getNDMRG( irrep_b );
2889                      #pragma omp for schedule(static)
2890                      for ( int combined = 0; combined < SIZE_iab; combined++ ){
2891                         const int i = combined % NOCC_i;
2892                         const int temp = combined / NOCC_i;
2893                         const int a = temp % NVIR_a;
2894                         const int b = temp / NVIR_a;
2895                         const double f_ii = fock->get( irrep_i, i, i );
2896                         const double f_aa = fock->get( irrep_a, N_OA_a + a, N_OA_a + a );
2897                         const double f_bb = fock->get( irrep_b, N_OA_b + b, N_OA_b + b );
2898                         const double beta = - f_dot_1dm + f_aa + f_bb - f_ii;
2899                         double * target = result + jump[ irrep + num_irreps * CHEMPS2_CASPT2_G_TRIPLET ] + SIZE * ( shift + combined );
2900                         #pragma omp simd
2901                         for ( int elem = 0; elem < SIZE; elem++ ){ target[ elem ] = 6 * ( FGG[ irrep ][ elem ] + beta ); }
2902                      }
2903                      shift += SIZE_iab;
2904                   }
2905                }
2906             }
2907          }
2908       }
2909 
2910       // FHH singlet and triplet
2911       {
2912          int shift = 0;
2913          for ( int irrep_ij = 0; irrep_ij < num_irreps; irrep_ij++ ){
2914             const int NOCC_ij = indices->getNOCC( irrep_ij );
2915             const int size_ij = ( NOCC_ij * ( NOCC_ij + 1 ) ) / 2;
2916             for ( int irrep_ab = 0; irrep_ab < num_irreps; irrep_ab++ ){
2917                const int NVIR_ab = indices->getNVIRT( irrep_ab );
2918                const int N_OA_ab = indices->getNOCC( irrep_ab ) + indices->getNDMRG( irrep_ab );
2919                double * target = result + jump[ num_irreps * CHEMPS2_CASPT2_H_SINGLET ] + shift;
2920                const int size_ijab = ( size_ij * NVIR_ab * ( NVIR_ab + 1 ) ) / 2;
2921                #pragma omp for schedule(static)
2922                for ( int combined = 0; combined < size_ijab; combined++ ){
2923                   Special::invert_triangle_two( combined % size_ij, triangle_idx );
2924                   const int i = triangle_idx[ 0 ];
2925                   const int j = triangle_idx[ 1 ];
2926                   Special::invert_triangle_two( combined / size_ij, triangle_idx );
2927                   const int a = triangle_idx[ 0 ];
2928                   const int b = triangle_idx[ 1 ];
2929                   const double f_ii = fock->get( irrep_ij, i, i );
2930                   const double f_jj = fock->get( irrep_ij, j, j );
2931                   const double f_aa = fock->get( irrep_ab, N_OA_ab + a, N_OA_ab + a );
2932                   const double f_bb = fock->get( irrep_ab, N_OA_ab + b, N_OA_ab + b );
2933                   const double term = f_aa + f_bb - f_ii - f_jj;
2934                   target[ combined ] = 4 * term;
2935                }
2936                shift += size_ijab;
2937             }
2938          }
2939       }
2940       {
2941          int shift = 0;
2942          for ( int irrep_ij = 0; irrep_ij < num_irreps; irrep_ij++ ){
2943             const int NOCC_ij = indices->getNOCC( irrep_ij );
2944             const int size_ij = ( NOCC_ij * ( NOCC_ij - 1 ) ) / 2;
2945             for ( int irrep_ab = 0; irrep_ab < num_irreps; irrep_ab++ ){
2946                const int NVIR_ab = indices->getNVIRT( irrep_ab );
2947                const int N_OA_ab = indices->getNOCC( irrep_ab ) + indices->getNDMRG( irrep_ab );
2948                double * target = result + jump[ num_irreps * CHEMPS2_CASPT2_H_TRIPLET ] + shift;
2949                const int size_ijab = ( size_ij * NVIR_ab * ( NVIR_ab - 1 ) ) / 2;
2950                #pragma omp for schedule(static)
2951                for ( int combined = 0; combined < size_ijab; combined++ ){
2952                   Special::invert_lower_triangle_two( combined % size_ij, triangle_idx );
2953                   const int i = triangle_idx[ 0 ];
2954                   const int j = triangle_idx[ 1 ];
2955                   Special::invert_lower_triangle_two( combined / size_ij, triangle_idx );
2956                   const int a = triangle_idx[ 0 ];
2957                   const int b = triangle_idx[ 1 ];
2958                   const double f_ii = fock->get( irrep_ij, i, i );
2959                   const double f_jj = fock->get( irrep_ij, j, j );
2960                   const double f_aa = fock->get( irrep_ab, N_OA_ab + a, N_OA_ab + a );
2961                   const double f_bb = fock->get( irrep_ab, N_OA_ab + b, N_OA_ab + b );
2962                   const double term = f_aa + f_bb - f_ii - f_jj;
2963                   target[ combined ] = 12 * term;
2964                }
2965                shift += size_ijab;
2966             }
2967          }
2968       }
2969       for ( int irrep = 1; irrep < num_irreps; irrep++ ){
2970          int shift = 0;
2971          for ( int irrep_i = 0; irrep_i < num_irreps; irrep_i++ ){
2972             const int irrep_j = Irreps::directProd( irrep, irrep_i );
2973             if ( irrep_i < irrep_j ){
2974                const int NOCC_i = indices->getNOCC( irrep_i );
2975                const int NOCC_j = indices->getNOCC( irrep_j );
2976                for ( int irrep_a = 0; irrep_a < num_irreps; irrep_a++ ){
2977                   const int irrep_b = Irreps::directProd( irrep, irrep_a );
2978                   if ( irrep_a < irrep_b ){
2979                      const int NVIR_a = indices->getNVIRT( irrep_a );
2980                      const int N_OA_a = indices->getNOCC( irrep_a ) + indices->getNDMRG( irrep_a );
2981                      const int N_OA_b = indices->getNOCC( irrep_b ) + indices->getNDMRG( irrep_b );
2982                      double * target_singlet = result + jump[ irrep + num_irreps * CHEMPS2_CASPT2_H_SINGLET ] + shift;
2983                      double * target_triplet = result + jump[ irrep + num_irreps * CHEMPS2_CASPT2_H_TRIPLET ] + shift;
2984                      const int size_ijab = NOCC_i * NOCC_j * NVIR_a * indices->getNVIRT( irrep_b );
2985                      #pragma omp for schedule(static)
2986                      for ( int combined = 0; combined < size_ijab; combined++ ){
2987                         const int i = combined % NOCC_i;
2988                         const int temp1 = combined / NOCC_i;
2989                         const int j = temp1 % NOCC_j;
2990                         const int temp2 = temp1 / NOCC_j;
2991                         const int a = temp2 % NVIR_a;
2992                         const int b = temp2 / NVIR_a;
2993                         const double f_ii = fock->get( irrep_i, i, i );
2994                         const double f_jj = fock->get( irrep_j, j, j );
2995                         const double f_aa = fock->get( irrep_a, N_OA_a + a, N_OA_a + a );
2996                         const double f_bb = fock->get( irrep_b, N_OA_b + b, N_OA_b + b );
2997                         const double term = f_aa + f_bb - f_ii - f_jj;
2998                         target_singlet[ combined ] =  4 * term;
2999                         target_triplet[ combined ] = 12 * term;
3000                      }
3001                      shift += size_ijab;
3002                   }
3003                }
3004             }
3005          }
3006       }
3007 
3008    }
3009 
3010 }
3011 
construct_rhs(const DMRGSCFmatrix * oei,const DMRGSCFintegrals * integrals)3012 void CheMPS2::CASPT2::construct_rhs( const DMRGSCFmatrix * oei, const DMRGSCFintegrals * integrals ){
3013 
3014    /*
3015       VA:  < H E_ti E_uv > = sum_w ( t_iw + sum_k [ 2 (iw|kk) - (ik|kw) ] ) [ 2 delta_tw Gamma_uv - Gamma_tuwv - delta_wu Gamma_tv ]
3016                            + sum_xzy (ix|zy) SAA[ Ii ][ xyztuv ]
3017 
3018       VB:  < H E_ti E_uj >
3019            < S_tiuj | H > = sum_xy (ix|jy) SBB_singlet[ It x Iu ][ xytu ] / sqrt( 1 + delta_ij )
3020                           = sum_{x<=y} [ (ix|jy) + (iy|jx) ] SBB_singlet[ It x Iu ][ xytu ] * (( x==y ) ? 0.5 : 1.0 ) / sqrt( 1 + delta_ij )
3021            < T_tiuj | H > = sum_xy (ix|jy) SBB_triplet[ It x Iu ][ xytu ]
3022                           = sum_{x<y}  [ (ix|jy) - (iy|jx) ] SBB_triplet[ It x Iu ][ xytu ]
3023 
3024       VC:  < H E_at E_uv > = sum_w ( t_wa + sum_k [ 2 (wa|kk) - (wk|ka) ] - sum_y (wy|ya) ) < E_wt E_uv >
3025                            + sum_wxy (xy|wa) < E_xy E_wt E_uv >
3026            < H E_at E_uv > = sum_w ( t_wa + sum_k [ 2 (wa|kk) - (wk|ka) ] - sum_y (wy|ya) ) [ Gamma_wutv + delta_ut Gamma_wv ]
3027                            + sum_zxy (zy|xa) SCC[ Ia ][ xyztuv ]
3028 
3029       VD1: < H E_ai E_tu > = ( t_ia + sum_k [ 2 (ia|kk) - (ik|ka) ] ) [ 2 Gamma_tu ]
3030                            + sum_xy [ (ia|yx) - 0.5 * (ix|ya) ] SD1D1[ It x Iu ][ xytu ]
3031            < H E_ai E_tu > = ( t_ia + sum_k [ 2 (ia|kk) - (ik|ka) ] ) [ 2 Gamma_tu ]
3032                            + sum_xy (ia|yx) SD1D1[ It x Iu ][ xytu ]
3033                            + sum_xy (ix|ya) SD2D1[ It x Iu ][ xytu ]
3034 
3035       VD2: < H E_ti E_au > = ( t_ia + sum_k [ 2 (ia|kk) - (ik|ka) ] ) [ - Gamma_tu ]
3036                            + sum_xy (ia|xy) [ - Gamma_xtyu ]
3037                            + sum_xy (iy|xa) [ - Gamma_txyu ]
3038                            + sum_y [ 2 (it|ya) - (ia|yt) ] Gamma_yu
3039            < H E_ti E_au > = ( t_ia + sum_k [ 2 (ia|kk) - (ik|ka) ] ) [ - Gamma_tu ]
3040                            + sum_xy (ia|yx) SD1D2[ It x Iu ][ xytu ]
3041                            + sum_xy (ix|ya) SD2D2[ It x Iu ][ xytu ]
3042 
3043       VE:  < H E_ti E_aj >
3044            < S_tiaj | H > = sum_w [ (aj|wi) + (ai|wj) ] * 1 * SEE[ It ][ wt ] / sqrt( 1 + delta_ij )
3045            < T_tiaj | H > = sum_w [ (aj|wi) - (ai|wj) ] * 3 * SEE[ It ][ wt ]
3046 
3047       VF:  < H E_at E_bu >
3048            < S_atbu | H > = sum_xy (ax|by) SFF_singlet[ It x Iu ][ xytu ] / sqrt( 1 + delta_ab )
3049                           = sum_{x<=y} [ (ax|by) + (ay|bx) ] SFF_singlet[ It x Iu ][ xytu ] * (( x==y ) ? 0.5 : 1.0 ) / sqrt( 1 + delta_ab )
3050            < T_atbu | H > = sum_xy (ax|by) SFF_triplet[ It x Iu ][ xytu ]
3051                           = sum_{x<y}  [ (ax|by) - (ay|bx) ] SFF_triplet[ It x Iu ][ xytu ]
3052 
3053       VG:  < H E_ai E_bt >
3054            < S_aibt | H > = sum_w [ (ai|bw) + (bi|aw) ] * 1 * SGG[ It ][ wt ] / sqrt( 1 + delta_ab )
3055            < T_aibt | H > = sum_w [ (ai|bw) - (bi|aw) ] * 3 * SGG[ It ][ wt ]
3056 
3057       VH:  < H E_ai E_bj >
3058            < S_aibj | H > = 2 * [ (ai|bj) + (aj|bi) ] / sqrt( ( 1 + delta_ij ) * ( 1 + delta_ab ) )
3059            < T_aibj | H > = 6 * [ (ai|bj) - (aj|bi) ]
3060    */
3061 
3062    const int LAS = indices->getDMRGcumulative( num_irreps );
3063 
3064    // First construct MAT[p,q] = ( t_pq + sum_k [ 2 (pq|kk) - (pk|kq) ] )
3065    DMRGSCFmatrix * MAT = new DMRGSCFmatrix( indices );
3066    for ( int irrep = 0; irrep < num_irreps; irrep++ ){
3067       const int NORB = indices->getNORB( irrep );
3068       for ( int row = 0; row < NORB; row++ ){
3069          for ( int col = row; col < NORB; col++ ){
3070             double value = oei->get( irrep, row, col );
3071             for ( int irrep_occ = 0; irrep_occ < num_irreps; irrep_occ++ ){
3072                const int NOCC = indices->getNOCC( irrep_occ );
3073                for ( int occ = 0; occ < NOCC; occ++ ){
3074                   value += ( 2 * integrals->FourIndexAPI( irrep, irrep_occ, irrep, irrep_occ, row, occ, col, occ )
3075                                - integrals->FourIndexAPI( irrep, irrep_occ, irrep_occ, irrep, row, occ, occ, col ) ); // Physics notation at 4-index
3076                }
3077             }
3078             MAT->set( irrep, row, col, value );
3079             MAT->set( irrep, col, row, value );
3080          }
3081       }
3082    }
3083 
3084    // and MAT2[p,q] = sum_y (py|yq)
3085    DMRGSCFmatrix * MAT2 = new DMRGSCFmatrix( indices );
3086    for ( int irrep = 0; irrep < num_irreps; irrep++ ){
3087       const int NORB = indices->getNORB( irrep );
3088       for ( int row = 0; row < NORB; row++ ){
3089          for ( int col = row; col < NORB; col++ ){
3090             double value = 0.0;
3091             for ( int irrep_act = 0; irrep_act < num_irreps; irrep_act++ ){
3092                const int NOCC = indices->getNOCC( irrep_act );
3093                const int NACT = indices->getNDMRG( irrep_act );
3094                for ( int act = 0; act < NACT; act++ ){
3095                   value += integrals->FourIndexAPI( irrep, irrep_act, irrep_act, irrep, row, NOCC + act, NOCC + act, col ); // Physics notation at 4-index
3096                }
3097             }
3098             MAT2->set( irrep, row, col, value );
3099             MAT2->set( irrep, col, row, value );
3100          }
3101       }
3102    }
3103 
3104    vector_rhs = new double[ jump[ num_irreps * CHEMPS2_CASPT2_NUM_CASES ] ];
3105 
3106    #pragma omp parallel
3107    {
3108       const int max_size = get_maxsize();
3109       double * workspace = new double[ max_size ];
3110 
3111       // VA
3112       for ( int irrep = 0; irrep < num_irreps; irrep++ ){
3113          if ( size_A[ irrep ] > 0 ){
3114             const int NOCC = indices->getNOCC( irrep );
3115             const int NACT = indices->getNDMRG( irrep );
3116             const int d_w  = indices->getDMRGcumulative( irrep );
3117             #pragma omp for schedule(static)
3118             for ( int count_i = 0; count_i < NOCC; count_i++ ){
3119 
3120                double * target = vector_rhs + jump[ irrep + num_irreps * CHEMPS2_CASPT2_A ] + size_A[ irrep ] * count_i;
3121 
3122                // Fill workspace[ xyz ] with (ix|zy)
3123                // Fill target[ tuv ] with sum_w MAT[i,w] [ 2 delta_tw Gamma_uv - Gamma_tuwv - delta_wu Gamma_tv ]
3124                //                           = 2 MAT[i,t] Gamma_uv - MAT[i,u] Gamma_tv - sum_w MAT[i,w] Gamma_tuwv
3125                int jump_xyz = 0;
3126                for ( int irrep_x = 0; irrep_x < num_irreps; irrep_x++ ){
3127                   const int occ_x = indices->getNOCC( irrep_x );
3128                   const int num_x = indices->getNDMRG( irrep_x );
3129                   const int d_x   = indices->getDMRGcumulative( irrep_x );
3130                   for ( int irrep_y = 0; irrep_y < num_irreps; irrep_y++ ){
3131                      const int irrep_z = Irreps::directProd( Irreps::directProd( irrep, irrep_x ), irrep_y );
3132                      const int occ_y = indices->getNOCC( irrep_y );
3133                      const int occ_z = indices->getNOCC( irrep_z );
3134                      const int num_y = indices->getNDMRG( irrep_y );
3135                      const int num_z = indices->getNDMRG( irrep_z );
3136                      const int d_y   = indices->getDMRGcumulative( irrep_y );
3137                      const int d_z   = indices->getDMRGcumulative( irrep_z );
3138 
3139                      // workspace[ xyz ] = (ix|zy)
3140                      for ( int z = 0; z < num_z; z++ ){
3141                         for ( int y = 0; y < num_y; y++ ){
3142                            for ( int x = 0; x < num_x; x++ ){
3143                               const double ix_zy = integrals->get_coulomb( irrep, irrep_x, irrep_z, irrep_y, count_i, occ_x + x, occ_z + z, occ_y + y );
3144                               workspace[ jump_xyz + x + num_x * ( y + num_y * z ) ] = ix_zy;
3145                            }
3146                         }
3147                      }
3148 
3149                      // target[ tuv ] = - sum_w MAT[i,w] Gamma_tuwv
3150                      for ( int z = 0; z < num_z; z++ ){
3151                         for ( int y = 0; y < num_y; y++ ){
3152                            for ( int x = 0; x < num_x; x++ ){
3153                               double value = 0.0;
3154                               for ( int w = 0; w < NACT; w++ ){
3155                                  value += MAT->get( irrep, count_i, NOCC + w ) * two_rdm[ d_x + x + LAS * ( d_y + y + LAS * ( d_w + w + LAS * ( d_z + z ) ) ) ];
3156                               }
3157                               target[ jump_xyz + x + num_x * ( y + num_y * z ) ] = - value;
3158                            }
3159                         }
3160                      }
3161 
3162                      // target[ tuv ] += 2 MAT[i,t] Gamma_uv
3163                      if ( irrep_x == irrep ){
3164                         for ( int z = 0; z < num_z; z++ ){
3165                            for ( int y = 0; y < num_y; y++ ){
3166                               for ( int x = 0; x < num_x; x++ ){
3167                                  target[ jump_xyz + x + num_x * ( y + num_y * z ) ] += 2 * MAT->get( irrep, count_i, occ_x + x ) * one_rdm[ d_y + y + LAS * ( d_z + z ) ];
3168                               }
3169                            }
3170                         }
3171                      }
3172 
3173                      // target[ tuv ] -= MAT[i,u] Gamma_tv
3174                      if ( irrep_y == irrep ){
3175                         for ( int z = 0; z < num_z; z++ ){
3176                            for ( int y = 0; y < num_y; y++ ){
3177                               for ( int x = 0; x < num_x; x++ ){
3178                                  target[ jump_xyz + x + num_x * ( y + num_y * z ) ] -= MAT->get( irrep, count_i, occ_y + y ) * one_rdm[ d_x + x + LAS * ( d_z + z ) ];
3179                               }
3180                            }
3181                         }
3182                      }
3183                      jump_xyz += num_x * num_y * num_z;
3184                   }
3185                }
3186                assert( jump_xyz == size_A[ irrep ] );
3187 
3188                // Perform target[ tuv ] += sum_xzy (ix|zy) SAA[ It x Iu x Iv ][ xyztuv ]
3189                char notrans = 'N';
3190                int inc1 = 1;
3191                double one = 1.0;
3192                dgemv_( &notrans, &jump_xyz, &jump_xyz, &one, SAA[ irrep ], &jump_xyz, workspace, &inc1, &one, target, &inc1 );
3193             }
3194             assert( NOCC * size_A[ irrep ] == jump[ 1 + irrep + num_irreps * CHEMPS2_CASPT2_A ] - jump[ irrep + num_irreps * CHEMPS2_CASPT2_A ] );
3195          }
3196       }
3197 
3198       // VC
3199       for ( int irrep = 0; irrep < num_irreps; irrep++ ){
3200          if ( size_C[ irrep ] > 0 ){
3201             const int NOCC = indices->getNOCC( irrep );
3202             const int NVIR = indices->getNVIRT( irrep );
3203             const int NACT = indices->getNDMRG( irrep );
3204             const int N_OA = NOCC + NACT;
3205             const int d_w  = indices->getDMRGcumulative( irrep );
3206             #pragma omp for schedule(static)
3207             for ( int count_a = 0; count_a < NVIR; count_a++ ){
3208 
3209                double * target = vector_rhs + jump[ irrep + num_irreps * CHEMPS2_CASPT2_C ] + size_C[ irrep ] * count_a;
3210 
3211                // Fill workspace[ xyz ] with (zy|xa)
3212                // Fill target[ tuv ] with sum_w (MAT[w,a] - MAT2[w,a]) [ Gamma_wutv + delta_ut Gamma_wv ]
3213                int jump_xyz = 0;
3214                for ( int irrep_x = 0; irrep_x < num_irreps; irrep_x++ ){
3215                   const int occ_x = indices->getNOCC( irrep_x );
3216                   const int num_x = indices->getNDMRG( irrep_x );
3217                   const int d_x   = indices->getDMRGcumulative( irrep_x );
3218                   for ( int irrep_y = 0; irrep_y < num_irreps; irrep_y++ ){
3219                      const int irrep_z = Irreps::directProd( Irreps::directProd( irrep, irrep_x ), irrep_y );
3220                      const int occ_y = indices->getNOCC( irrep_y );
3221                      const int occ_z = indices->getNOCC( irrep_z );
3222                      const int num_y = indices->getNDMRG( irrep_y );
3223                      const int num_z = indices->getNDMRG( irrep_z );
3224                      const int d_y   = indices->getDMRGcumulative( irrep_y );
3225                      const int d_z   = indices->getDMRGcumulative( irrep_z );
3226 
3227                      // workspace[ xyz ] = (zy|xa)
3228                      for ( int z = 0; z < num_z; z++ ){
3229                         for ( int y = 0; y < num_y; y++ ){
3230                            for ( int x = 0; x < num_x; x++ ){
3231                               const double zy_xa = integrals->get_coulomb( irrep_z, irrep_y, irrep_x, irrep, occ_z + z, occ_y + y, occ_x + x, N_OA + count_a );
3232                               workspace[ jump_xyz + x + num_x * ( y + num_y * z ) ] = zy_xa;
3233                            }
3234                         }
3235                      }
3236 
3237                      // target[ tuv ] = sum_w ( MAT[w,a] - MAT2[w,a] ) Gamma_wutv
3238                      for ( int z = 0; z < num_z; z++ ){
3239                         for ( int y = 0; y < num_y; y++ ){
3240                            for ( int x = 0; x < num_x; x++ ){
3241                               double value = 0.0;
3242                               for ( int w = 0; w < NACT; w++ ){
3243                                  value += ( ( MAT->get( irrep, NOCC + w, N_OA + count_a ) - MAT2->get( irrep, NOCC + w, N_OA + count_a ) )
3244                                           * two_rdm[ d_w + w + LAS * ( d_y + y + LAS * ( d_x + x + LAS * ( d_z + z ) ) ) ] );
3245                               }
3246                               target[ jump_xyz + x + num_x * ( y + num_y * z ) ] = value;
3247                            }
3248                         }
3249                      }
3250 
3251                      // target[ tuv ] += sum_w ( MAT[w,a] - MAT2[w,a] ) delta_ut Gamma_wv
3252                      if (( irrep_z == irrep ) && ( irrep_x == irrep_y )){
3253                         for ( int z = 0; z < num_z; z++ ){ // v
3254                            double value = 0.0;
3255                            for ( int w = 0; w < NACT; w++ ){
3256                               value += ( ( MAT->get( irrep, NOCC + w, N_OA + count_a ) - MAT2->get( irrep, NOCC + w, N_OA + count_a ) )
3257                                        * one_rdm[ d_w + w + LAS * ( d_z + z ) ] );
3258                            }
3259                            for ( int xy = 0; xy < num_x; xy++ ){ // tu
3260                               target[ jump_xyz + xy + num_x * ( xy + num_x * z ) ] += value;
3261                            }
3262                         }
3263                      }
3264                      jump_xyz += num_x * num_y * num_z;
3265                   }
3266                }
3267                assert( jump_xyz == size_C[ irrep ] );
3268 
3269                // Perform target[ tuv ] += sum_zxy (zy|xa) SCC[ It x Iu x Iv ][ xyztuv ]
3270                char notrans = 'N';
3271                int inc1 = 1;
3272                double one = 1.0;
3273                dgemv_( &notrans, &jump_xyz, &jump_xyz, &one, SCC[ irrep ], &jump_xyz, workspace, &inc1, &one, target, &inc1 );
3274             }
3275             assert( NVIR * size_C[ irrep ] == jump[ 1 + irrep + num_irreps * CHEMPS2_CASPT2_C ] - jump[ irrep + num_irreps * CHEMPS2_CASPT2_C ] );
3276          }
3277       }
3278       #pragma omp single
3279       {
3280          delete MAT2;
3281       }
3282 
3283       // VD1 and VD2
3284       for ( int irrep = 0; irrep < num_irreps; irrep++ ){
3285          if ( size_D[ irrep ] > 0 ){
3286             int shift = 0;
3287             const int D2JUMP = size_D[ irrep ] / 2;
3288             for ( int irrep_i = 0; irrep_i < num_irreps; irrep_i++ ){
3289                const int irrep_a = Irreps::directProd( irrep_i, irrep );
3290                assert( shift == shift_D_nonactive( indices, irrep_i, irrep_a ) );
3291                const int NOCC_i = indices->getNOCC( irrep_i );
3292                const int N_OA_a = indices->getNOCC( irrep_a ) + indices->getNDMRG( irrep_a );
3293                const int NVIR_a = indices->getNVIRT( irrep_a );
3294                const int loopsize = NOCC_i * NVIR_a;
3295                #pragma omp for schedule(static)
3296                for ( int combined = 0; combined < loopsize; combined++ ){
3297 
3298                   const int count_i = combined % NOCC_i;
3299                   const int count_a = combined / NOCC_i;
3300                   double * target = vector_rhs + jump[ irrep + num_irreps * CHEMPS2_CASPT2_D ] + size_D[ irrep ] * ( shift + combined );
3301                   const double MAT_ia = ( ( irrep_i == irrep_a ) ? MAT->get( irrep_i, count_i, N_OA_a + count_a ) : 0.0 );
3302 
3303                   /* Fill workspace[          xy ] with (ia|yx)
3304                      Fill workspace[ D2JUMP + xy ] with (ix|ya)
3305                           then ( workspace * SDD )_D1 = sum_xy (ia|yx) SD1D1[ It x Iu ][ xytu ]
3306                                                       + sum_xy (ix|ya) SD2D1[ It x Iu ][ xytu ]
3307                            and ( workspace * SDD )_D2 = sum_xy (ia|yx) SD1D2[ It x Iu ][ xytu ]
3308                                                       + sum_xy (ix|ya) SD2D2[ It x Iu ][ xytu ]
3309                      Fill target[          tu ] = 2 MAT[i,a] Gamma_tu
3310                      Fill target[ D2JUMP + tu ] = - MAT[i,a] Gamma_tu */
3311                   int jump_xy = 0;
3312                   for ( int irrep_x = 0; irrep_x < num_irreps; irrep_x++ ){
3313                      const int irrep_y = Irreps::directProd( irrep, irrep_x );
3314                      const int occ_x = indices->getNOCC( irrep_x );
3315                      const int occ_y = indices->getNOCC( irrep_y );
3316                      const int num_x = indices->getNDMRG( irrep_x );
3317                      const int num_y = indices->getNDMRG( irrep_y );
3318 
3319                      // workspace[ xy ] = (ia|yx)
3320                      for ( int y = 0; y < num_y; y++ ){
3321                         for ( int x = 0; x < num_x; x++ ){
3322                            const double ia_yx = integrals->get_coulomb( irrep_y, irrep_x, irrep_i, irrep_a, occ_y + y, occ_x + x, count_i, N_OA_a + count_a );
3323                            workspace[ jump_xy + x + num_x * y ] = ia_yx;
3324                         }
3325                      }
3326 
3327                      // workspace[ D2JUMP + xy ] = (ix|ya)
3328                      for ( int y = 0; y < num_y; y++ ){
3329                         for ( int x = 0; x < num_x; x++ ){
3330                            const double ix_ya = integrals->get_coulomb( irrep_i, irrep_x, irrep_y, irrep_a, count_i, occ_x + x, occ_y + y, N_OA_a + count_a );
3331                            workspace[ D2JUMP + jump_xy + x + num_x * y ] = ix_ya;
3332                         }
3333                      }
3334 
3335                      // target[ tu          ] = 2 MAT[i,a] Gamma_tu
3336                      // target[ D2JUMP + tu ] = - MAT[i,a] Gamma_tu
3337                      if ( irrep == 0 ){
3338                         const int d_xy = indices->getDMRGcumulative( irrep_x );
3339                         for ( int y = 0; y < num_y; y++ ){
3340                            for ( int x = 0; x < num_x; x++ ){
3341                               const double value = MAT_ia * one_rdm[ d_xy + x + LAS * ( d_xy + y ) ];
3342                               target[          jump_xy + x + num_x * y ] = 2 * value;
3343                               target[ D2JUMP + jump_xy + x + num_x * y ] =   - value;
3344                            }
3345                         }
3346                      } else {
3347                         for ( int y = 0; y < num_y; y++ ){
3348                            for ( int x = 0; x < num_x; x++ ){
3349                               target[          jump_xy + x + num_x * y ] = 0.0;
3350                               target[ D2JUMP + jump_xy + x + num_x * y ] = 0.0;
3351                            }
3352                         }
3353                      }
3354                      jump_xy += num_x * num_y;
3355                   }
3356                   jump_xy = 2 * jump_xy;
3357                   assert( jump_xy == size_D[ irrep ] );
3358 
3359                   // Perform target += workspace * SDD[ irrep ]
3360                   char notrans = 'N';
3361                   int inc1 = 1;
3362                   double one = 1.0;
3363                   dgemv_( &notrans, &jump_xy, &jump_xy, &one, SDD[ irrep ], &jump_xy, workspace, &inc1, &one, target, &inc1 );
3364                }
3365                shift += loopsize;
3366             }
3367             assert( shift * size_D[ irrep ] == jump[ 1 + irrep + num_irreps * CHEMPS2_CASPT2_D ] - jump[ irrep + num_irreps * CHEMPS2_CASPT2_D ] );
3368          }
3369       }
3370       #pragma omp single
3371       {
3372          delete MAT;
3373       }
3374       const double SQRT_0p5 = sqrt( 0.5 );
3375       int triangle_idx[] = { 0, 0 };
3376 
3377       // VB singlet and triplet
3378       if ( size_B_singlet[ 0 ] > 0 ){ // First do irrep == Ii x Ij == Ix x Iy == It x Iu == 0
3379          int shift = 0; // First do SINGLET
3380          for ( int irrep_ij = 0; irrep_ij < num_irreps; irrep_ij++ ){
3381             assert( shift == shift_B_nonactive( indices, irrep_ij, irrep_ij, +1 ) );
3382             const int NOCC_ij = indices->getNOCC( irrep_ij );
3383             const int loopsize = ( NOCC_ij * ( NOCC_ij + 1 ) ) / 2;
3384             #pragma omp for schedule(static)
3385             for ( int combined = 0; combined < loopsize; combined++ ){
3386 
3387                Special::invert_triangle_two( combined, triangle_idx );
3388                const int i = triangle_idx[ 0 ];
3389                const int j = triangle_idx[ 1 ];
3390                assert( combined == i + ( j * ( j + 1 ) ) / 2 );
3391 
3392                // Fill workspace[ xy ] with [ (ix|jy) + (iy|jx) ] * (( x==y ) ? 0.5 : 1.0 )
3393                int jump_xy = 0;
3394                for ( int irrep_xy = 0; irrep_xy < num_irreps; irrep_xy++ ){
3395                   const int occ_xy = indices->getNOCC( irrep_xy );
3396                   const int num_xy = indices->getNDMRG( irrep_xy );
3397 
3398                   for ( int x = 0; x < num_xy; x++ ){
3399                      for ( int y = x; y < num_xy; y++ ){ // 0 <= x <= y < num_xy
3400                         const double ix_jy = integrals->get_coulomb( irrep_ij, irrep_xy, irrep_ij, irrep_xy, i, occ_xy + x, j, occ_xy + y );
3401                         const double iy_jx = integrals->get_coulomb( irrep_ij, irrep_xy, irrep_ij, irrep_xy, i, occ_xy + y, j, occ_xy + x );
3402                         workspace[ jump_xy + x + (y*(y+1))/2 ] = ( ix_jy + iy_jx ) * (( x==y ) ? 0.5 : 1.0 );
3403                      }
3404                   }
3405                   jump_xy += ( num_xy * ( num_xy + 1 ) ) / 2;
3406                }
3407                assert( jump_xy == size_B_singlet[ 0 ] );
3408 
3409                // Perform target[ tu ] = sum_{x<=y} [ (ix|jy) + (iy|jx) ] SBB_singlet[ It x Iu ][ xytu ] * (( x==y ) ? 0.5 : 1.0 ) / sqrt( 1 + delta_ij )
3410                char notrans = 'N';
3411                int inc1 = 1;
3412                double alpha = (( i == j ) ? SQRT_0p5 : 1.0 );
3413                double set = 0.0;
3414                double * target = vector_rhs + jump[ num_irreps * CHEMPS2_CASPT2_B_SINGLET ] + size_B_singlet[ 0 ] * ( shift + combined );
3415                dgemv_( &notrans, &jump_xy, &jump_xy, &alpha, SBB_singlet[ 0 ], &jump_xy, workspace, &inc1, &set, target, &inc1 );
3416             }
3417             shift += loopsize;
3418          }
3419          assert( shift * size_B_singlet[ 0 ] == jump[ 1 + num_irreps * CHEMPS2_CASPT2_B_SINGLET ] - jump[ num_irreps * CHEMPS2_CASPT2_B_SINGLET ] );
3420       }
3421       if ( size_B_triplet[ 0 ] > 0 ){ // Then do TRIPLET
3422          int shift = 0;
3423          for ( int irrep_ij = 0; irrep_ij < num_irreps; irrep_ij++ ){
3424             assert( shift == shift_B_nonactive( indices, irrep_ij, irrep_ij, -1 ) );
3425             const int NOCC_ij = indices->getNOCC( irrep_ij );
3426             const int loopsize = ( NOCC_ij * ( NOCC_ij - 1 ) ) / 2;
3427             #pragma omp for schedule(static)
3428             for ( int combined = 0; combined < loopsize; combined++ ){
3429 
3430                Special::invert_lower_triangle_two( combined, triangle_idx );
3431                const int i = triangle_idx[ 0 ];
3432                const int j = triangle_idx[ 1 ];
3433                assert( combined == i + ( j * ( j - 1 ) ) / 2 );
3434 
3435                // Fill workspace[ xy ] with [ (ix|jy) - (iy|jx) ]
3436                int jump_xy = 0;
3437                for ( int irrep_xy = 0; irrep_xy < num_irreps; irrep_xy++ ){
3438                   const int occ_xy = indices->getNOCC( irrep_xy );
3439                   const int num_xy = indices->getNDMRG( irrep_xy );
3440 
3441                   for ( int x = 0; x < num_xy; x++ ){
3442                      for ( int y = x+1; y < num_xy; y++ ){ // 0 <= x < y < num_xy
3443                         const double ix_jy = integrals->get_coulomb( irrep_ij, irrep_xy, irrep_ij, irrep_xy, i, occ_xy + x, j, occ_xy + y );
3444                         const double iy_jx = integrals->get_coulomb( irrep_ij, irrep_xy, irrep_ij, irrep_xy, i, occ_xy + y, j, occ_xy + x );
3445                         workspace[ jump_xy + x + (y*(y-1))/2 ] = ix_jy - iy_jx;
3446                      }
3447                   }
3448                   jump_xy += ( num_xy * ( num_xy - 1 ) ) / 2;
3449                }
3450                assert( jump_xy == size_B_triplet[ 0 ] );
3451 
3452                // Perform target[ tu ] = sum_{x<y} [ (ix|jy) - (iy|jx) ] SBB_triplet[ It x Iu ][ xytu ]
3453                char notrans = 'N';
3454                int inc1 = 1;
3455                double alpha = 1.0;
3456                double set = 0.0;
3457                double * target = vector_rhs + jump[ num_irreps * CHEMPS2_CASPT2_B_TRIPLET ] + size_B_triplet[ 0 ] * ( shift + combined );
3458                dgemv_( &notrans, &jump_xy, &jump_xy, &alpha, SBB_triplet[ 0 ], &jump_xy, workspace, &inc1, &set, target, &inc1 );
3459             }
3460             shift += loopsize;
3461          }
3462          assert( shift * size_B_triplet[ 0 ] == jump[ 1 + num_irreps * CHEMPS2_CASPT2_B_TRIPLET ] - jump[ num_irreps * CHEMPS2_CASPT2_B_TRIPLET ] );
3463       }
3464       for ( int irrep = 1; irrep < num_irreps; irrep++ ){
3465          assert( size_B_singlet[ irrep ] == size_B_triplet[ irrep ] );
3466          if ( size_B_singlet[ irrep ] > 0 ){
3467             int shift = 0;
3468             for ( int irrep_i = 0; irrep_i < num_irreps; irrep_i++ ){
3469                const int irrep_j = Irreps::directProd( irrep, irrep_i );
3470                if ( irrep_i < irrep_j ){
3471                   assert( shift == shift_B_nonactive( indices, irrep_i, irrep_j, 0 ) );
3472                   const int NOCC_i = indices->getNOCC( irrep_i );
3473                   const int NOCC_j = indices->getNOCC( irrep_j );
3474                   const int loopsize = NOCC_i * NOCC_j;
3475                   #pragma omp for schedule(static)
3476                   for ( int combined = 0; combined < loopsize; combined++ ){
3477 
3478                      const int i = combined % NOCC_i;
3479                      const int j = combined / NOCC_i;
3480 
3481                      // Fill workspace[ xy ] with [ (ix|jy) + (iy|jx) ]
3482                      int jump_xy = 0;
3483                      for ( int irrep_x = 0; irrep_x < num_irreps; irrep_x++ ){
3484                         const int irrep_y = Irreps::directProd( irrep, irrep_x );
3485                         if ( irrep_x < irrep_y ){
3486                            const int occ_x = indices->getNOCC( irrep_x );
3487                            const int occ_y = indices->getNOCC( irrep_y );
3488                            const int num_x = indices->getNDMRG( irrep_x );
3489                            const int num_y = indices->getNDMRG( irrep_y );
3490 
3491                            for ( int y = 0; y < num_y; y++ ){
3492                               for ( int x = 0; x < num_x; x++ ){
3493                                  const double ix_jy = integrals->get_coulomb( irrep_i, irrep_x, irrep_j, irrep_y, i, occ_x + x, j, occ_y + y );
3494                                  const double iy_jx = integrals->get_coulomb( irrep_i, irrep_y, irrep_j, irrep_x, i, occ_y + y, j, occ_x + x );
3495                                  workspace[ jump_xy + x + num_x * y ] = ix_jy + iy_jx;
3496                               }
3497                            }
3498                            jump_xy += num_x * num_y;
3499                         }
3500                      }
3501                      assert( jump_xy == size_B_singlet[ irrep ] );
3502 
3503                      // Perform target[ tu ] = sum_{x<y} [ (ix|jy) + (iy|jx) ] SBB_singlet[ It x Iu ][ xytu ]
3504                      char notrans = 'N';
3505                      int inc1 = 1;
3506                      double alpha = 1.0;
3507                      double set = 0.0;
3508                      double * target = vector_rhs + jump[ irrep + num_irreps * CHEMPS2_CASPT2_B_SINGLET ] + size_B_singlet[ irrep ] * ( shift + combined );
3509                      dgemv_( &notrans, &jump_xy, &jump_xy, &alpha, SBB_singlet[ irrep ], &jump_xy, workspace, &inc1, &set, target, &inc1 );
3510 
3511                      // Fill workspace[ xy ] with [ (ix|jy) - (iy|jx) ]
3512                      jump_xy = 0;
3513                      for ( int irrep_x = 0; irrep_x < num_irreps; irrep_x++ ){
3514                         const int irrep_y = Irreps::directProd( irrep, irrep_x );
3515                         if ( irrep_x < irrep_y ){
3516                            const int occ_x = indices->getNOCC( irrep_x );
3517                            const int occ_y = indices->getNOCC( irrep_y );
3518                            const int num_x = indices->getNDMRG( irrep_x );
3519                            const int num_y = indices->getNDMRG( irrep_y );
3520 
3521                            for ( int y = 0; y < num_y; y++ ){
3522                               for ( int x = 0; x < num_x; x++ ){
3523                                  const double ix_jy = integrals->get_coulomb( irrep_i, irrep_x, irrep_j, irrep_y, i, occ_x + x, j, occ_y + y );
3524                                  const double iy_jx = integrals->get_coulomb( irrep_i, irrep_y, irrep_j, irrep_x, i, occ_y + y, j, occ_x + x );
3525                                  workspace[ jump_xy + x + num_x * y ] = ix_jy - iy_jx;
3526                               }
3527                            }
3528                            jump_xy += num_x * num_y;
3529                         }
3530                      }
3531                      assert( jump_xy == size_B_triplet[ irrep ] );
3532 
3533                      // Perform target[ tu ] = sum_{x<y} [ (ix|jy) - (iy|jx) ] SBB_triplet[ It x Iu ][ xytu ]
3534                      target = vector_rhs + jump[ irrep + num_irreps * CHEMPS2_CASPT2_B_TRIPLET ] + size_B_triplet[ irrep ] * ( shift + combined );
3535                      dgemv_( &notrans, &jump_xy, &jump_xy, &alpha, SBB_triplet[ irrep ], &jump_xy, workspace, &inc1, &set, target, &inc1 );
3536                   }
3537                   shift += loopsize;
3538                }
3539             }
3540             assert( shift * size_B_singlet[ irrep ] == jump[ 1 + irrep + num_irreps * CHEMPS2_CASPT2_B_SINGLET ] - jump[ irrep + num_irreps * CHEMPS2_CASPT2_B_SINGLET ] );
3541             assert( shift * size_B_triplet[ irrep ] == jump[ 1 + irrep + num_irreps * CHEMPS2_CASPT2_B_TRIPLET ] - jump[ irrep + num_irreps * CHEMPS2_CASPT2_B_TRIPLET ] );
3542          }
3543       }
3544 
3545       // VF singlet and triplet
3546       if ( size_F_singlet[ 0 ] > 0 ){ // First do irrep == Ii x Ij == Ix x Iy == It x Iu == 0
3547          int shift = 0; // First do SINGLET
3548          for ( int irrep_ab = 0; irrep_ab < num_irreps; irrep_ab++ ){
3549             assert( shift == shift_F_nonactive( indices, irrep_ab, irrep_ab, +1 ) );
3550             const int N_OA_ab = indices->getNOCC( irrep_ab ) + indices->getNDMRG( irrep_ab );
3551             const int NVIR_ab = indices->getNVIRT( irrep_ab );
3552             const int loopsize = ( NVIR_ab * ( NVIR_ab + 1 ) ) / 2;
3553             #pragma omp for schedule(static)
3554             for ( int combined = 0; combined < loopsize; combined++ ){
3555 
3556                Special::invert_triangle_two( combined, triangle_idx );
3557                const int a = triangle_idx[ 0 ];
3558                const int b = triangle_idx[ 1 ];
3559                assert( combined == a + ( b * ( b + 1 ) ) / 2 );
3560 
3561                // Fill workspace[ xy ] with [ (ax|by) + (ay|bx) ] * (( x==y ) ? 0.5 : 1.0 )
3562                int jump_xy = 0;
3563                for ( int irrep_xy = 0; irrep_xy < num_irreps; irrep_xy++ ){
3564                   const int occ_xy = indices->getNOCC( irrep_xy );
3565                   const int num_xy = indices->getNDMRG( irrep_xy );
3566 
3567                   for ( int x = 0; x < num_xy; x++ ){
3568                      for ( int y = x; y < num_xy; y++ ){ // 0 <= x <= y < num_xy
3569                         const double ax_by = integrals->get_exchange( irrep_xy, irrep_xy, irrep_ab, irrep_ab, occ_xy + x, occ_xy + y, N_OA_ab + a, N_OA_ab + b );
3570                         const double ay_bx = integrals->get_exchange( irrep_xy, irrep_xy, irrep_ab, irrep_ab, occ_xy + y, occ_xy + x, N_OA_ab + a, N_OA_ab + b );
3571                         workspace[ jump_xy + x + (y*(y+1))/2 ] = ( ax_by + ay_bx ) * (( x==y ) ? 0.5 : 1.0 );
3572                      }
3573                   }
3574                   jump_xy += ( num_xy * ( num_xy + 1 ) ) / 2;
3575                }
3576                assert( jump_xy == size_F_singlet[ 0 ] );
3577 
3578                // Perform target[ tu ] = sum_{x<=y} [ (ax|by) + (ay|bx) ] SFF_singlet[ It x Iu ][ xytu ] * (( x==y ) ? 0.5 : 1.0 ) / sqrt( 1 + delta_ab )
3579                char notrans = 'N';
3580                int inc1 = 1;
3581                double alpha = (( a == b ) ? SQRT_0p5 : 1.0 );
3582                double set = 0.0;
3583                double * target = vector_rhs + jump[ 0 + num_irreps * CHEMPS2_CASPT2_F_SINGLET ] + size_F_singlet[ 0 ] * ( shift + combined );
3584                dgemv_( &notrans, &jump_xy, &jump_xy, &alpha, SFF_singlet[ 0 ], &jump_xy, workspace, &inc1, &set, target, &inc1 );
3585             }
3586             shift += loopsize;
3587          }
3588          assert( shift * size_F_singlet[ 0 ] == jump[ 1 + num_irreps * CHEMPS2_CASPT2_F_SINGLET ] - jump[ num_irreps * CHEMPS2_CASPT2_F_SINGLET ] );
3589       }
3590       if ( size_F_triplet[ 0 ] > 0 ){ // Then do TRIPLET
3591          int shift = 0;
3592          for ( int irrep_ab = 0; irrep_ab < num_irreps; irrep_ab++ ){
3593             assert( shift == shift_F_nonactive( indices, irrep_ab, irrep_ab, -1 ) );
3594             const int N_OA_ab = indices->getNOCC( irrep_ab ) + indices->getNDMRG( irrep_ab );
3595             const int NVIR_ab = indices->getNVIRT( irrep_ab );
3596             const int loopsize = ( NVIR_ab * ( NVIR_ab - 1 ) ) / 2;
3597             #pragma omp for schedule(static)
3598             for ( int combined = 0; combined < loopsize; combined++ ){
3599 
3600                Special::invert_lower_triangle_two( combined, triangle_idx );
3601                const int a = triangle_idx[ 0 ];
3602                const int b = triangle_idx[ 1 ];
3603                assert( combined == a + ( b * ( b - 1 ) ) / 2 );
3604 
3605                // Fill workspace[ xy ] with [ (ax|by) - (ay|bx) ]
3606                int jump_xy = 0;
3607                for ( int irrep_xy = 0; irrep_xy < num_irreps; irrep_xy++ ){
3608                   const int occ_xy = indices->getNOCC( irrep_xy );
3609                   const int num_xy = indices->getNDMRG( irrep_xy );
3610 
3611                   for ( int x = 0; x < num_xy; x++ ){
3612                      for ( int y = x+1; y < num_xy; y++ ){ // 0 <= x < y < num_xy
3613                         const double ax_by = integrals->get_exchange( irrep_xy, irrep_xy, irrep_ab, irrep_ab, occ_xy + x, occ_xy + y, N_OA_ab + a, N_OA_ab + b );
3614                         const double ay_bx = integrals->get_exchange( irrep_xy, irrep_xy, irrep_ab, irrep_ab, occ_xy + y, occ_xy + x, N_OA_ab + a, N_OA_ab + b );
3615                         workspace[ jump_xy + x + (y*(y-1))/2 ] = ax_by - ay_bx;
3616                      }
3617                   }
3618                   jump_xy += ( num_xy * ( num_xy - 1 ) ) / 2;
3619                }
3620                assert( jump_xy == size_F_triplet[ 0 ] );
3621 
3622                // Perform target[ tu ] = sum_{x<y} [ (ax|by) - (ay|bx) ] SFF_triplet[ It x Iu ][ xytu ]
3623                char notrans = 'N';
3624                int inc1 = 1;
3625                double alpha = 1.0;
3626                double set = 0.0;
3627                double * target = vector_rhs + jump[ 0 + num_irreps * CHEMPS2_CASPT2_F_TRIPLET ] + size_F_triplet[ 0 ] * ( shift + combined );
3628                dgemv_( &notrans, &jump_xy, &jump_xy, &alpha, SFF_triplet[ 0 ], &jump_xy, workspace, &inc1, &set, target, &inc1 );
3629             }
3630             shift += loopsize;
3631          }
3632          assert( shift * size_F_triplet[ 0 ] == jump[ 1 + num_irreps * CHEMPS2_CASPT2_F_TRIPLET ] - jump[ num_irreps * CHEMPS2_CASPT2_F_TRIPLET ] );
3633       }
3634       for ( int irrep = 1; irrep < num_irreps; irrep++ ){
3635          assert( size_F_singlet[ irrep ] == size_F_triplet[ irrep ] );
3636          if ( size_F_singlet[ irrep ] > 0 ){
3637             int shift = 0;
3638             for ( int irrep_a = 0; irrep_a < num_irreps; irrep_a++ ){
3639                const int irrep_b = Irreps::directProd( irrep, irrep_a );
3640                if ( irrep_a < irrep_b ){
3641                   assert( shift == shift_F_nonactive( indices, irrep_a, irrep_b, 0 ) );
3642                   const int N_OA_a = indices->getNOCC( irrep_a ) + indices->getNDMRG( irrep_a );
3643                   const int N_OA_b = indices->getNOCC( irrep_b ) + indices->getNDMRG( irrep_b );
3644                   const int NVIR_a = indices->getNVIRT( irrep_a );
3645                   const int NVIR_b = indices->getNVIRT( irrep_b );
3646                   const int loopsize = NVIR_a * NVIR_b;
3647                   #pragma omp for schedule(static)
3648                   for ( int combined = 0; combined < loopsize; combined++ ){
3649 
3650                      const int a = combined % NVIR_a;
3651                      const int b = combined / NVIR_a;
3652 
3653                      // Fill workspace[ xy ] with [ (ax|by) + (ay|bx) ]
3654                      int jump_xy = 0;
3655                      for ( int irrep_x = 0; irrep_x < num_irreps; irrep_x++ ){
3656                         const int irrep_y = Irreps::directProd( irrep, irrep_x );
3657                         if ( irrep_x < irrep_y ){
3658                            const int occ_x = indices->getNOCC( irrep_x );
3659                            const int occ_y = indices->getNOCC( irrep_y );
3660                            const int num_x = indices->getNDMRG( irrep_x );
3661                            const int num_y = indices->getNDMRG( irrep_y );
3662 
3663                            for ( int y = 0; y < num_y; y++ ){
3664                               for ( int x = 0; x < num_x; x++ ){
3665                                  const double ax_by = integrals->get_exchange( irrep_x, irrep_y, irrep_a, irrep_b, occ_x + x, occ_y + y, N_OA_a + a, N_OA_b + b );
3666                                  const double ay_bx = integrals->get_exchange( irrep_y, irrep_x, irrep_a, irrep_b, occ_y + y, occ_x + x, N_OA_a + a, N_OA_b + b );
3667                                  workspace[ jump_xy + x + num_x * y ] = ax_by + ay_bx;
3668                               }
3669                            }
3670                            jump_xy += num_x * num_y;
3671                         }
3672                      }
3673                      assert( jump_xy == size_F_singlet[ irrep ] );
3674 
3675                      // Perform target[ tu ] = sum_{x<y} [ (ax|by) + (ay|bx) ] SFF_singlet[ It x Iu ][ xytu ]
3676                      char notrans = 'N';
3677                      int inc1 = 1;
3678                      double alpha = 1.0;
3679                      double set = 0.0;
3680                      double * target = vector_rhs + jump[ irrep + num_irreps * CHEMPS2_CASPT2_F_SINGLET ] + size_F_singlet[ irrep ] * ( shift + combined );
3681                      dgemv_( &notrans, &jump_xy, &jump_xy, &alpha, SFF_singlet[ irrep ], &jump_xy, workspace, &inc1, &set, target, &inc1 );
3682 
3683                      // Fill workspace[ xy ] with [ (ax|by) - (ay|bx) ]
3684                      jump_xy = 0;
3685                      for ( int irrep_x = 0; irrep_x < num_irreps; irrep_x++ ){
3686                         const int irrep_y = Irreps::directProd( irrep, irrep_x );
3687                         if ( irrep_x < irrep_y ){
3688                            const int occ_x = indices->getNOCC( irrep_x );
3689                            const int occ_y = indices->getNOCC( irrep_y );
3690                            const int num_x = indices->getNDMRG( irrep_x );
3691                            const int num_y = indices->getNDMRG( irrep_y );
3692 
3693                            for ( int y = 0; y < num_y; y++ ){
3694                               for ( int x = 0; x < num_x; x++ ){
3695                                  const double ax_by = integrals->get_exchange( irrep_x, irrep_y, irrep_a, irrep_b, occ_x + x, occ_y + y, N_OA_a + a, N_OA_b + b );
3696                                  const double ay_bx = integrals->get_exchange( irrep_y, irrep_x, irrep_a, irrep_b, occ_y + y, occ_x + x, N_OA_a + a, N_OA_b + b );
3697                                  workspace[ jump_xy + x + num_x * y ] = ax_by - ay_bx;
3698                               }
3699                            }
3700                            jump_xy += num_x * num_y;
3701                         }
3702                      }
3703                      assert( jump_xy == size_F_triplet[ irrep ] );
3704 
3705                      // Perform target[ tu ] = sum_{x<y} [ (ax|by) - (ay|bx) ] SFF_triplet[ It x Iu ][ xytu ]
3706                      target = vector_rhs + jump[ irrep + num_irreps * CHEMPS2_CASPT2_F_TRIPLET ] + size_F_triplet[ irrep ] * ( shift + combined );
3707                      dgemv_( &notrans, &jump_xy, &jump_xy, &alpha, SFF_triplet[ irrep ], &jump_xy, workspace, &inc1, &set, target, &inc1 );
3708                   }
3709                   shift += loopsize;
3710                }
3711             }
3712             assert( shift * size_F_singlet[ irrep ] == jump[ 1 + irrep + num_irreps * CHEMPS2_CASPT2_F_SINGLET ] - jump[ irrep + num_irreps * CHEMPS2_CASPT2_F_SINGLET ] );
3713             assert( shift * size_F_triplet[ irrep ] == jump[ 1 + irrep + num_irreps * CHEMPS2_CASPT2_F_TRIPLET ] - jump[ irrep + num_irreps * CHEMPS2_CASPT2_F_TRIPLET ] );
3714          }
3715       }
3716       delete [] workspace;
3717 
3718       // VE singlet and triplet
3719       for ( int irrep = 0; irrep < num_irreps; irrep++ ){
3720          const int occ_t = indices->getNOCC( irrep );
3721          const int num_t = indices->getNDMRG( irrep );
3722          if ( num_t > 0 ){
3723             double * target_singlet = vector_rhs + jump[ irrep + num_irreps * CHEMPS2_CASPT2_E_SINGLET ];
3724             double * target_triplet = vector_rhs + jump[ irrep + num_irreps * CHEMPS2_CASPT2_E_TRIPLET ];
3725             int shift_singlet = 0;
3726             int shift_triplet = 0;
3727             for ( int irrep_a = 0; irrep_a < num_irreps; irrep_a++ ){
3728                const int NVIR_a = indices->getNVIRT( irrep_a );
3729                const int N_OA_a = indices->getNOCC( irrep_a ) + indices->getNDMRG( irrep_a );
3730                const int irrep_occ = Irreps::directProd( irrep_a, irrep );
3731                if ( irrep_occ == 0 ){
3732                   for ( int irrep_ij = 0; irrep_ij < num_irreps; irrep_ij++ ){
3733                      assert( shift_singlet == shift_E_nonactive( indices, irrep_a, irrep_ij, irrep_ij, +1 ) );
3734                      assert( shift_triplet == shift_E_nonactive( indices, irrep_a, irrep_ij, irrep_ij, -1 ) );
3735                      const int NOCC_ij = indices->getNOCC( irrep_ij );
3736                      const int loopsize_singlet = ( NVIR_a * NOCC_ij * ( NOCC_ij + 1 ) ) / 2;
3737                      const int loopsize_triplet = ( NVIR_a * NOCC_ij * ( NOCC_ij - 1 ) ) / 2;
3738                      #pragma omp for schedule(static)
3739                      for ( int combined_singlet = 0; combined_singlet < loopsize_singlet; combined_singlet++ ){
3740 
3741                         const int a = combined_singlet % NVIR_a;
3742                         const int remainder = combined_singlet / NVIR_a;
3743                         Special::invert_triangle_two( remainder, triangle_idx );
3744                         const int i = triangle_idx[ 0 ];
3745                         const int j = triangle_idx[ 1 ];
3746                         const int combined_triplet = a + NVIR_a * ( i + ( j * ( j - 1 ) ) / 2 );
3747                         const double ij_factor = (( i == j ) ? SQRT_0p5 : 1.0 );
3748 
3749                         const int count_aij_singlet = shift_singlet + combined_singlet;
3750                         const int count_aij_triplet = shift_triplet + combined_triplet;
3751                         for ( int t = 0; t < num_t; t++ ){
3752                            double value_singlet = 0.0;
3753                            double value_triplet = 0.0;
3754                            for ( int w = 0; w < num_t; w++ ){
3755                               const double SEE_wt = SEE[ irrep ][ w + num_t * t ];
3756                               const double aj_wi  = integrals->get_coulomb( irrep_ij, irrep, irrep_ij, irrep_a, i, occ_t + w, j, N_OA_a + a );
3757                               const double ai_wj  = integrals->get_coulomb( irrep_ij, irrep, irrep_ij, irrep_a, j, occ_t + w, i, N_OA_a + a );
3758                               value_singlet +=     SEE_wt * ( aj_wi + ai_wj );
3759                               value_triplet += 3 * SEE_wt * ( aj_wi - ai_wj );
3760                            }
3761                            target_singlet[ t + num_t * count_aij_singlet ] = value_singlet * ij_factor;
3762              if ( j > i ){ target_triplet[ t + num_t * count_aij_triplet ] = value_triplet; }
3763                         }
3764                      }
3765                      shift_singlet += loopsize_singlet;
3766                      shift_triplet += loopsize_triplet;
3767                   }
3768                } else {
3769                   for ( int irrep_i = 0; irrep_i < num_irreps; irrep_i++ ){
3770                      const int irrep_j = Irreps::directProd( irrep_i, irrep_occ );
3771                      if ( irrep_i < irrep_j ){
3772                         assert( shift_singlet == shift_E_nonactive( indices, irrep_a, irrep_i, irrep_j, +1 ) );
3773                         assert( shift_triplet == shift_E_nonactive( indices, irrep_a, irrep_i, irrep_j, -1 ) );
3774                         const int NOCC_i = indices->getNOCC( irrep_i );
3775                         const int NOCC_j = indices->getNOCC( irrep_j );
3776                         const int loopsize = NOCC_i * NOCC_j * NVIR_a;
3777                         #pragma omp for schedule(static)
3778                         for ( int combined = 0; combined < loopsize; combined++ ){
3779 
3780                            const int a = combined % NVIR_a;
3781                            const int remainder = combined / NVIR_a;
3782                            const int i = remainder % NOCC_i;
3783                            const int j = remainder / NOCC_i;
3784 
3785                            const int count_aij_singlet = shift_singlet + combined;
3786                            const int count_aij_triplet = shift_triplet + combined;
3787                            for ( int t = 0; t < num_t; t++ ){
3788                               double value_singlet = 0.0;
3789                               double value_triplet = 0.0;
3790                               for ( int w = 0; w < num_t; w++ ){
3791                                  const double SEE_wt = SEE[ irrep ][ w + num_t * t ];
3792                                  const double aj_wi  = integrals->get_coulomb( irrep_i, irrep, irrep_j, irrep_a, i, occ_t + w, j, N_OA_a + a );
3793                                  const double ai_wj  = integrals->get_coulomb( irrep_j, irrep, irrep_i, irrep_a, j, occ_t + w, i, N_OA_a + a );
3794                                  value_singlet +=     SEE_wt * ( aj_wi + ai_wj );
3795                                  value_triplet += 3 * SEE_wt * ( aj_wi - ai_wj );
3796                               }
3797                               target_singlet[ t + num_t * count_aij_singlet ] = value_singlet;
3798                               target_triplet[ t + num_t * count_aij_triplet ] = value_triplet;
3799                            }
3800                         }
3801                         shift_singlet += loopsize;
3802                         shift_triplet += loopsize;
3803                      }
3804                   }
3805                }
3806             }
3807             assert( shift_singlet * num_t == jump[ 1 + irrep + num_irreps * CHEMPS2_CASPT2_E_SINGLET ] - jump[ irrep + num_irreps * CHEMPS2_CASPT2_E_SINGLET ] );
3808             assert( shift_triplet * num_t == jump[ 1 + irrep + num_irreps * CHEMPS2_CASPT2_E_TRIPLET ] - jump[ irrep + num_irreps * CHEMPS2_CASPT2_E_TRIPLET ] );
3809          }
3810       }
3811 
3812       // VG singlet and triplet
3813       for ( int irrep = 0; irrep < num_irreps; irrep++ ){
3814          const int occ_t = indices->getNOCC( irrep );
3815          const int num_t = indices->getNDMRG( irrep );
3816          if ( num_t > 0 ){
3817             double * target_singlet = vector_rhs + jump[ irrep + num_irreps * CHEMPS2_CASPT2_G_SINGLET ];
3818             double * target_triplet = vector_rhs + jump[ irrep + num_irreps * CHEMPS2_CASPT2_G_TRIPLET ];
3819             int shift_singlet = 0;
3820             int shift_triplet = 0;
3821             for ( int irrep_i = 0; irrep_i < num_irreps; irrep_i++ ){
3822                const int NOCC_i = indices->getNOCC( irrep_i );
3823                const int irrep_virt = Irreps::directProd( irrep_i, irrep );
3824                if ( irrep_virt == 0 ){
3825                   for ( int irrep_ab = 0; irrep_ab < num_irreps; irrep_ab++ ){
3826                      assert( shift_singlet == shift_G_nonactive( indices, irrep_i, irrep_ab, irrep_ab, +1 ) );
3827                      assert( shift_triplet == shift_G_nonactive( indices, irrep_i, irrep_ab, irrep_ab, -1 ) );
3828                      const int N_OA_ab = indices->getNOCC( irrep_ab ) + indices->getNDMRG( irrep_ab );
3829                      const int NVIR_ab = indices->getNVIRT( irrep_ab );
3830                      const int loopsize_singlet = ( NOCC_i * NVIR_ab * ( NVIR_ab + 1 ) ) / 2;
3831                      const int loopsize_triplet = ( NOCC_i * NVIR_ab * ( NVIR_ab - 1 ) ) / 2;
3832                      #pragma omp for schedule(static)
3833                      for ( int combined_singlet = 0; combined_singlet < loopsize_singlet; combined_singlet++ ){
3834 
3835                         const int i = combined_singlet % NOCC_i;
3836                         const int remainder = combined_singlet / NOCC_i;
3837                         Special::invert_triangle_two( remainder, triangle_idx );
3838                         const int a = triangle_idx[ 0 ];
3839                         const int b = triangle_idx[ 1 ];
3840                         const int combined_triplet = i + NOCC_i * ( a + ( b * ( b - 1 ) ) / 2 );
3841                         const double ab_factor = (( a == b ) ? SQRT_0p5 : 1.0 );
3842 
3843                         const int count_abi_singlet = shift_singlet + combined_singlet;
3844                         const int count_abi_triplet = shift_triplet + combined_triplet;
3845                         for ( int t = 0; t < num_t; t++ ){
3846                            double value_singlet = 0.0;
3847                            double value_triplet = 0.0;
3848                            for ( int u = 0; u < num_t; u++ ){
3849                               const double SGG_ut = SGG[ irrep ][ u + num_t * t ];
3850                               const double ai_bu  = integrals->get_exchange( irrep_i, irrep, irrep_ab, irrep_ab, i, occ_t + u, N_OA_ab + a, N_OA_ab + b );
3851                               const double bi_au  = integrals->get_exchange( irrep_i, irrep, irrep_ab, irrep_ab, i, occ_t + u, N_OA_ab + b, N_OA_ab + a );
3852                               value_singlet +=     SGG_ut * ( ai_bu + bi_au );
3853                               value_triplet += 3 * SGG_ut * ( ai_bu - bi_au );
3854                            }
3855                            target_singlet[ t + num_t * count_abi_singlet ] = value_singlet * ab_factor;
3856              if ( b > a ){ target_triplet[ t + num_t * count_abi_triplet ] = value_triplet; }
3857                         }
3858                      }
3859                      shift_singlet += loopsize_singlet;
3860                      shift_triplet += loopsize_triplet;
3861                   }
3862                } else {
3863                   for ( int irrep_a = 0; irrep_a < num_irreps; irrep_a++ ){
3864                      const int irrep_b = Irreps::directProd( irrep_a, irrep_virt );
3865                      if ( irrep_a < irrep_b ){
3866                         assert( shift_singlet == shift_G_nonactive( indices, irrep_i, irrep_a, irrep_b, +1 ) );
3867                         assert( shift_triplet == shift_G_nonactive( indices, irrep_i, irrep_a, irrep_b, -1 ) );
3868                         const int N_OA_a = indices->getNOCC( irrep_a ) + indices->getNDMRG( irrep_a );
3869                         const int N_OA_b = indices->getNOCC( irrep_b ) + indices->getNDMRG( irrep_b );
3870                         const int NVIR_a = indices->getNVIRT( irrep_a );
3871                         const int NVIR_b = indices->getNVIRT( irrep_b );
3872                         const int loopsize = NOCC_i * NVIR_a * NVIR_b;
3873                         #pragma omp for schedule(static)
3874                         for ( int combined = 0; combined < loopsize; combined++ ){
3875 
3876                            const int i = combined % NOCC_i;
3877                            const int remainder = combined / NOCC_i;
3878                            const int a = remainder % NVIR_a;
3879                            const int b = remainder / NVIR_a;
3880 
3881                            const int count_abi_singlet = shift_singlet + combined;
3882                            const int count_abi_triplet = shift_triplet + combined;
3883                            for ( int t = 0; t < num_t; t++ ){
3884                               double value_singlet = 0.0;
3885                               double value_triplet = 0.0;
3886                               for ( int u = 0; u < num_t; u++ ){
3887                                  const double SGG_ut = SGG[ irrep ][ u + num_t * t ];
3888                                  const double ai_bu  = integrals->get_exchange( irrep_i, irrep, irrep_a, irrep_b, i, occ_t + u, N_OA_a + a, N_OA_b + b );
3889                                  const double bi_au  = integrals->get_exchange( irrep_i, irrep, irrep_b, irrep_a, i, occ_t + u, N_OA_b + b, N_OA_a + a );
3890                                  value_singlet +=     SGG_ut * ( ai_bu + bi_au );
3891                                  value_triplet += 3 * SGG_ut * ( ai_bu - bi_au );
3892                               }
3893                               target_singlet[ t + num_t * count_abi_singlet ] = value_singlet;
3894                               target_triplet[ t + num_t * count_abi_triplet ] = value_triplet;
3895                            }
3896                         }
3897                         shift_singlet += loopsize;
3898                         shift_triplet += loopsize;
3899                      }
3900                   }
3901                }
3902             }
3903             assert( shift_singlet * num_t == jump[ 1 + irrep + num_irreps * CHEMPS2_CASPT2_G_SINGLET ] - jump[ irrep + num_irreps * CHEMPS2_CASPT2_G_SINGLET ] );
3904             assert( shift_triplet * num_t == jump[ 1 + irrep + num_irreps * CHEMPS2_CASPT2_G_TRIPLET ] - jump[ irrep + num_irreps * CHEMPS2_CASPT2_G_TRIPLET ] );
3905          }
3906       }
3907 
3908       // VH singlet and triplet
3909       for ( int irrep = 0; irrep < num_irreps; irrep++ ){
3910          int shift_singlet = 0;
3911          int shift_triplet = 0;
3912          double * target_singlet = vector_rhs + jump[ irrep + num_irreps * CHEMPS2_CASPT2_H_SINGLET ];
3913          double * target_triplet = vector_rhs + jump[ irrep + num_irreps * CHEMPS2_CASPT2_H_TRIPLET ];
3914          if ( irrep == 0 ){ // irrep_i == irrep_j  and  irrep_a == irrep_b
3915             for ( int irrep_ij = 0; irrep_ij < num_irreps; irrep_ij++ ){
3916                const int nocc_ij = indices->getNOCC( irrep_ij );
3917                const int linsize_ij_singlet = ( nocc_ij * ( nocc_ij + 1 ) ) / 2;
3918                const int linsize_ij_triplet = ( nocc_ij * ( nocc_ij - 1 ) ) / 2;
3919                for ( int irrep_ab = 0; irrep_ab < num_irreps; irrep_ab++ ){
3920                   assert( shift_singlet == shift_H_nonactive( indices, irrep_ij, irrep_ij, irrep_ab, irrep_ab, +1 ) );
3921                   assert( shift_triplet == shift_H_nonactive( indices, irrep_ij, irrep_ij, irrep_ab, irrep_ab, -1 ) );
3922                   const int nvirt_ab = indices->getNVIRT( irrep_ab );
3923                   const int noa_ab = indices->getNOCC( irrep_ab ) + indices->getNDMRG( irrep_ab );
3924                   const int linsize_ab_singlet = ( nvirt_ab * ( nvirt_ab + 1 ) ) / 2;
3925                   const int linsize_ab_triplet = ( nvirt_ab * ( nvirt_ab - 1 ) ) / 2;
3926                   #pragma omp for schedule(static)
3927                   for ( int combined_ab_singlet = 0; combined_ab_singlet < linsize_ab_singlet; combined_ab_singlet++ ){
3928 
3929                      Special::invert_triangle_two( combined_ab_singlet, triangle_idx );
3930                      const int a = triangle_idx[ 0 ];
3931                      const int b = triangle_idx[ 1 ];
3932                      const double ab_factor = (( a == b ) ? SQRT_0p5 : 1.0 );
3933                      const int combined_ab_triplet = a + ( b * ( b - 1 ) ) / 2;
3934 
3935                      for ( int i = 0; i < nocc_ij; i++ ){
3936                         for ( int j = i; j < nocc_ij; j++ ){
3937                            const double ij_factor = (( i == j ) ? SQRT_0p5 : 1.0 );
3938                            const int counter_singlet = shift_singlet + i + ( j * ( j + 1 ) ) / 2 + linsize_ij_singlet * combined_ab_singlet;
3939                            const int counter_triplet = shift_triplet + i + ( j * ( j - 1 ) ) / 2 + linsize_ij_triplet * combined_ab_triplet;
3940 
3941                            const double ai_bj = integrals->get_exchange( irrep_ij, irrep_ij, irrep_ab, irrep_ab, i, j, noa_ab + a, noa_ab + b );
3942                            const double aj_bi = integrals->get_exchange( irrep_ij, irrep_ij, irrep_ab, irrep_ab, j, i, noa_ab + a, noa_ab + b );
3943                            target_singlet[ counter_singlet ] = 2 * ( ai_bj + aj_bi ) * ij_factor * ab_factor;
3944       if ((a<b) && (i<j)){ target_triplet[ counter_triplet ] = 6 * ( ai_bj - aj_bi ); }
3945                         }
3946                      }
3947                   }
3948                   shift_singlet += linsize_ij_singlet * linsize_ab_singlet;
3949                   shift_triplet += linsize_ij_triplet * linsize_ab_triplet;
3950                }
3951             }
3952          } else { // irrep_i < irrep_j = irrep_i x irrep   and   irrep_a < irrep_b = irrep_a x irrep
3953             for ( int irrep_i = 0; irrep_i < num_irreps; irrep_i++ ){
3954                const int irrep_j = Irreps::directProd( irrep, irrep_i );
3955                if ( irrep_i < irrep_j ){
3956                   const int nocc_i = indices->getNOCC( irrep_i );
3957                   const int nocc_j = indices->getNOCC( irrep_j );
3958                   const int linsize_ij = nocc_i * nocc_j;
3959                   for ( int irrep_a = 0; irrep_a < num_irreps; irrep_a++ ){
3960                      const int irrep_b = Irreps::directProd( irrep, irrep_a );
3961                      if ( irrep_a < irrep_b ){
3962                         assert( shift_singlet == shift_H_nonactive( indices, irrep_i, irrep_j, irrep_a, irrep_b, +1 ) );
3963                         assert( shift_triplet == shift_H_nonactive( indices, irrep_i, irrep_j, irrep_a, irrep_b, -1 ) );
3964                         const int nvir_a = indices->getNVIRT( irrep_a );
3965                         const int nvir_b = indices->getNVIRT( irrep_b );
3966                         const int noa_a  = indices->getNOCC( irrep_a ) + indices->getNDMRG( irrep_a );
3967                         const int noa_b  = indices->getNOCC( irrep_b ) + indices->getNDMRG( irrep_b );
3968                         const int linsize_ab = nvir_a * nvir_b;
3969                         #pragma omp for schedule(static)
3970                         for ( int combined_ab = 0; combined_ab < linsize_ab; combined_ab++ ){
3971 
3972                            const int a = combined_ab % nvir_a;
3973                            const int b = combined_ab / nvir_a;
3974 
3975                            for ( int j = 0; j < nocc_j; j++ ){
3976                               for ( int i = 0; i < nocc_i; i++ ){
3977                                  const int count_singlet = shift_singlet + i + nocc_i * ( j + nocc_j * combined_ab );
3978                                  const int count_triplet = shift_triplet + i + nocc_i * ( j + nocc_j * combined_ab );
3979                                  const double ai_bj = integrals->get_exchange( irrep_i, irrep_j, irrep_a, irrep_b, i, j, noa_a + a, noa_b + b );
3980                                  const double aj_bi = integrals->get_exchange( irrep_j, irrep_i, irrep_a, irrep_b, j, i, noa_a + a, noa_b + b );
3981                                  target_singlet[ count_singlet ] = 2 * ( ai_bj + aj_bi );
3982                                  target_triplet[ count_triplet ] = 6 * ( ai_bj - aj_bi );
3983                               }
3984                            }
3985                         }
3986                         shift_singlet += linsize_ij * linsize_ab;
3987                         shift_triplet += linsize_ij * linsize_ab;
3988                      }
3989                   }
3990                }
3991             }
3992          }
3993          assert( shift_singlet == jump[ 1 + irrep + num_irreps * CHEMPS2_CASPT2_H_SINGLET ] - jump[ irrep + num_irreps * CHEMPS2_CASPT2_H_SINGLET ] );
3994          assert( shift_triplet == jump[ 1 + irrep + num_irreps * CHEMPS2_CASPT2_H_TRIPLET ] - jump[ irrep + num_irreps * CHEMPS2_CASPT2_H_TRIPLET ] );
3995       }
3996    }// End of #pragma omp parallel
3997 
3998 }
3999 
jump_AC_active(const DMRGSCFindices * idx,const int irrep_t,const int irrep_u,const int irrep_v)4000 int CheMPS2::CASPT2::jump_AC_active( const DMRGSCFindices * idx, const int irrep_t, const int irrep_u, const int irrep_v ){
4001 
4002    const int irrep_tuv = Irreps::directProd( irrep_t, Irreps::directProd( irrep_u, irrep_v ) );
4003    const int n_irreps = idx->getNirreps();
4004 
4005    int jump_ac = 0;
4006    for ( int It = 0; It < n_irreps; It++ ){
4007       for ( int Iu = 0; Iu < n_irreps; Iu++ ){
4008          const int Iv = Irreps::directProd( irrep_tuv, Irreps::directProd( It, Iu ) );
4009          if (( It == irrep_t ) && ( Iu == irrep_u ) && ( Iv == irrep_v )){
4010             It = n_irreps;
4011             Iu = n_irreps;
4012          } else {
4013             jump_ac += idx->getNDMRG( It ) * idx->getNDMRG( Iu ) * idx->getNDMRG( Iv );
4014          }
4015       }
4016    }
4017    return jump_ac;
4018 
4019 }
4020 
jump_BF_active(const DMRGSCFindices * idx,const int irrep_t,const int irrep_u,const int ST)4021 int CheMPS2::CASPT2::jump_BF_active( const DMRGSCFindices * idx, const int irrep_t, const int irrep_u, const int ST ){
4022 
4023    assert( irrep_t <= irrep_u );
4024    const int irrep_tu = Irreps::directProd( irrep_t, irrep_u );
4025    const int n_irreps = idx->getNirreps();
4026 
4027    int jump_bf = 0;
4028    if ( irrep_tu == 0 ){
4029       for ( int Itu = 0; Itu < n_irreps; Itu++ ){
4030          if (( Itu == irrep_t ) && ( Itu == irrep_u )){
4031             Itu = n_irreps;
4032          } else {
4033             jump_bf += ( idx->getNDMRG( Itu ) * ( idx->getNDMRG( Itu ) + ST ) ) / 2;
4034          }
4035       }
4036    } else {
4037       for ( int It = 0; It < n_irreps; It++ ){
4038          const int Iu = Irreps::directProd( irrep_tu, It );
4039          if ( It < Iu ){
4040             if (( It == irrep_t ) && ( Iu == irrep_u )){
4041                It = n_irreps;
4042             } else {
4043                jump_bf += idx->getNDMRG( It ) * idx->getNDMRG( Iu );
4044             }
4045          }
4046       }
4047    }
4048    return jump_bf;
4049 
4050 }
4051 
make_FDE_FDG()4052 void CheMPS2::CASPT2::make_FDE_FDG(){
4053 
4054    /*
4055       FD1E singlet: < E_yx E_lb E_kw SE_tiaj > = 1 delta_ab ( delta_ik delta_jl + delta_il delta_jk ) / sqrt( 1 + delta_ij ) FD1E_singlet[ Ib x Il ][ It ][ w ][ xy, t ]
4056       FD2E singlet: < E_yb E_lx E_kw SE_tiaj > = 1 delta_ab ( delta_ik delta_jl + delta_il delta_jk ) / sqrt( 1 + delta_ij ) FD2E_singlet[ Ib x Il ][ It ][ w ][ xy, t ]
4057       FD1E triplet: < E_yx E_lb E_kw TE_tiaj > = 3 delta_ab ( delta_ik delta_jl - delta_il delta_jk ) / sqrt( 1 + delta_ij ) FD1E_triplet[ Ib x Il ][ It ][ w ][ xy, t ]
4058       FD2E triplet: < E_yb E_lx E_kw TE_tiaj > = 3 delta_ab ( delta_ik delta_jl - delta_il delta_jk ) / sqrt( 1 + delta_ij ) FD2E_triplet[ Ib x Il ][ It ][ w ][ xy, t ]
4059 
4060             FD1E_singlet[ Ib x Il ][ It ][ w ][ xy, t ] = ( + 2 delta_tw Gamma_yx
4061                                                             - delta_tx Gamma_yw
4062                                                             - Gamma_ytxw
4063                                                           )
4064 
4065             FD2E_singlet[ Ib x Il ][ It ][ w ][ xy, t ] = ( + Gamma_ytxw
4066                                                             + Gamma_ytwx
4067                                                             - delta_tx Gamma_yw
4068                                                             - delta_tw Gamma_yx
4069                                                           )
4070 
4071             FD1E_triplet[ Ib x Il ][ It ][ w ][ xy, t ] = ( + 2 delta_tw Gamma_yx
4072                                                             - delta_tx Gamma_yw
4073                                                             - Gamma_ytxw
4074                                                           )
4075 
4076             FD2E_triplet[ Ib x Il ][ It ][ w ][ xy, t ] = ( + Gamma_ytxw / 3
4077                                                             - Gamma_ytwx / 3
4078                                                             + delta_tx Gamma_yw
4079                                                             - delta_tw Gamma_yx
4080                                                           )
4081 
4082       FD1G singlet: < E_yx E_jd E_wc SG_aibt > = 1 delta_ij ( delta_ac delta_bd + delta_ad delta_bc ) / sqrt( 1 + delta_ab ) FD1G_singlet[ Ij x Id ][ It ][ w ][ xy, t ]
4083       FD2G singlet: < E_yd E_jx E_wc SG_aibt > = 1 delta_ij ( delta_ac delta_bd + delta_ad delta_bc ) / sqrt( 1 + delta_ab ) FD2G_singlet[ Ij x Id ][ It ][ w ][ xy, t ]
4084       FD1G triplet: < E_yx E_jd E_wc TG_aibt > = 3 delta_ij ( delta_ac delta_bd - delta_ad delta_bc ) / sqrt( 1 + delta_ab ) FD1G_triplet[ Ij x Id ][ It ][ w ][ xy, t ]
4085       FD2G triplet: < E_yd E_jx E_wc TG_aibt > = 3 delta_ij ( delta_ac delta_bd - delta_ad delta_bc ) / sqrt( 1 + delta_ab ) FD2G_triplet[ Ij x Id ][ It ][ w ][ xy, t ]
4086 
4087             FD1G_singlet[ Ij x Id ][ It ][ w ][ xy, t ] = ( + Gamma_ywxt
4088                                                             + delta_wx Gamma_yt
4089                                                           )
4090 
4091             FD2G_singlet[ Ij x Id ][ It ][ w ][ xy, t ] = ( + delta_xw Gamma_yt
4092                                                             - Gamma_ywtx
4093                                                             - Gamma_ywxt
4094                                                           )
4095 
4096             FD1G_triplet[ Ij x Id ][ It ][ w ][ xy, t ] = ( - Gamma_ywxt
4097                                                             - delta_wx Gamma_yt
4098                                                           )
4099 
4100             FD2G_triplet[ Ij x Id ][ It ][ w ][ xy, t ] = ( + delta_xw Gamma_yt
4101                                                             - Gamma ywtx / 3
4102                                                             + Gamma ywxt / 3
4103                                                           )
4104    */
4105 
4106    FDE_singlet = new double***[ num_irreps ];
4107    FDE_triplet = new double***[ num_irreps ];
4108    FDG_singlet = new double***[ num_irreps ];
4109    FDG_triplet = new double***[ num_irreps ];
4110 
4111    const int LAS = indices->getDMRGcumulative( num_irreps );
4112 
4113    for ( int irrep_left = 0; irrep_left < num_irreps; irrep_left++ ){
4114 
4115       const int SIZE_left = size_D[ irrep_left ];
4116       const int D2JUMP    = SIZE_left / 2;
4117       FDE_singlet[ irrep_left ] = new double**[ num_irreps ];
4118       FDE_triplet[ irrep_left ] = new double**[ num_irreps ];
4119       FDG_singlet[ irrep_left ] = new double**[ num_irreps ];
4120       FDG_triplet[ irrep_left ] = new double**[ num_irreps ];
4121 
4122       for ( int irrep_t = 0; irrep_t < num_irreps; irrep_t++ ){
4123 
4124          const int SIZE_right = indices->getNDMRG( irrep_t );
4125          const int d_t     = indices->getDMRGcumulative( irrep_t );
4126          const int num_t   = indices->getNDMRG( irrep_t );
4127          const int irrep_w = Irreps::directProd( irrep_left, irrep_t );
4128          const int d_w     = indices->getDMRGcumulative( irrep_w );
4129          const int num_w   = indices->getNDMRG( irrep_w );
4130          FDE_singlet[ irrep_left ][ irrep_t ] = new double*[ num_w ];
4131          FDE_triplet[ irrep_left ][ irrep_t ] = new double*[ num_w ];
4132          FDG_singlet[ irrep_left ][ irrep_t ] = new double*[ num_w ];
4133          FDG_triplet[ irrep_left ][ irrep_t ] = new double*[ num_w ];
4134 
4135          for ( int w = 0; w < num_w; w++ ){
4136 
4137             FDE_singlet[ irrep_left ][ irrep_t ][ w ] = new double[ SIZE_left * SIZE_right ];
4138             FDE_triplet[ irrep_left ][ irrep_t ][ w ] = new double[ SIZE_left * SIZE_right ];
4139             FDG_singlet[ irrep_left ][ irrep_t ][ w ] = new double[ SIZE_left * SIZE_right ];
4140             FDG_triplet[ irrep_left ][ irrep_t ][ w ] = new double[ SIZE_left * SIZE_right ];
4141 
4142             double * FDE_sing = FDE_singlet[ irrep_left ][ irrep_t ][ w ];
4143             double * FDE_trip = FDE_triplet[ irrep_left ][ irrep_t ][ w ];
4144             double * FDG_sing = FDG_singlet[ irrep_left ][ irrep_t ][ w ];
4145             double * FDG_trip = FDG_triplet[ irrep_left ][ irrep_t ][ w ];
4146 
4147             int jump_row = 0;
4148             for ( int irrep_x = 0; irrep_x < num_irreps; irrep_x++ ){
4149                const int d_x     = indices->getDMRGcumulative( irrep_x );
4150                const int num_x   = indices->getNDMRG( irrep_x );
4151                const int irrep_y = Irreps::directProd( irrep_left, irrep_x );
4152                const int d_y     = indices->getDMRGcumulative( irrep_y );
4153                const int num_y   = indices->getNDMRG( irrep_y );
4154 
4155                for ( int t = 0; t < num_t; t++ ){
4156                   for ( int y = 0; y < num_y; y++ ){
4157                      for ( int x = 0; x < num_x; x++ ){
4158                         const double gamma_ytxw = two_rdm[ d_y + y + LAS * ( d_t + t + LAS * ( d_x + x + LAS * ( d_w + w ))) ];
4159                         const double gamma_ytwx = two_rdm[ d_y + y + LAS * ( d_t + t + LAS * ( d_w + w + LAS * ( d_x + x ))) ];
4160                         FDE_sing[          jump_row + x + num_x * y + SIZE_left * t ] = - gamma_ytxw;
4161                         FDE_sing[ D2JUMP + jump_row + x + num_x * y + SIZE_left * t ] = + gamma_ytxw + gamma_ytwx;
4162                         FDE_trip[          jump_row + x + num_x * y + SIZE_left * t ] = - gamma_ytxw;
4163                         FDE_trip[ D2JUMP + jump_row + x + num_x * y + SIZE_left * t ] = ( gamma_ytxw - gamma_ytwx ) / 3.0;
4164                         const double gamma_ywtx = two_rdm[ d_y + y + LAS * ( d_w + w + LAS * ( d_t + t + LAS * ( d_x + x ))) ];
4165                         const double gamma_ywxt = two_rdm[ d_y + y + LAS * ( d_w + w + LAS * ( d_x + x + LAS * ( d_t + t ))) ];
4166                         FDG_sing[          jump_row + x + num_x * y + SIZE_left * t ] = + gamma_ywxt;
4167                         FDG_sing[ D2JUMP + jump_row + x + num_x * y + SIZE_left * t ] = - gamma_ywxt - gamma_ywtx;
4168                         FDG_trip[          jump_row + x + num_x * y + SIZE_left * t ] = - gamma_ywxt;
4169                         FDG_trip[ D2JUMP + jump_row + x + num_x * y + SIZE_left * t ] = ( gamma_ywxt - gamma_ywtx ) / 3.0;
4170                      }
4171                   }
4172                }
4173 
4174                if (( irrep_t == irrep_w ) && ( irrep_x == irrep_y )){
4175                   for ( int x = 0; x < num_x; x++ ){
4176                      for ( int y = 0; y < num_x; y++ ){
4177                         const double gamma_yx = one_rdm[ d_y + y + LAS * ( d_x + x ) ];
4178                         FDE_sing[          jump_row + x + num_x * y + SIZE_left * w ] += 2 * gamma_yx;
4179                         FDE_sing[ D2JUMP + jump_row + x + num_x * y + SIZE_left * w ] -= gamma_yx;
4180                         FDE_trip[          jump_row + x + num_x * y + SIZE_left * w ] += 2 * gamma_yx;
4181                         FDE_trip[ D2JUMP + jump_row + x + num_x * y + SIZE_left * w ] -= gamma_yx;
4182                      }
4183                   }
4184                }
4185 
4186                if (( irrep_t == irrep_x ) && ( irrep_y == irrep_w )){
4187                   for ( int y = 0; y < num_y; y++ ){
4188                      const double gamma_yw = one_rdm[ d_y + y + LAS * ( d_w + w ) ];
4189                      for ( int tx = 0; tx < num_x; tx++ ){
4190                         FDE_sing[          jump_row + tx + num_x * y + SIZE_left * tx ] -= gamma_yw;
4191                         FDE_sing[ D2JUMP + jump_row + tx + num_x * y + SIZE_left * tx ] -= gamma_yw;
4192                         FDE_trip[          jump_row + tx + num_x * y + SIZE_left * tx ] -= gamma_yw;
4193                         FDE_trip[ D2JUMP + jump_row + tx + num_x * y + SIZE_left * tx ] += gamma_yw;
4194                      }
4195                   }
4196                }
4197 
4198                if (( irrep_t == irrep_y ) && ( irrep_w == irrep_x )){
4199                   for ( int t = 0; t < num_t; t++ ){
4200                      for ( int y = 0; y < num_y; y++ ){
4201                         const double gamma_yt = one_rdm[ d_y + y + LAS * ( d_t + t ) ];
4202                         FDG_sing[          jump_row + w + num_x * y + SIZE_left * t ] += gamma_yt;
4203                         FDG_sing[ D2JUMP + jump_row + w + num_x * y + SIZE_left * t ] += gamma_yt;
4204                         FDG_trip[          jump_row + w + num_x * y + SIZE_left * t ] -= gamma_yt;
4205                         FDG_trip[ D2JUMP + jump_row + w + num_x * y + SIZE_left * t ] += gamma_yt;
4206                      }
4207                   }
4208                }
4209                jump_row += num_x * num_y;
4210             }
4211          }
4212       }
4213    }
4214 
4215 }
4216 
make_FEH_FGH()4217 void CheMPS2::CASPT2::make_FEH_FGH(){
4218 
4219    /*
4220       FEH singlet: < SE_xkdl E_wc SH_aibj > = 2 delta_ik delta_jl ( delta_ac delta_bd + delta_ad delta_bc ) / sqrt( 1 + delta_ab ) FEH[ Ix ][ w ][ x ]
4221       FEH triplet: < TE_xkdl E_wc TH_aibj > = 6 delta_ik delta_jl ( delta_ac delta_bd - delta_ad delta_bc ) / sqrt( 1 + delta_ab ) FEH[ Ix ][ w ][ x ]
4222       FGH singlet: < SG_cldx E_kw SH_aibj > = 2 delta_ac delta_bd ( delta_il delta_jk + delta_ik delta_jl ) / sqrt( 1 + delta_ij ) FGH[ Ix ][ w ][ x ]
4223       FGH triplet: < TG_cldx E_kw TH_aibj > = 6 delta_ac delta_bd ( delta_il delta_jk - delta_ik delta_jl ) / sqrt( 1 + delta_ij ) FGH[ Ix ][ w ][ x ]
4224 
4225             FEH[ Ix ][ w ][ x ] = + SEE[ Ix ][ xw ]
4226             FGH[ Ix ][ w ][ x ] = - SGG[ Ix ][ xw ]
4227    */
4228 
4229    FEH = new double**[ num_irreps ];
4230    FGH = new double**[ num_irreps ];
4231 
4232    for ( int irrep_x = 0; irrep_x < num_irreps; irrep_x++ ){
4233 
4234       const int NACT = indices->getNDMRG( irrep_x );
4235       FEH[ irrep_x ] = new double*[ NACT ];
4236       FGH[ irrep_x ] = new double*[ NACT ];
4237 
4238       for ( int w = 0; w < NACT; w++ ){
4239          FEH[ irrep_x ][ w ] = new double[ NACT ];
4240          FGH[ irrep_x ][ w ] = new double[ NACT ];
4241 
4242          for ( int x = 0; x < NACT; x++ ){
4243             FEH[ irrep_x ][ w ][ x ] = + SEE[ irrep_x ][ x + NACT * w ];
4244             FGH[ irrep_x ][ w ][ x ] = - SGG[ irrep_x ][ x + NACT * w ];
4245          }
4246       }
4247    }
4248 
4249 }
4250 
make_FBE_FFG_singlet()4251 void CheMPS2::CASPT2::make_FBE_FFG_singlet(){
4252 
4253    /*
4254       FBE singlet : < SB_xkyl E_wc SE_tiaj > = 2 delta_ac delta_ik delta_jl FBE_singlet[ Ik x Il ][ It ][ w ][ xy, t ]
4255       FFG singlet : < SF_cxdy E_kw SG_aibt > = 2 delta_ac delta_bd delta_ik FFG_singlet[ Ic x Id ][ It ][ w ][ xy, t ]
4256 
4257             FBE_singlet[ Ik x Il ][ It ][ w ][ xy, t ] = + SBB_singlet[ Ik x Il ][ xy, tw ]
4258             FFG_singlet[ Ic x Id ][ It ][ w ][ xy, t ] = - SFF_singlet[ Ic x Id ][ xy, tw ]
4259    */
4260 
4261    FBE_singlet = new double***[ num_irreps ];
4262    FFG_singlet = new double***[ num_irreps ];
4263 
4264    for ( int irrep_left = 0; irrep_left < num_irreps; irrep_left++ ){
4265 
4266       assert( size_B_singlet[ irrep_left ] == size_F_singlet[ irrep_left ] ); // At construction
4267       const int SIZE_left = size_B_singlet[ irrep_left ];
4268       FBE_singlet[ irrep_left ] = new double**[ num_irreps ];
4269       FFG_singlet[ irrep_left ] = new double**[ num_irreps ];
4270 
4271       for ( int irrep_t = 0; irrep_t < num_irreps; irrep_t++ ){
4272 
4273          const int SIZE_right = indices->getNDMRG( irrep_t );
4274          const int num_t   = indices->getNDMRG( irrep_t );
4275          const int irrep_w = Irreps::directProd( irrep_left, irrep_t );
4276          const int num_w   = indices->getNDMRG( irrep_w );
4277          FBE_singlet[ irrep_left ][ irrep_t ] = new double*[ num_w ];
4278          FFG_singlet[ irrep_left ][ irrep_t ] = new double*[ num_w ];
4279 
4280          for ( int w = 0; w < num_w; w++ ){
4281 
4282             FBE_singlet[ irrep_left ][ irrep_t ][ w ] = new double[ SIZE_left * SIZE_right ];
4283             FFG_singlet[ irrep_left ][ irrep_t ][ w ] = new double[ SIZE_left * SIZE_right ];
4284 
4285             double * BEptr = FBE_singlet[ irrep_left ][ irrep_t ][ w ];
4286             double * FGptr = FFG_singlet[ irrep_left ][ irrep_t ][ w ];
4287 
4288             if ( irrep_left == 0 ){ // irrep_t == irrep_w
4289                const int jump_singlet = jump_BF_active( indices, irrep_t, irrep_w, +1 );
4290                for ( int t = 0; t < w; t++ ){
4291                   for ( int xy = 0; xy < SIZE_left; xy++ ){
4292                      BEptr[ xy + SIZE_left * t ] = + SBB_singlet[ irrep_left ][ xy + SIZE_left * ( jump_singlet + t + ( w * ( w + 1 ) ) / 2 ) ];
4293                      FGptr[ xy + SIZE_left * t ] = - SFF_singlet[ irrep_left ][ xy + SIZE_left * ( jump_singlet + t + ( w * ( w + 1 ) ) / 2 ) ];
4294                   }
4295                }
4296                for ( int t = w; t < num_t; t++ ){
4297                   for ( int xy = 0; xy < SIZE_left; xy++ ){
4298                      BEptr[ xy + SIZE_left * t ] = + SBB_singlet[ irrep_left ][ xy + SIZE_left * ( jump_singlet + w + ( t * ( t + 1 ) ) / 2 ) ];
4299                      FGptr[ xy + SIZE_left * t ] = - SFF_singlet[ irrep_left ][ xy + SIZE_left * ( jump_singlet + w + ( t * ( t + 1 ) ) / 2 ) ];
4300                   }
4301                }
4302             } else { // irrep_t != irrep_w
4303                if ( irrep_t < irrep_w ){
4304                   const int jump_singlet = jump_BF_active( indices, irrep_t, irrep_w, +1 );
4305                   for ( int t = 0; t < num_t; t++ ){
4306                      for ( int xy = 0; xy < SIZE_left; xy++ ){
4307                         BEptr[ xy + SIZE_left * t ] = + SBB_singlet[ irrep_left ][ xy + SIZE_left * ( jump_singlet + t + num_t * w ) ];
4308                         FGptr[ xy + SIZE_left * t ] = - SFF_singlet[ irrep_left ][ xy + SIZE_left * ( jump_singlet + t + num_t * w ) ];
4309                      }
4310                   }
4311                } else {
4312                   const int jump_singlet = jump_BF_active( indices, irrep_w, irrep_t, +1 );
4313                   for ( int t = 0; t < num_t; t++ ){
4314                      for ( int xy = 0; xy < SIZE_left; xy++ ){
4315                         BEptr[ xy + SIZE_left * t ] = + SBB_singlet[ irrep_left ][ xy + SIZE_left * ( jump_singlet + w + num_w * t ) ];
4316                         FGptr[ xy + SIZE_left * t ] = - SFF_singlet[ irrep_left ][ xy + SIZE_left * ( jump_singlet + w + num_w * t ) ];
4317                      }
4318                   }
4319                }
4320             }
4321          }
4322       }
4323    }
4324 
4325 }
4326 
make_FBE_FFG_triplet()4327 void CheMPS2::CASPT2::make_FBE_FFG_triplet(){
4328 
4329    /*
4330       FBE triplet : < TB_xkyl E_wc TE_tiaj > = 2 delta_ac delta_ik delta_jl FBE_triplet[ Ik x Il ][ It ][ w ][ xy, t ]
4331       FFG triplet : < TF_cxdy E_kw TG_aibt > = 2 delta_ac delta_bd delta_ik FFG_triplet[ Ic x Id ][ It ][ w ][ xy, t ]
4332 
4333             FBE_triplet[ Ik x Il ][ It ][ w ][ xy, t ] = + SBB_triplet[ Ik x Il ][ xy, tw ]
4334             FFG_triplet[ Ic x Id ][ It ][ w ][ xy, t ] = + SFF_triplet[ Ic x Id ][ xy, tw ]
4335    */
4336 
4337    FBE_triplet = new double***[ num_irreps ];
4338    FFG_triplet = new double***[ num_irreps ];
4339 
4340    for ( int irrep_left = 0; irrep_left < num_irreps; irrep_left++ ){
4341 
4342       assert( size_B_triplet[ irrep_left ] == size_F_triplet[ irrep_left ] ); // At construction
4343       const int SIZE_left = size_B_triplet[ irrep_left ];
4344       FBE_triplet[ irrep_left ] = new double**[ num_irreps ];
4345       FFG_triplet[ irrep_left ] = new double**[ num_irreps ];
4346 
4347       for ( int irrep_t = 0; irrep_t < num_irreps; irrep_t++ ){
4348 
4349          const int SIZE_right = indices->getNDMRG( irrep_t );
4350          const int num_t   = indices->getNDMRG( irrep_t );
4351          const int irrep_w = Irreps::directProd( irrep_left, irrep_t );
4352          const int num_w   = indices->getNDMRG( irrep_w );
4353          FBE_triplet[ irrep_left ][ irrep_t ] = new double*[ num_w ];
4354          FFG_triplet[ irrep_left ][ irrep_t ] = new double*[ num_w ];
4355 
4356          for ( int w = 0; w < num_w; w++ ){
4357 
4358             FBE_triplet[ irrep_left ][ irrep_t ][ w ] = new double[ SIZE_left * SIZE_right ];
4359             FFG_triplet[ irrep_left ][ irrep_t ][ w ] = new double[ SIZE_left * SIZE_right ];
4360 
4361             double * BEptr = FBE_triplet[ irrep_left ][ irrep_t ][ w ];
4362             double * FGptr = FFG_triplet[ irrep_left ][ irrep_t ][ w ];
4363 
4364             if ( irrep_left == 0 ){ // irrep_t == irrep_w
4365                const int jump_triplet = jump_BF_active( indices, irrep_t, irrep_w, -1 );
4366                for ( int t = 0; t < w; t++ ){
4367                   for ( int xy = 0; xy < SIZE_left; xy++ ){
4368                      BEptr[ xy + SIZE_left * t ] = + SBB_triplet[ irrep_left ][ xy + SIZE_left * ( jump_triplet + t + ( w * ( w - 1 ) ) / 2 ) ];
4369                      FGptr[ xy + SIZE_left * t ] = + SFF_triplet[ irrep_left ][ xy + SIZE_left * ( jump_triplet + t + ( w * ( w - 1 ) ) / 2 ) ];
4370                   }
4371                }
4372                for ( int xy = 0; xy < SIZE_left; xy++ ){
4373                   BEptr[ xy + SIZE_left * w ] = 0.0;
4374                   FGptr[ xy + SIZE_left * w ] = 0.0;
4375                }
4376                for ( int t = w+1; t < num_t; t++ ){
4377                   for ( int xy = 0; xy < SIZE_left; xy++ ){
4378                      BEptr[ xy + SIZE_left * t ] = - SBB_triplet[ irrep_left ][ xy + SIZE_left * ( jump_triplet + w + ( t * ( t - 1 ) ) / 2 ) ];
4379                      FGptr[ xy + SIZE_left * t ] = - SFF_triplet[ irrep_left ][ xy + SIZE_left * ( jump_triplet + w + ( t * ( t - 1 ) ) / 2 ) ];
4380                   }
4381                }
4382             } else { // irrep_t != irrep_w
4383                if ( irrep_t < irrep_w ){
4384                   const int jump_triplet = jump_BF_active( indices, irrep_t, irrep_w, -1 );
4385                   for ( int t = 0; t < num_t; t++ ){
4386                      for ( int xy = 0; xy < SIZE_left; xy++ ){
4387                         BEptr[ xy + SIZE_left * t ] = + SBB_triplet[ irrep_left ][ xy + SIZE_left * ( jump_triplet + t + num_t * w ) ];
4388                         FGptr[ xy + SIZE_left * t ] = + SFF_triplet[ irrep_left ][ xy + SIZE_left * ( jump_triplet + t + num_t * w ) ];
4389                      }
4390                   }
4391                } else {
4392                   const int jump_triplet = jump_BF_active( indices, irrep_w, irrep_t, -1 );
4393                   for ( int t = 0; t < num_t; t++ ){
4394                      for ( int xy = 0; xy < SIZE_left; xy++ ){
4395                         BEptr[ xy + SIZE_left * t ] = - SBB_triplet[ irrep_left ][ xy + SIZE_left * ( jump_triplet + w + num_w * t ) ];
4396                         FGptr[ xy + SIZE_left * t ] = - SFF_triplet[ irrep_left ][ xy + SIZE_left * ( jump_triplet + w + num_w * t ) ];
4397                      }
4398                   }
4399                }
4400             }
4401          }
4402       }
4403    }
4404 
4405 }
4406 
make_FAB_FCF_singlet()4407 void CheMPS2::CASPT2::make_FAB_FCF_singlet(){
4408 
4409    /*
4410       FAB singlet : < E_zy E_lx | E_kw | SB_tiuj > = ( delta_ik delta_jl + delta_jk delta_il ) / sqrt( 1 + delta_ij ) * FAB_singlet[ Il ][ Ii x Ij ][ w ][ xyz, tu ]
4411 
4412       FCF singlet : < E_zy E_xd | E_wc | SF_atbu > = ( delta_ac delta_bd + delta_ad delta_bc ) / sqrt( 1 + delta_ab ) * FCF_singlet[ Id ][ Ia x Ib ][ w ][ xyz, tu ]
4413 
4414             FAB_singlet[ Il ][ Ii x Ij ][ w ][ xyz, tu ] = ( + 2 delta_tw delta_ux Gamma_zy
4415                                                              + 2 delta_uw delta_tx Gamma_zy
4416                                                              - delta_tw Gamma_zuyx
4417                                                              - delta_tw delta_uy Gamma_zx
4418                                                              - delta_uw Gamma_ztyx
4419                                                              - delta_uw delta_ty Gamma_zx
4420                                                              - SAA[ Il ][ xyz, utw ]
4421                                                              - SAA[ Il ][ xyz, tuw ]
4422                                                            )
4423 
4424             FCF_singlet[ Id ][ Ia x Ib ][ w ][ xyz, tu ] = ( + SCC[ Id ][ xyz, uwt ]
4425                                                              + SCC[ Id ][ xyz, twu ]
4426                                                              - delta_uw Gamma_zxyt
4427                                                              - delta_uw delta_xy Gamma_zt
4428                                                              - delta_tw Gamma_zxyu
4429                                                              - delta_tw delta_xy Gamma_zu
4430                                                            )
4431    */
4432 
4433    FAB_singlet = new double***[ num_irreps ];
4434    FCF_singlet = new double***[ num_irreps ];
4435 
4436    const int LAS = indices->getDMRGcumulative( num_irreps );
4437 
4438    for ( int irrep_left = 0; irrep_left < num_irreps; irrep_left++ ){
4439 
4440       assert( size_A[ irrep_left ] == size_C[ irrep_left ] ); // At construction
4441       const int SIZE_left = size_A[ irrep_left ];
4442       FAB_singlet[ irrep_left ] = new double**[ num_irreps ];
4443       FCF_singlet[ irrep_left ] = new double**[ num_irreps ];
4444 
4445       { // irrep_right == 0
4446 
4447          assert( size_B_singlet[ 0 ] == size_F_singlet[ 0 ] ); // At construction
4448          const int SIZE_right = size_B_singlet[ 0 ];
4449          const int num_w = indices->getNDMRG( irrep_left ); // irrep_w == irrep_left x irrep_right = irrep_left
4450          FAB_singlet[ irrep_left ][ 0 ] = new double*[ num_w ];
4451          FCF_singlet[ irrep_left ][ 0 ] = new double*[ num_w ];
4452 
4453          for ( int w = 0; w < num_w; w++ ){
4454 
4455             FAB_singlet[ irrep_left ][ 0 ][ w ] = new double[ SIZE_left * SIZE_right ];
4456             FCF_singlet[ irrep_left ][ 0 ][ w ] = new double[ SIZE_left * SIZE_right ];
4457 
4458             double * ABptr = FAB_singlet[ irrep_left ][ 0 ][ w ];
4459             double * CFptr = FCF_singlet[ irrep_left ][ 0 ][ w ];
4460 
4461             int jump_col = 0;
4462             for ( int irrep_ut = 0; irrep_ut < num_irreps; irrep_ut++ ){
4463                const int d_ut   = indices->getDMRGcumulative( irrep_ut );
4464                const int num_ut = indices->getNDMRG( irrep_ut );
4465                assert( jump_col == jump_BF_active( indices, irrep_ut, irrep_ut, +1 ) );
4466 
4467                const int jump_AB = jump_AC_active( indices, irrep_ut, irrep_ut, irrep_left );
4468                const int jump_CF = jump_AC_active( indices, irrep_ut, irrep_left, irrep_ut );
4469 
4470                for ( int t = 0; t < num_ut; t++ ){
4471                   for ( int u = t; u < num_ut; u++ ){ // 0 <= t <= u < num_ut
4472                      for ( int xyz = 0; xyz < SIZE_left; xyz++ ){
4473                         ABptr[ xyz + SIZE_left * ( jump_col + t + ( u * ( u + 1 ) ) / 2 ) ] = ( - SAA[ irrep_left ][ xyz + SIZE_left * ( jump_AB + u + num_ut * ( t + num_ut * w )) ]
4474                                                                                                 - SAA[ irrep_left ][ xyz + SIZE_left * ( jump_AB + t + num_ut * ( u + num_ut * w )) ] );
4475                         CFptr[ xyz + SIZE_left * ( jump_col + t + ( u * ( u + 1 ) ) / 2 ) ] = ( + SCC[ irrep_left ][ xyz + SIZE_left * ( jump_CF + u + num_ut * ( w + num_w  * t )) ]
4476                                                                                                 + SCC[ irrep_left ][ xyz + SIZE_left * ( jump_CF + t + num_ut * ( w + num_w  * u )) ] );
4477                      }
4478                   }
4479                }
4480 
4481                int jump_row = 0;
4482                for ( int irrep_x = 0; irrep_x < num_irreps; irrep_x++ ){
4483                   const int d_x   = indices->getDMRGcumulative( irrep_x );
4484                   const int num_x = indices->getNDMRG( irrep_x );
4485                   for ( int irrep_y = 0; irrep_y < num_irreps; irrep_y++ ){
4486                      const int irrep_z = Irreps::directProd( Irreps::directProd( irrep_left, irrep_x ), irrep_y );
4487                      const int d_y   = indices->getDMRGcumulative( irrep_y );
4488                      const int num_y = indices->getNDMRG( irrep_y );
4489                      const int d_z   = indices->getDMRGcumulative( irrep_z );
4490                      const int num_z = indices->getNDMRG( irrep_z );
4491                      assert( jump_row == jump_AC_active( indices, irrep_x, irrep_y, irrep_z ) );
4492 
4493                      if ( irrep_ut == irrep_left ){
4494 
4495                         // FAB_singlet[ xyz,tu ] -= delta_tw Gamma_zuyx
4496                         for ( int u = w; u < num_ut; u++ ){
4497                            for ( int x = 0; x < num_x; x++ ){
4498                               for ( int y = 0; y < num_y; y++ ){
4499                                  for ( int z = 0; z < num_z; z++ ){
4500                                     const double gamma_zuyx = two_rdm[ d_z + z + LAS * ( d_ut + u + LAS * ( d_y + y + LAS * ( d_x + x ))) ];
4501                                     ABptr[ jump_row + x + num_x * ( y + num_y * z ) + SIZE_left * ( jump_col + w + ( u * ( u + 1 ) ) / 2 ) ] -= gamma_zuyx;
4502                                  }
4503                               }
4504                            }
4505                         }
4506 
4507                         // FAB_singlet[ xyz,tu ] -= delta_tw delta_uy Gamma_zx
4508                         if ( irrep_ut == irrep_y ){
4509                            for ( int x = 0; x < num_x; x++ ){
4510                               for ( int z = 0; z < num_z; z++ ){
4511                                  const double gamma_zx = one_rdm[ d_z + z + LAS * ( d_x + x ) ];
4512                                  for ( int uy = w; uy < num_y; uy++ ){
4513                                     ABptr[ jump_row + x + num_x * ( uy + num_y * z ) + SIZE_left * ( jump_col + w + ( uy * ( uy + 1 ) ) / 2 ) ] -= gamma_zx;
4514                                  }
4515                               }
4516                            }
4517                         }
4518 
4519                         // FAB_singlet[ xyz,tu ] += 2 delta_tw delta_ux Gamma_zy
4520                         if ( irrep_ut == irrep_x ){
4521                            for ( int y = 0; y < num_y; y++ ){
4522                               for ( int z = 0; z < num_z; z++ ){
4523                                  const double gamma_zy = one_rdm[ d_z + z + LAS * ( d_y + y ) ];
4524                                  for ( int ux = w; ux < num_x; ux++ ){
4525                                     ABptr[ jump_row + ux + num_x * ( y + num_y * z ) + SIZE_left * ( jump_col + w + ( ux * ( ux + 1 ) ) / 2 ) ] += 2 * gamma_zy;
4526                                  }
4527                               }
4528                            }
4529                         }
4530 
4531                         // FCF_singlet[ xyz,tu ] -= delta_tw Gamma_zxyu
4532                         for ( int u = w; u < num_ut; u++ ){
4533                            for ( int x = 0; x < num_x; x++ ){
4534                               for ( int y = 0; y < num_y; y++ ){
4535                                  for ( int z = 0; z < num_z; z++ ){
4536                                     const double gamma_zxyu = two_rdm[ d_z + z + LAS * ( d_x + x + LAS * ( d_y + y + LAS * ( d_ut + u ))) ];
4537                                     CFptr[ jump_row + x + num_x * ( y + num_y * z ) + SIZE_left * ( jump_col + w + ( u * ( u + 1 ) ) / 2 ) ] -= gamma_zxyu;
4538                                  }
4539                               }
4540                            }
4541                         }
4542 
4543                         // FCF_singlet[ xyz,tu ] -= delta_tw delta_xy Gamma_zu
4544                         if ( irrep_x == irrep_y ){
4545                            for ( int u = w; u < num_ut; u++ ){
4546                               for ( int z = 0; z < num_z; z++ ){
4547                                  const double gamma_zu = one_rdm[ d_z + z + LAS * ( d_ut + u ) ];
4548                                  for ( int xy = 0; xy < num_x; xy++ ){
4549                                     CFptr[ jump_row + xy + num_x * ( xy + num_x * z ) + SIZE_left * ( jump_col + w + ( u * ( u + 1 ) ) / 2 ) ] -= gamma_zu;
4550                                  }
4551                               }
4552                            }
4553                         }
4554 
4555                         // FAB_singlet[ xyz,tu ] -= delta_uw Gamma_ztyx
4556                         for ( int t = 0; t <= w; t++ ){
4557                            for ( int x = 0; x < num_x; x++ ){
4558                               for ( int y = 0; y < num_y; y++ ){
4559                                  for ( int z = 0; z < num_z; z++ ){
4560                                     const double gamma_ztyx = two_rdm[ d_z + z + LAS * ( d_ut + t + LAS * ( d_y + y + LAS * ( d_x + x ))) ];
4561                                     ABptr[ jump_row + x + num_x * ( y + num_y * z ) + SIZE_left * ( jump_col + t + ( w * ( w + 1 ) ) / 2 ) ] -= gamma_ztyx;
4562                                  }
4563                               }
4564                            }
4565                         }
4566 
4567                         // FAB_singlet[ xyz,tu ] -= delta_uw delta_ty Gamma_zx
4568                         if ( irrep_ut == irrep_y ){
4569                            for ( int x = 0; x < num_x; x++ ){
4570                               for ( int z = 0; z < num_z; z++ ){
4571                                  const double gamma_zx = one_rdm[ d_z + z + LAS * ( d_x + x ) ];
4572                                  for ( int ty = 0; ty <= w; ty++ ){
4573                                     ABptr[ jump_row + x + num_x * ( ty + num_y * z ) + SIZE_left * ( jump_col + ty + ( w * ( w + 1 ) ) / 2 ) ] -= gamma_zx;
4574                                  }
4575                               }
4576                            }
4577                         }
4578 
4579                         // FAB_singlet[ xyz,tu ] += 2 delta_uw delta_tx Gamma_zy
4580                         if ( irrep_ut == irrep_x ){
4581                            for ( int y = 0; y < num_y; y++ ){
4582                               for ( int z = 0; z < num_z; z++ ){
4583                                  const double gamma_zy = one_rdm[ d_z + z + LAS * ( d_y + y ) ];
4584                                  for ( int tx = 0; tx <= w; tx++ ){
4585                                     ABptr[ jump_row + tx + num_x * ( y + num_y * z ) + SIZE_left * ( jump_col + tx + ( w * ( w + 1 ) ) / 2 ) ] += 2 * gamma_zy;
4586                                  }
4587                               }
4588                            }
4589                         }
4590 
4591                         // FCF_singlet[ xyz,tu ] -= delta_uw Gamma_zxyt
4592                         for ( int t = 0; t <= w; t++ ){
4593                            for ( int x = 0; x < num_x; x++ ){
4594                               for ( int y = 0; y < num_y; y++ ){
4595                                  for ( int z = 0; z < num_z; z++ ){
4596                                     const double gamma_zxyt = two_rdm[ d_z + z + LAS * ( d_x + x + LAS * ( d_y + y + LAS * ( d_ut + t ))) ];
4597                                     CFptr[ jump_row + x + num_x * ( y + num_y * z ) + SIZE_left * ( jump_col + t + ( w * ( w + 1 ) ) / 2 ) ] -= gamma_zxyt;
4598                                  }
4599                               }
4600                            }
4601                         }
4602 
4603                         // FCF_singlet[ xyz,tu ] -= delta_uw delta_xy Gamma_zt
4604                         if ( irrep_x == irrep_y ){
4605                            for ( int t = 0; t <= w; t++ ){
4606                               for ( int z = 0; z < num_z; z++ ){
4607                                  const double gamma_zt = one_rdm[ d_z + z + LAS * ( d_ut + t ) ];
4608                                  for ( int xy = 0; xy < num_x; xy++ ){
4609                                     CFptr[ jump_row + xy + num_x * ( xy + num_x * z ) + SIZE_left * ( jump_col + t + ( w * ( w + 1 ) ) / 2 ) ] -= gamma_zt;
4610                                  }
4611                               }
4612                            }
4613                         }
4614                      }
4615                      jump_row += num_x * num_y * num_z;
4616                   }
4617                }
4618                jump_col += ( num_ut * ( num_ut + 1 ) ) / 2;
4619             }
4620          }
4621       }
4622 
4623       for ( int irrep_right = 1; irrep_right < num_irreps; irrep_right++ ){
4624 
4625          assert( size_B_singlet[ irrep_right ] == size_F_singlet[ irrep_right ] ); // At construction
4626          const int SIZE_right = size_B_singlet[ irrep_right ];
4627          const int irrep_w = Irreps::directProd( irrep_left, irrep_right );
4628          const int num_w   = indices->getNDMRG( irrep_w );
4629          FAB_singlet[ irrep_left ][ irrep_right ] = new double*[ num_w ];
4630          FCF_singlet[ irrep_left ][ irrep_right ] = new double*[ num_w ];
4631 
4632          for ( int w = 0; w < num_w; w++ ){
4633 
4634             FAB_singlet[ irrep_left ][ irrep_right ][ w ] = new double[ SIZE_left * SIZE_right ];
4635             FCF_singlet[ irrep_left ][ irrep_right ][ w ] = new double[ SIZE_left * SIZE_right ];
4636 
4637             double * ABptr = FAB_singlet[ irrep_left ][ irrep_right ][ w ];
4638             double * CFptr = FCF_singlet[ irrep_left ][ irrep_right ][ w ];
4639 
4640             int jump_col = 0;
4641             for ( int irrep_t = 0; irrep_t < num_irreps; irrep_t++ ){
4642                const int irrep_u = Irreps::directProd( irrep_right, irrep_t );
4643                if ( irrep_t < irrep_u ){
4644                   const int d_t   = indices->getDMRGcumulative( irrep_t );
4645                   const int num_t = indices->getNDMRG( irrep_t );
4646                   const int d_u   = indices->getDMRGcumulative( irrep_u );
4647                   const int num_u = indices->getNDMRG( irrep_u );
4648                   assert( jump_col == jump_BF_active( indices, irrep_t, irrep_u, +1 ) );
4649 
4650                   const int jump_AB1 = jump_AC_active( indices, irrep_u, irrep_t, irrep_w );
4651                   const int jump_AB2 = jump_AC_active( indices, irrep_t, irrep_u, irrep_w );
4652                   const int jump_CF1 = jump_AC_active( indices, irrep_u, irrep_w, irrep_t );
4653                   const int jump_CF2 = jump_AC_active( indices, irrep_t, irrep_w, irrep_u );
4654 
4655                   for ( int t = 0; t < num_t; t++ ){
4656                      for ( int u = 0; u < num_u; u++ ){
4657                         for ( int xyz = 0; xyz < SIZE_left; xyz++ ){
4658                            ABptr[ xyz + SIZE_left * ( jump_col + t + num_t * u ) ] = ( - SAA[ irrep_left ][ xyz + SIZE_left * ( jump_AB1 + u + num_u * ( t + num_t * w )) ]
4659                                                                                        - SAA[ irrep_left ][ xyz + SIZE_left * ( jump_AB2 + t + num_t * ( u + num_u * w )) ] );
4660                            CFptr[ xyz + SIZE_left * ( jump_col + t + num_t * u ) ] = ( + SCC[ irrep_left ][ xyz + SIZE_left * ( jump_CF1 + u + num_u * ( w + num_w * t )) ]
4661                                                                                        + SCC[ irrep_left ][ xyz + SIZE_left * ( jump_CF2 + t + num_t * ( w + num_w * u )) ] );
4662                         }
4663                      }
4664                   }
4665 
4666                   int jump_row = 0;
4667                   for ( int irrep_x = 0; irrep_x < num_irreps; irrep_x++ ){
4668                      const int d_x   = indices->getDMRGcumulative( irrep_x );
4669                      const int num_x = indices->getNDMRG( irrep_x );
4670                      for ( int irrep_y = 0; irrep_y < num_irreps; irrep_y++ ){
4671                         const int irrep_z = Irreps::directProd( Irreps::directProd( irrep_left, irrep_x ), irrep_y );
4672                         const int d_y   = indices->getDMRGcumulative( irrep_y );
4673                         const int num_y = indices->getNDMRG( irrep_y );
4674                         const int d_z   = indices->getDMRGcumulative( irrep_z );
4675                         const int num_z = indices->getNDMRG( irrep_z );
4676                         assert( jump_row == jump_AC_active( indices, irrep_x, irrep_y, irrep_z ) );
4677 
4678                         if ( irrep_t == irrep_w ){
4679 
4680                            // FAB_singlet[ xyz,tu ] -= delta_tw Gamma_zuyx
4681                            for ( int u = 0; u < num_u; u++ ){
4682                               for ( int x = 0; x < num_x; x++ ){
4683                                  for ( int y = 0; y < num_y; y++ ){
4684                                     for ( int z = 0; z < num_z; z++ ){
4685                                        const double gamma_zuyx = two_rdm[ d_z + z + LAS * ( d_u + u + LAS * ( d_y + y + LAS * ( d_x + x ))) ];
4686                                        ABptr[ jump_row + x + num_x * ( y + num_y * z ) + SIZE_left * ( jump_col + w + num_w * u ) ] -= gamma_zuyx;
4687                                     }
4688                                  }
4689                               }
4690                            }
4691 
4692                            // FAB_singlet[ xyz,tu ] -= delta_tw delta_uy Gamma_zx
4693                            if ( irrep_u == irrep_y ){
4694                               for ( int x = 0; x < num_x; x++ ){
4695                                  for ( int z = 0; z < num_z; z++ ){
4696                                     const double gamma_zx = one_rdm[ d_z + z + LAS * ( d_x + x ) ];
4697                                     for ( int uy = 0; uy < num_y; uy++ ){
4698                                        ABptr[ jump_row + x + num_x * ( uy + num_y * z ) + SIZE_left * ( jump_col + w + num_w * uy ) ] -= gamma_zx;
4699                                     }
4700                                  }
4701                               }
4702                            }
4703 
4704                            // FAB_singlet[ xyz,tu ] += 2 delta_tw delta_ux Gamma_zy
4705                            if ( irrep_u == irrep_x ){
4706                               for ( int y = 0; y < num_y; y++ ){
4707                                  for ( int z = 0; z < num_z; z++ ){
4708                                     const double gamma_zy = one_rdm[ d_z + z + LAS * ( d_y + y ) ];
4709                                     for ( int ux = 0; ux < num_x; ux++ ){
4710                                        ABptr[ jump_row + ux + num_x * ( y + num_y * z ) + SIZE_left * ( jump_col + w + num_w * ux ) ] += 2 * gamma_zy;
4711                                     }
4712                                  }
4713                               }
4714                            }
4715 
4716                            // FCF_singlet[ xyz,tu ] -= delta_tw Gamma_zxyu
4717                            for ( int u = 0; u < num_u; u++ ){
4718                               for ( int x = 0; x < num_x; x++ ){
4719                                  for ( int y = 0; y < num_y; y++ ){
4720                                     for ( int z = 0; z < num_z; z++ ){
4721                                        const double gamma_zxyu = two_rdm[ d_z + z + LAS * ( d_x + x + LAS * ( d_y + y + LAS * ( d_u + u ))) ];
4722                                        CFptr[ jump_row + x + num_x * ( y + num_y * z ) + SIZE_left * ( jump_col + w + num_w * u ) ] -= gamma_zxyu;
4723                                     }
4724                                  }
4725                               }
4726                            }
4727 
4728                            // FCF_singlet[ xyz,tu ] -= delta_tw delta_xy Gamma_zu
4729                            if ( irrep_x == irrep_y ){
4730                               for ( int u = 0; u < num_u; u++ ){
4731                                  for ( int z = 0; z < num_z; z++ ){
4732                                     const double gamma_zu = one_rdm[ d_z + z + LAS * ( d_u + u ) ];
4733                                     for ( int xy = 0; xy < num_x; xy++ ){
4734                                        CFptr[ jump_row + xy + num_x * ( xy + num_x * z ) + SIZE_left * ( jump_col + w + num_w * u ) ] -= gamma_zu;
4735                                     }
4736                                  }
4737                               }
4738                            }
4739                         }
4740 
4741                         if ( irrep_u == irrep_w ){
4742 
4743                            // FAB_singlet[ xyz,tu ] -= delta_uw Gamma_ztyx
4744                            for ( int t = 0; t < num_t; t++ ){
4745                               for ( int x = 0; x < num_x; x++ ){
4746                                  for ( int y = 0; y < num_y; y++ ){
4747                                     for ( int z = 0; z < num_z; z++ ){
4748                                        const double gamma_ztyx = two_rdm[ d_z + z + LAS * ( d_t + t + LAS * ( d_y + y + LAS * ( d_x + x ))) ];
4749                                        ABptr[ jump_row + x + num_x * ( y + num_y * z ) + SIZE_left * ( jump_col + t + num_t * w ) ] -= gamma_ztyx;
4750                                     }
4751                                  }
4752                               }
4753                            }
4754 
4755                            // FAB_singlet[ xyz,tu ] -= delta_uw delta_ty Gamma_zx
4756                            if ( irrep_t == irrep_y ){
4757                               for ( int x = 0; x < num_x; x++ ){
4758                                  for ( int z = 0; z < num_z; z++ ){
4759                                     const double gamma_zx = one_rdm[ d_z + z + LAS * ( d_x + x ) ];
4760                                     for ( int ty = 0; ty < num_y; ty++ ){
4761                                        ABptr[ jump_row + x + num_x * ( ty + num_y * z ) + SIZE_left * ( jump_col + ty + num_y * w ) ] -= gamma_zx;
4762                                     }
4763                                  }
4764                               }
4765                            }
4766 
4767                            // FAB_singlet[ xyz,tu ] += 2 delta_uw delta_tx Gamma_zy
4768                            if ( irrep_t == irrep_x ){
4769                               for ( int y = 0; y < num_y; y++ ){
4770                                  for ( int z = 0; z < num_z; z++ ){
4771                                     const double gamma_zy = one_rdm[ d_z + z + LAS * ( d_y + y ) ];
4772                                     for ( int tx = 0; tx < num_x; tx++ ){
4773                                        ABptr[ jump_row + tx + num_x * ( y + num_y * z ) + SIZE_left * ( jump_col + tx + num_x * w ) ] += 2 * gamma_zy;
4774                                     }
4775                                  }
4776                               }
4777                            }
4778 
4779                            // FCF_singlet[ xyz,tu ] -= delta_uw Gamma_zxyt
4780                            for ( int t = 0; t < num_t; t++ ){
4781                               for ( int x = 0; x < num_x; x++ ){
4782                                  for ( int y = 0; y < num_y; y++ ){
4783                                     for ( int z = 0; z < num_z; z++ ){
4784                                        const double gamma_zxyt = two_rdm[ d_z + z + LAS * ( d_x + x + LAS * ( d_y + y + LAS * ( d_t + t ))) ];
4785                                        CFptr[ jump_row + x + num_x * ( y + num_y * z ) + SIZE_left * ( jump_col + t + num_t * w ) ] -= gamma_zxyt;
4786                                     }
4787                                  }
4788                               }
4789                            }
4790 
4791                            // FCF_singlet[ xyz,tu ] -= delta_uw delta_xy Gamma_zt
4792                            if ( irrep_x == irrep_y ){
4793                               for ( int t = 0; t < num_t; t++ ){
4794                                  for ( int z = 0; z < num_z; z++ ){
4795                                     const double gamma_zt = one_rdm[ d_z + z + LAS * ( d_t + t ) ];
4796                                     for ( int xy = 0; xy < num_x; xy++ ){
4797                                        CFptr[ jump_row + xy + num_x * ( xy + num_x * z ) + SIZE_left * ( jump_col + t + num_t * w ) ] -= gamma_zt;
4798                                     }
4799                                  }
4800                               }
4801                            }
4802                         }
4803                         jump_row += num_x * num_y * num_z;
4804                      }
4805                   }
4806                   jump_col += num_t * num_u;
4807                }
4808             }
4809          }
4810       }
4811    }
4812 
4813 }
4814 
make_FAB_FCF_triplet()4815 void CheMPS2::CASPT2::make_FAB_FCF_triplet(){
4816 
4817    /*
4818       FAB triplet : < E_zy E_lx | E_kw | TB_tiuj > = ( delta_ik delta_jl - delta_jk delta_il ) / sqrt( 1 + delta_ij ) * FAB_triplet[ Il ][ Ii x Ij ][ w ][ xyz, tu ]
4819 
4820       FCF triplet : < E_zy E_xd | E_wc | TF_atbu > = ( delta_ac delta_bd - delta_ad delta_bc ) / sqrt( 1 + delta_ab ) * FCF_triplet[ Id ][ Ia x Ib ][ w ][ xyz, tu ]
4821 
4822             FAB_triplet[ Il ][ Ii x Ij ][ w ][ xyz, tu ] = ( + 6 delta_tw delta_ux Gamma_zy
4823                                                              - 6 delta_uw delta_tx Gamma_zy
4824                                                              - 3 delta_tw Gamma_zuyx
4825                                                              - 3 delta_tw delta_uy Gamma_zx
4826                                                              + 3 delta_uw Gamma_ztyx
4827                                                              + 3 delta_uw delta_ty Gamma_zx
4828                                                              - SAA[ Il ][ xyz, utw ]
4829                                                              + SAA[ Il ][ xyz, tuw ]
4830                                                            )
4831 
4832             FCF_triplet[ Id ][ Ia x Ib ][ w ][ xyz, tu ] = ( + SCC[ Id ][ xyz, uwt ]
4833                                                              - SCC[ Id ][ xyz, twu ]
4834                                                              - delta_uw Gamma_zxyt
4835                                                              - delta_uw delta_xy Gamma_zt
4836                                                              + delta_tw Gamma_zxyu
4837                                                              + delta_tw delta_xy Gamma_zu
4838                                                            )
4839    */
4840 
4841    FAB_triplet = new double***[ num_irreps ];
4842    FCF_triplet = new double***[ num_irreps ];
4843 
4844    const int LAS = indices->getDMRGcumulative( num_irreps );
4845 
4846    for ( int irrep_left = 0; irrep_left < num_irreps; irrep_left++ ){
4847 
4848       assert( size_A[ irrep_left ] == size_C[ irrep_left ] ); // At construction
4849       const int SIZE_left = size_A[ irrep_left ];
4850       FAB_triplet[ irrep_left ] = new double**[ num_irreps ];
4851       FCF_triplet[ irrep_left ] = new double**[ num_irreps ];
4852 
4853       { // irrep_right == 0
4854 
4855          assert( size_B_triplet[ 0 ] == size_F_triplet[ 0 ] ); // At construction
4856          const int SIZE_right = size_B_triplet[ 0 ];
4857          const int num_w = indices->getNDMRG( irrep_left );
4858          FAB_triplet[ irrep_left ][ 0 ] = new double*[ num_w ];
4859          FCF_triplet[ irrep_left ][ 0 ] = new double*[ num_w ];
4860 
4861          for ( int w = 0; w < num_w; w++ ){
4862 
4863             FAB_triplet[ irrep_left ][ 0 ][ w ] = new double[ SIZE_left * SIZE_right ];
4864             FCF_triplet[ irrep_left ][ 0 ][ w ] = new double[ SIZE_left * SIZE_right ];
4865 
4866             double * ABptr = FAB_triplet[ irrep_left ][ 0 ][ w ];
4867             double * CFptr = FCF_triplet[ irrep_left ][ 0 ][ w ];
4868 
4869             int jump_col = 0;
4870             for ( int irrep_ut = 0; irrep_ut < num_irreps; irrep_ut++ ){
4871                const int d_ut   = indices->getDMRGcumulative( irrep_ut );
4872                const int num_ut = indices->getNDMRG( irrep_ut );
4873                assert( jump_col == jump_BF_active( indices, irrep_ut, irrep_ut, -1 ) );
4874 
4875                const int jump_AB = jump_AC_active( indices, irrep_ut, irrep_ut, irrep_left );
4876                const int jump_CF = jump_AC_active( indices, irrep_ut, irrep_left, irrep_ut );
4877 
4878                for ( int t = 0; t < num_ut; t++ ){
4879                   for ( int u = t+1; u < num_ut; u++ ){ // 0 <= t < u < num_ut
4880                      for ( int xyz = 0; xyz < SIZE_left; xyz++ ){
4881                         ABptr[ xyz + SIZE_left * ( jump_col + t + ( u * ( u - 1 ) ) / 2 ) ] = ( - SAA[ irrep_left ][ xyz + SIZE_left * ( jump_AB + u + num_ut * ( t + num_ut * w )) ]
4882                                                                                                 + SAA[ irrep_left ][ xyz + SIZE_left * ( jump_AB + t + num_ut * ( u + num_ut * w )) ] );
4883                         CFptr[ xyz + SIZE_left * ( jump_col + t + ( u * ( u - 1 ) ) / 2 ) ] = ( + SCC[ irrep_left ][ xyz + SIZE_left * ( jump_CF + u + num_ut * ( w + num_w  * t )) ]
4884                                                                                                 - SCC[ irrep_left ][ xyz + SIZE_left * ( jump_CF + t + num_ut * ( w + num_w  * u )) ] );
4885                      }
4886                   }
4887                }
4888 
4889                int jump_row = 0;
4890                for ( int irrep_x = 0; irrep_x < num_irreps; irrep_x++ ){
4891                   const int d_x   = indices->getDMRGcumulative( irrep_x );
4892                   const int num_x = indices->getNDMRG( irrep_x );
4893                   for ( int irrep_y = 0; irrep_y < num_irreps; irrep_y++ ){
4894                      const int irrep_z = Irreps::directProd( Irreps::directProd( irrep_left, irrep_x ), irrep_y );
4895                      const int d_y   = indices->getDMRGcumulative( irrep_y );
4896                      const int num_y = indices->getNDMRG( irrep_y );
4897                      const int d_z   = indices->getDMRGcumulative( irrep_z );
4898                      const int num_z = indices->getNDMRG( irrep_z );
4899                      assert( jump_row == jump_AC_active( indices, irrep_x, irrep_y, irrep_z ) );
4900 
4901                      if ( irrep_ut == irrep_left ){
4902 
4903                         // FAB_triplet[ xyz,tu ] -= 3 delta_tw Gamma_zuyx
4904                         for ( int u = w+1; u < num_ut; u++ ){
4905                            for ( int x = 0; x < num_x; x++ ){
4906                               for ( int y = 0; y < num_y; y++ ){
4907                                  for ( int z = 0; z < num_z; z++ ){
4908                                     const double gamma_zuyx = two_rdm[ d_z + z + LAS * ( d_ut + u + LAS * ( d_y + y + LAS * ( d_x + x ))) ];
4909                                     ABptr[ jump_row + x + num_x * ( y + num_y * z ) + SIZE_left * ( jump_col + w + ( u * ( u - 1 ) ) / 2 ) ] -= 3 * gamma_zuyx;
4910                                  }
4911                               }
4912                            }
4913                         }
4914 
4915                         // FAB_triplet[ xyz,tu ] -= 3 delta_tw delta_uy Gamma_zx
4916                         if ( irrep_ut == irrep_y ){
4917                            for ( int x = 0; x < num_x; x++ ){
4918                               for ( int z = 0; z < num_z; z++ ){
4919                                  const double gamma_zx = one_rdm[ d_z + z + LAS * ( d_x + x ) ];
4920                                  for ( int uy = w+1; uy < num_y; uy++ ){
4921                                     ABptr[ jump_row + x + num_x * ( uy + num_y * z ) + SIZE_left * ( jump_col + w + ( uy * ( uy - 1 ) ) / 2 ) ] -= 3 * gamma_zx;
4922                                  }
4923                               }
4924                            }
4925                         }
4926 
4927                         // FAB_triplet[ xyz,tu ] += 6 delta_tw delta_ux Gamma_zy
4928                         if ( irrep_ut == irrep_x ){
4929                            for ( int y = 0; y < num_y; y++ ){
4930                               for ( int z = 0; z < num_z; z++ ){
4931                                  const double gamma_zy = one_rdm[ d_z + z + LAS * ( d_y + y ) ];
4932                                  for ( int ux = w+1; ux < num_x; ux++ ){
4933                                     ABptr[ jump_row + ux + num_x * ( y + num_y * z ) + SIZE_left * ( jump_col + w + ( ux * ( ux - 1 ) ) / 2 ) ] += 6 * gamma_zy;
4934                                  }
4935                               }
4936                            }
4937                         }
4938 
4939                         // FCF_triplet[ xyz,tu ] += delta_tw Gamma_zxyu
4940                         for ( int u = w+1; u < num_ut; u++ ){
4941                            for ( int x = 0; x < num_x; x++ ){
4942                               for ( int y = 0; y < num_y; y++ ){
4943                                  for ( int z = 0; z < num_z; z++ ){
4944                                     const double gamma_zxyu = two_rdm[ d_z + z + LAS * ( d_x + x + LAS * ( d_y + y + LAS * ( d_ut + u ))) ];
4945                                     CFptr[ jump_row + x + num_x * ( y + num_y * z ) + SIZE_left * ( jump_col + w + ( u * ( u - 1 ) ) / 2 ) ] += gamma_zxyu;
4946                                  }
4947                               }
4948                            }
4949                         }
4950 
4951                         // FCF_triplet[ xyz,tu ] += delta_tw delta_xy Gamma_zu
4952                         if ( irrep_x == irrep_y ){
4953                            for ( int u = w+1; u < num_ut; u++ ){
4954                               for ( int z = 0; z < num_z; z++ ){
4955                                  const double gamma_zu = one_rdm[ d_z + z + LAS * ( d_ut + u ) ];
4956                                  for ( int xy = 0; xy < num_x; xy++ ){
4957                                     CFptr[ jump_row + xy + num_x * ( xy + num_x * z ) + SIZE_left * ( jump_col + w + ( u * ( u - 1 ) ) / 2 ) ] += gamma_zu;
4958                                  }
4959                               }
4960                            }
4961                         }
4962 
4963                         // FAB_triplet[ xyz,tu ] += 3 delta_uw Gamma_ztyx
4964                         for ( int t = 0; t < w; t++ ){
4965                            for ( int x = 0; x < num_x; x++ ){
4966                               for ( int y = 0; y < num_y; y++ ){
4967                                  for ( int z = 0; z < num_z; z++ ){
4968                                     const double gamma_ztyx = two_rdm[ d_z + z + LAS * ( d_ut + t + LAS * ( d_y + y + LAS * ( d_x + x ))) ];
4969                                     ABptr[ jump_row + x + num_x * ( y + num_y * z ) + SIZE_left * ( jump_col + t + ( w * ( w - 1 ) ) / 2 ) ] += 3 * gamma_ztyx;
4970                                  }
4971                               }
4972                            }
4973                         }
4974 
4975                         // FAB_triplet[ xyz,tu ] += 3 delta_uw delta_ty Gamma_zx
4976                         if ( irrep_ut == irrep_y ){
4977                            for ( int x = 0; x < num_x; x++ ){
4978                               for ( int z = 0; z < num_z; z++ ){
4979                                  const double gamma_zx = one_rdm[ d_z + z + LAS * ( d_x + x ) ];
4980                                  for ( int ty = 0; ty < w; ty++ ){
4981                                     ABptr[ jump_row + x + num_x * ( ty + num_y * z ) + SIZE_left * ( jump_col + ty + ( w * ( w - 1 ) ) / 2 ) ] += 3 * gamma_zx;
4982                                  }
4983                               }
4984                            }
4985                         }
4986 
4987                         // FAB_triplet[ xyz,tu ] -= 6 delta_uw delta_tx Gamma_zy
4988                         if ( irrep_ut == irrep_x ){
4989                            for ( int y = 0; y < num_y; y++ ){
4990                               for ( int z = 0; z < num_z; z++ ){
4991                                  const double gamma_zy = one_rdm[ d_z + z + LAS * ( d_y + y ) ];
4992                                  for ( int tx = 0; tx < w; tx++ ){
4993                                     ABptr[ jump_row + tx + num_x * ( y + num_y * z ) + SIZE_left * ( jump_col + tx + ( w * ( w - 1 ) ) / 2 ) ] -= 6 * gamma_zy;
4994                                  }
4995                               }
4996                            }
4997                         }
4998 
4999                         // FCF_triplet[ xyz,tu ] -= delta_uw Gamma_zxyt
5000                         for ( int t = 0; t < w; t++ ){
5001                            for ( int x = 0; x < num_x; x++ ){
5002                               for ( int y = 0; y < num_y; y++ ){
5003                                  for ( int z = 0; z < num_z; z++ ){
5004                                     const double gamma_zxyt = two_rdm[ d_z + z + LAS * ( d_x + x + LAS * ( d_y + y + LAS * ( d_ut + t ))) ];
5005                                     CFptr[ jump_row + x + num_x * ( y + num_y * z ) + SIZE_left * ( jump_col + t + ( w * ( w - 1 ) ) / 2 ) ] -= gamma_zxyt;
5006                                  }
5007                               }
5008                            }
5009                         }
5010 
5011                         // FCF_triplet[ xyz,tu ] -= delta_uw delta_xy Gamma_zt
5012                         if ( irrep_x == irrep_y ){
5013                            for ( int t = 0; t < w; t++ ){
5014                               for ( int z = 0; z < num_z; z++ ){
5015                                  const double gamma_zt = one_rdm[ d_z + z + LAS * ( d_ut + t ) ];
5016                                  for ( int xy = 0; xy < num_x; xy++ ){
5017                                     CFptr[ jump_row + xy + num_x * ( xy + num_x * z ) + SIZE_left * ( jump_col + t + ( w * ( w - 1 ) ) / 2 ) ] -= gamma_zt;
5018                                  }
5019                               }
5020                            }
5021                         }
5022                      }
5023                      jump_row += num_x * num_y * num_z;
5024                   }
5025                }
5026                jump_col += ( num_ut * ( num_ut - 1 ) ) / 2;
5027             }
5028          }
5029       }
5030 
5031       for ( int irrep_right = 1; irrep_right < num_irreps; irrep_right++ ){
5032 
5033          assert( size_B_triplet[ irrep_right ] == size_F_triplet[ irrep_right ] ); // At construction
5034          const int SIZE_right = size_B_triplet[ irrep_right ];
5035          const int irrep_w = Irreps::directProd( irrep_left, irrep_right );
5036          const int num_w   = indices->getNDMRG( irrep_w );
5037          FAB_triplet[ irrep_left ][ irrep_right ] = new double*[ num_w ];
5038          FCF_triplet[ irrep_left ][ irrep_right ] = new double*[ num_w ];
5039 
5040          for ( int w = 0; w < num_w; w++ ){
5041 
5042             FAB_triplet[ irrep_left ][ irrep_right ][ w ] = new double[ SIZE_left * SIZE_right ];
5043             FCF_triplet[ irrep_left ][ irrep_right ][ w ] = new double[ SIZE_left * SIZE_right ];
5044 
5045             double * ABptr = FAB_triplet[ irrep_left ][ irrep_right ][ w ];
5046             double * CFptr = FCF_triplet[ irrep_left ][ irrep_right ][ w ];
5047 
5048             int jump_col = 0;
5049             for ( int irrep_t = 0; irrep_t < num_irreps; irrep_t++ ){
5050                const int irrep_u = Irreps::directProd( irrep_right, irrep_t );
5051                if ( irrep_t < irrep_u ){
5052                   const int d_t   = indices->getDMRGcumulative( irrep_t );
5053                   const int num_t = indices->getNDMRG( irrep_t );
5054                   const int d_u   = indices->getDMRGcumulative( irrep_u );
5055                   const int num_u = indices->getNDMRG( irrep_u );
5056                   assert( jump_col == jump_BF_active( indices, irrep_t, irrep_u, -1 ) );
5057 
5058                   const int jump_AB1 = jump_AC_active( indices, irrep_u, irrep_t, irrep_w );
5059                   const int jump_AB2 = jump_AC_active( indices, irrep_t, irrep_u, irrep_w );
5060                   const int jump_CF1 = jump_AC_active( indices, irrep_u, irrep_w, irrep_t );
5061                   const int jump_CF2 = jump_AC_active( indices, irrep_t, irrep_w, irrep_u );
5062 
5063                   for ( int t = 0; t < num_t; t++ ){
5064                      for ( int u = 0; u < num_u; u++ ){
5065                         for ( int xyz = 0; xyz < SIZE_left; xyz++ ){
5066                            ABptr[ xyz + SIZE_left * ( jump_col + t + num_t * u ) ] = ( - SAA[ irrep_left ][ xyz + SIZE_left * ( jump_AB1 + u + num_u * ( t + num_t * w )) ]
5067                                                                                        + SAA[ irrep_left ][ xyz + SIZE_left * ( jump_AB2 + t + num_t * ( u + num_u * w )) ] );
5068                            CFptr[ xyz + SIZE_left * ( jump_col + t + num_t * u ) ] = ( + SCC[ irrep_left ][ xyz + SIZE_left * ( jump_CF1 + u + num_u * ( w + num_w * t )) ]
5069                                                                                        - SCC[ irrep_left ][ xyz + SIZE_left * ( jump_CF2 + t + num_t * ( w + num_w * u )) ] );
5070                         }
5071                      }
5072                   }
5073 
5074                   int jump_row = 0;
5075                   for ( int irrep_x = 0; irrep_x < num_irreps; irrep_x++ ){
5076                      const int d_x   = indices->getDMRGcumulative( irrep_x );
5077                      const int num_x = indices->getNDMRG( irrep_x );
5078                      for ( int irrep_y = 0; irrep_y < num_irreps; irrep_y++ ){
5079                         const int irrep_z = Irreps::directProd( Irreps::directProd( irrep_left, irrep_x ), irrep_y );
5080                         const int d_y   = indices->getDMRGcumulative( irrep_y );
5081                         const int num_y = indices->getNDMRG( irrep_y );
5082                         const int d_z   = indices->getDMRGcumulative( irrep_z );
5083                         const int num_z = indices->getNDMRG( irrep_z );
5084                         assert( jump_row == jump_AC_active( indices, irrep_x, irrep_y, irrep_z ) );
5085 
5086                         if ( irrep_t == irrep_w ){
5087 
5088                            // FAB_triplet[ xyz,tu ] -= 3 delta_tw Gamma_zuyx
5089                            for ( int u = 0; u < num_u; u++ ){
5090                               for ( int x = 0; x < num_x; x++ ){
5091                                  for ( int y = 0; y < num_y; y++ ){
5092                                     for ( int z = 0; z < num_z; z++ ){
5093                                        const double gamma_zuyx = two_rdm[ d_z + z + LAS * ( d_u + u + LAS * ( d_y + y + LAS * ( d_x + x ))) ];
5094                                        ABptr[ jump_row + x + num_x * ( y + num_y * z ) + SIZE_left * ( jump_col + w + num_w * u ) ] -= 3 * gamma_zuyx;
5095                                     }
5096                                  }
5097                               }
5098                            }
5099 
5100                            // FAB_triplet[ xyz,tu ] -= 3 delta_tw delta_uy Gamma_zx
5101                            if ( irrep_u == irrep_y ){
5102                               for ( int x = 0; x < num_x; x++ ){
5103                                  for ( int z = 0; z < num_z; z++ ){
5104                                     const double gamma_zx = one_rdm[ d_z + z + LAS * ( d_x + x ) ];
5105                                     for ( int uy = 0; uy < num_y; uy++ ){
5106                                        ABptr[ jump_row + x + num_x * ( uy + num_y * z ) + SIZE_left * ( jump_col + w + num_w * uy ) ] -= 3 * gamma_zx;
5107                                     }
5108                                  }
5109                               }
5110                            }
5111 
5112                            // FAB_triplet[ xyz,tu ] += 6 delta_tw delta_ux Gamma_zy
5113                            if ( irrep_u == irrep_x ){
5114                               for ( int y = 0; y < num_y; y++ ){
5115                                  for ( int z = 0; z < num_z; z++ ){
5116                                     const double gamma_zy = one_rdm[ d_z + z + LAS * ( d_y + y ) ];
5117                                     for ( int ux = 0; ux < num_x; ux++ ){
5118                                        ABptr[ jump_row + ux + num_x * ( y + num_y * z ) + SIZE_left * ( jump_col + w + num_w * ux ) ] += 6 * gamma_zy;
5119                                     }
5120                                  }
5121                               }
5122                            }
5123 
5124                            // FCF_triplet[ xyz,tu ] += delta_tw Gamma_zxyu
5125                            for ( int u = 0; u < num_u; u++ ){
5126                               for ( int x = 0; x < num_x; x++ ){
5127                                  for ( int y = 0; y < num_y; y++ ){
5128                                     for ( int z = 0; z < num_z; z++ ){
5129                                        const double gamma_zxyu = two_rdm[ d_z + z + LAS * ( d_x + x + LAS * ( d_y + y + LAS * ( d_u + u ))) ];
5130                                        CFptr[ jump_row + x + num_x * ( y + num_y * z ) + SIZE_left * ( jump_col + w + num_w * u ) ] += gamma_zxyu;
5131                                     }
5132                                  }
5133                               }
5134                            }
5135 
5136                            // FCF_triplet[ xyz,tu ] += delta_tw delta_xy Gamma_zu
5137                            if ( irrep_x == irrep_y ){
5138                               for ( int u = 0; u < num_u; u++ ){
5139                                  for ( int z = 0; z < num_z; z++ ){
5140                                     const double gamma_zu = one_rdm[ d_z + z + LAS * ( d_u + u ) ];
5141                                     for ( int xy = 0; xy < num_x; xy++ ){
5142                                        CFptr[ jump_row + xy + num_x * ( xy + num_x * z ) + SIZE_left * ( jump_col + w + num_w * u ) ] += gamma_zu;
5143                                     }
5144                                  }
5145                               }
5146                            }
5147                         }
5148 
5149                         if ( irrep_u == irrep_w ){
5150 
5151                            // FAB_triplet[ xyz,tu ] += 3 * delta_uw Gamma_ztyx
5152                            for ( int t = 0; t < num_t; t++ ){
5153                               for ( int x = 0; x < num_x; x++ ){
5154                                  for ( int y = 0; y < num_y; y++ ){
5155                                     for ( int z = 0; z < num_z; z++ ){
5156                                        const double gamma_ztyx = two_rdm[ d_z + z + LAS * ( d_t + t + LAS * ( d_y + y + LAS * ( d_x + x ))) ];
5157                                        ABptr[ jump_row + x + num_x * ( y + num_y * z ) + SIZE_left * ( jump_col + t + num_t * w ) ] += 3 * gamma_ztyx;
5158                                     }
5159                                  }
5160                               }
5161                            }
5162 
5163                            // FAB_triplet[ xyz,tu ] += 3 * delta_uw delta_ty Gamma_zx
5164                            if ( irrep_t == irrep_y ){
5165                               for ( int x = 0; x < num_x; x++ ){
5166                                  for ( int z = 0; z < num_z; z++ ){
5167                                     const double gamma_zx = one_rdm[ d_z + z + LAS * ( d_x + x ) ];
5168                                     for ( int ty = 0; ty < num_y; ty++ ){
5169                                        ABptr[ jump_row + x + num_x * ( ty + num_y * z ) + SIZE_left * ( jump_col + ty + num_y * w ) ] += 3 * gamma_zx;
5170                                     }
5171                                  }
5172                               }
5173                            }
5174 
5175                            // FAB_triplet[ xyz,tu ] -= 6 delta_uw delta_tx Gamma_zy
5176                            if ( irrep_t == irrep_x ){
5177                               for ( int y = 0; y < num_y; y++ ){
5178                                  for ( int z = 0; z < num_z; z++ ){
5179                                     const double gamma_zy = one_rdm[ d_z + z + LAS * ( d_y + y ) ];
5180                                     for ( int tx = 0; tx < num_x; tx++ ){
5181                                        ABptr[ jump_row + tx + num_x * ( y + num_y * z ) + SIZE_left * ( jump_col + tx + num_x * w ) ] -= 6 * gamma_zy;
5182                                     }
5183                                  }
5184                               }
5185                            }
5186 
5187                            // FCF_triplet[ xyz,tu ] -= delta_uw Gamma_zxyt
5188                            for ( int t = 0; t < num_t; t++ ){
5189                               for ( int x = 0; x < num_x; x++ ){
5190                                  for ( int y = 0; y < num_y; y++ ){
5191                                     for ( int z = 0; z < num_z; z++ ){
5192                                        const double gamma_zxyt = two_rdm[ d_z + z + LAS * ( d_x + x + LAS * ( d_y + y + LAS * ( d_t + t ))) ];
5193                                        CFptr[ jump_row + x + num_x * ( y + num_y * z ) + SIZE_left * ( jump_col + t + num_t * w ) ] -= gamma_zxyt;
5194                                     }
5195                                  }
5196                               }
5197                            }
5198 
5199                            // FCF_triplet[ xyz,tu ] -= delta_uw delta_xy Gamma_zt
5200                            if ( irrep_x == irrep_y ){
5201                               for ( int t = 0; t < num_t; t++ ){
5202                                  for ( int z = 0; z < num_z; z++ ){
5203                                     const double gamma_zt = one_rdm[ d_z + z + LAS * ( d_t + t ) ];
5204                                     for ( int xy = 0; xy < num_x; xy++ ){
5205                                        CFptr[ jump_row + xy + num_x * ( xy + num_x * z ) + SIZE_left * ( jump_col + t + num_t * w ) ] -= gamma_zt;
5206                                     }
5207                                  }
5208                               }
5209                            }
5210                         }
5211                         jump_row += num_x * num_y * num_z;
5212                      }
5213                   }
5214                   jump_col += num_t * num_u;
5215                }
5216             }
5217          }
5218       }
5219    }
5220 
5221 }
5222 
make_FAD_FCD()5223 void CheMPS2::CASPT2::make_FAD_FCD(){
5224 
5225    /*
5226       FAD1: < E_zy E_jx E_wc E_ai E_tu > = delta_ac delta_ij FAD1[ Ij ][ Ii x Ia ][ w ][ xyz, tu ]
5227       FAD2: < E_zy E_jx E_wc E_ti E_au > = delta_ac delta_ij FAD2[ Ij ][ Ii x Ia ][ w ][ xyz, tu ]
5228       FCD1: < E_zy E_xb E_kw E_ai E_tu > = delta_ik delta_ab FCD1[ Ib ][ Ii x Ia ][ w ][ xyz, tu ]
5229       FCD2: < E_zy E_xb E_kw E_ti E_au > = delta_ik delta_ab FCD2[ Ib ][ Ii x Ia ][ w ][ xyz, tu ]
5230 
5231             FAD1[ Ij ][ Ii x Ia ][ w ][ xyz, tu ] = + SAA[ Ij ][ xyzwtu ]
5232 
5233             FAD2[ Ij ][ Ii x Ia ][ w ][ xyz, tu ] = + SAA[ Ij ][ xyztwu ]
5234 
5235             FCD1[ Ib ][ Ii x Ia ][ w ][ xyz, tu ] = - SCC[ Ib ][ xyzwtu ]
5236 
5237             FCD2[ Ib ][ Ii x Ia ][ w ][ xyz, tu ] = ( + 2 delta_tw Gamma_zxyu
5238                                                       + 2 delta_tw delta_xy Gamma_zu
5239                                                       + delta_tu Gamma_zxyw
5240                                                       + delta_tu delta_xy Gamma_zw
5241                                                       - SCC[ Ib ][ xyzutw ]
5242                                                     )
5243    */
5244 
5245    FAD = new double***[ num_irreps ];
5246    FCD = new double***[ num_irreps ];
5247 
5248    const int LAS = indices->getDMRGcumulative( num_irreps );
5249 
5250    for ( int irrep_left = 0; irrep_left < num_irreps; irrep_left++ ){
5251 
5252       assert( size_A[ irrep_left ] == size_C[ irrep_left ] ); // At construction
5253       const int SIZE_left = size_A[ irrep_left ];
5254       FAD[ irrep_left ] = new double**[ num_irreps ];
5255       FCD[ irrep_left ] = new double**[ num_irreps ];
5256 
5257       for ( int irrep_right = 0; irrep_right < num_irreps; irrep_right++ ){
5258 
5259          const int SIZE_right = size_D[ irrep_right ];
5260          const int D2JUMP  = SIZE_right / 2;
5261          const int irrep_w = Irreps::directProd( irrep_left, irrep_right );
5262          const int d_w     = indices->getDMRGcumulative( irrep_w );
5263          const int num_w   = indices->getNDMRG( irrep_w );
5264          FAD[ irrep_left ][ irrep_right ] = new double*[ num_w ];
5265          FCD[ irrep_left ][ irrep_right ] = new double*[ num_w ];
5266 
5267          for ( int w = 0; w < num_w; w++ ){
5268 
5269             FAD[ irrep_left ][ irrep_right ][ w ] = new double[ SIZE_left * SIZE_right ];
5270             FCD[ irrep_left ][ irrep_right ][ w ] = new double[ SIZE_left * SIZE_right ];
5271 
5272             double * AD1ptr = FAD[ irrep_left ][ irrep_right ][ w ];
5273             double * AD2ptr = FAD[ irrep_left ][ irrep_right ][ w ] + SIZE_left * D2JUMP;
5274             double * CD1ptr = FCD[ irrep_left ][ irrep_right ][ w ];
5275             double * CD2ptr = FCD[ irrep_left ][ irrep_right ][ w ] + SIZE_left * D2JUMP;
5276 
5277             int jump_col = 0;
5278             for ( int irrep_t = 0; irrep_t < num_irreps; irrep_t++ ){
5279                const int irrep_u = Irreps::directProd( irrep_right, irrep_t );
5280                const int num_t = indices->getNDMRG( irrep_t );
5281                const int d_u   = indices->getDMRGcumulative( irrep_u );
5282                const int num_u = indices->getNDMRG( irrep_u );
5283 
5284                const int jump_AD1 = jump_AC_active( indices, irrep_w, irrep_t, irrep_u );
5285                const int jump_AD2 = jump_AC_active( indices, irrep_t, irrep_w, irrep_u );
5286                const int jump_CD2 = jump_AC_active( indices, irrep_u, irrep_t, irrep_w );
5287 
5288                for ( int t = 0; t < num_t; t++ ){
5289                   for ( int u = 0; u < num_u; u++ ){
5290                      for ( int xyz = 0; xyz < SIZE_left; xyz++ ){
5291                         AD1ptr[ xyz + SIZE_left * ( jump_col + t + num_t * u ) ] = + SAA[ irrep_left ][ xyz + SIZE_left * ( jump_AD1 + w + num_w * ( t + num_t * u )) ];
5292                         AD2ptr[ xyz + SIZE_left * ( jump_col + t + num_t * u ) ] = + SAA[ irrep_left ][ xyz + SIZE_left * ( jump_AD2 + t + num_t * ( w + num_w * u )) ];
5293                         CD1ptr[ xyz + SIZE_left * ( jump_col + t + num_t * u ) ] = - SCC[ irrep_left ][ xyz + SIZE_left * ( jump_AD1 + w + num_w * ( t + num_t * u )) ];
5294                         CD2ptr[ xyz + SIZE_left * ( jump_col + t + num_t * u ) ] = - SCC[ irrep_left ][ xyz + SIZE_left * ( jump_CD2 + u + num_u * ( t + num_t * w )) ];
5295                      }
5296                   }
5297                }
5298 
5299                int jump_row = 0;
5300                for ( int irrep_x = 0; irrep_x < num_irreps; irrep_x++ ){
5301                   const int d_x   = indices->getDMRGcumulative( irrep_x );
5302                   const int num_x = indices->getNDMRG( irrep_x );
5303                   for ( int irrep_y = 0; irrep_y < num_irreps; irrep_y++ ){
5304                      const int irrep_z = Irreps::directProd( Irreps::directProd( irrep_left, irrep_x ), irrep_y );
5305                      const int d_y   = indices->getDMRGcumulative( irrep_y );
5306                      const int num_y = indices->getNDMRG( irrep_y );
5307                      const int d_z   = indices->getDMRGcumulative( irrep_z );
5308                      const int num_z = indices->getNDMRG( irrep_z );
5309                      assert( jump_row == jump_AC_active( indices, irrep_x, irrep_y, irrep_z ) );
5310 
5311                      if ( irrep_t == irrep_w ){
5312                         // + 2 delta_tw Gamma_zxyu
5313                         for ( int u = 0; u < num_u; u++ ){
5314                            for ( int x = 0; x < num_x; x++ ){
5315                               for ( int y = 0; y < num_y; y++ ){
5316                                  for ( int z = 0; z < num_z; z++ ){
5317                                     const double gamma_zxyu = two_rdm[ d_z + z + LAS * ( d_x + x + LAS * ( d_y + y + LAS * ( d_u + u ))) ];
5318                                     CD2ptr[ jump_row + x + num_x * ( y + num_y * z ) + SIZE_left * ( jump_col + w + num_t * u ) ] += 2 * gamma_zxyu;
5319                                  }
5320                               }
5321                            }
5322                         }
5323 
5324                         // + 2 delta_tw delta_xy Gamma_zu
5325                         if ( irrep_x == irrep_y ){
5326                            for ( int u = 0; u < num_u; u++ ){
5327                               for ( int z = 0; z < num_z; z++ ){
5328                                  const double gamma_zu = one_rdm[ d_z + z + LAS * ( d_u + u ) ];
5329                                  for ( int xy = 0; xy < num_x; xy++ ){
5330                                     CD2ptr[ jump_row + xy + num_x * ( xy + num_x * z ) + SIZE_left * ( jump_col + w + num_t * u ) ] += 2 * gamma_zu;
5331                                  }
5332                               }
5333                            }
5334                         }
5335                      }
5336 
5337                      if ( irrep_t == irrep_u ){
5338                         // + delta_tu Gamma_zxyw
5339                         for ( int tu = 0; tu < num_u; tu++ ){
5340                            for ( int x = 0; x < num_x; x++ ){
5341                               for ( int y = 0; y < num_y; y++ ){
5342                                  for ( int z = 0; z < num_z; z++ ){
5343                                     const double gamma_zxyw = two_rdm[ d_z + z + LAS * ( d_x + x + LAS * ( d_y + y + LAS * ( d_w + w ))) ];
5344                                     CD2ptr[ jump_row + x + num_x * ( y + num_y * z ) + SIZE_left * ( jump_col + tu + num_u * tu ) ] += gamma_zxyw;
5345                                  }
5346                               }
5347                            }
5348                         }
5349 
5350                         // + delta_tu delta_xy Gamma_zw
5351                         if ( irrep_x == irrep_y ){
5352                            for ( int z = 0; z < num_z; z++ ){
5353                               const double gamma_zw = one_rdm[ d_z + z + LAS * ( d_w + w ) ];
5354                               for ( int tu = 0; tu < num_u; tu++ ){
5355                                  for ( int xy = 0; xy < num_x; xy++ ){
5356                                     CD2ptr[ jump_row + xy + num_x * ( xy + num_x * z ) + SIZE_left * ( jump_col + tu + num_u * tu ) ] += gamma_zw;
5357                                  }
5358                               }
5359                            }
5360                         }
5361                      }
5362                      jump_row += num_x * num_y * num_z;
5363                   }
5364                }
5365                jump_col += num_t * num_u;
5366             }
5367          }
5368       }
5369    }
5370 
5371 }
5372 
make_AA_CC(const bool OVLP,const double IPEA)5373 void CheMPS2::CASPT2::make_AA_CC( const bool OVLP, const double IPEA ){
5374 
5375    /*
5376       SAA: < E_zy E_jx | 1 | E_ti E_uv > = delta_ij ( SAA[ Ii ][ xyztuv ] )
5377       FAA: < E_zy E_jx | F | E_ti E_uv > = delta_ij ( FAA[ Ii ][ xyztuv ] + ( 2 sum_k f_kk - f_ii ) SAA[ Ii ][ xyztuv ] )
5378 
5379          SAA[ Ii ][ xyztuv ] = ( - Gamma_{ztuyxv}
5380                                  + 2 delta_tx Gamma_{zuyv}
5381                                  -   delta_uy Gamma_{tzxv}
5382                                  -   delta_ty Gamma_{zuxv}
5383                                  -   delta_ux Gamma_{ztyv}
5384                                  + ( 2 delta_tx delta_uy - delta_ux delta_ty ) Gamma_{zv}
5385                                )
5386 
5387          FAA[ Ii ][ xyztuv ] = ( - f_dot_4dm[ ztuyxv ] + ( f_tt + f_uu + f_xx + f_yy ) SAA[ Ii ][ xyztuv ]
5388                                  + 2 delta_tx ( f_dot_3dm[ zuyv ] - f_tt Gamma_{zuyv} )
5389                                  -   delta_uy ( f_dot_3dm[ ztvx ] - f_uu Gamma_{ztvx} )
5390                                  -   delta_ty ( f_dot_3dm[ zuxv ] - f_tt Gamma_{zuxv} )
5391                                  -   delta_ux ( f_dot_3dm[ ztyv ] - f_uu Gamma_{ztyv} )
5392                                  + ( 2 delta_tx delta_uy - delta_ux delta_ty ) ( f_dot_2dm[ zv ] - ( f_tt + f_uu ) Gamma_{zv} )
5393                                )
5394 
5395       SCC: < E_zy E_xb | 1 | E_at E_uv > = delta_ab ( SCC[ Ia ][ xyztuv ] )
5396       FCC: < E_zy E_xb | F | E_at E_uv > = delta_ab ( FCC[ Ia ][ xyztuv ] + ( 2 sum_k f_kk + f_aa ) SCC[ Ia ][ xyztuv ] )
5397 
5398          SCC[ Ia ][ xyztuv ] = ( + Gamma_{zxuytv}
5399                                  + delta_uy Gamma_{xztv}
5400                                  + delta_xy Gamma_{zutv}
5401                                  + delta_ut Gamma_{zxyv}
5402                                  + delta_ut delta_xy Gamma_{zv}
5403                                )
5404 
5405          FCC[ Ia ][ xyztuv ] = ( + f_dot_4dm[ zxuytv ] + ( f_uu + f_yy ) SCC[ Ia ][ xyztuv ]
5406                                  + delta_uy ( f_dot_3dm[ xztv ] - f_uu Gamma_{xztv} )
5407                                  + delta_xy ( f_dot_3dm[ zutv ] - f_xx Gamma_{zutv} )
5408                                  + delta_ut ( f_dot_3dm[ zxyv ] - f_tt Gamma_{zxyv} )
5409                                  + delta_ut delta_xy ( f_dot_2dm[ zv ] - ( f_tt + f_xx ) Gamma_{zv} )
5410                                )
5411    */
5412 
5413    if ( OVLP ){ SAA = new double*[ num_irreps ];
5414                 SCC = new double*[ num_irreps ]; }
5415    else {       FAA = new double*[ num_irreps ];
5416                 FCC = new double*[ num_irreps ]; }
5417 
5418    const int LAS = indices->getDMRGcumulative( num_irreps );
5419 
5420    for ( int irrep = 0; irrep < num_irreps; irrep++ ){
5421 
5422       assert( size_A[ irrep ] == size_C[ irrep ] ); // At construction
5423       const int SIZE = size_A[ irrep ];
5424       if ( OVLP ){ SAA[ irrep ] = new double[ SIZE * SIZE ];
5425                    SCC[ irrep ] = new double[ SIZE * SIZE ]; }
5426       else {       FAA[ irrep ] = new double[ SIZE * SIZE ];
5427                    FCC[ irrep ] = new double[ SIZE * SIZE ]; }
5428 
5429       int jump_col = 0;
5430       for ( int irrep_t = 0; irrep_t < num_irreps; irrep_t++ ){
5431          const int d_t    = indices->getDMRGcumulative( irrep_t );
5432          const int num_t  = indices->getNDMRG( irrep_t );
5433          const int nocc_t = indices->getNOCC( irrep_t );
5434          for ( int irrep_u = 0; irrep_u < num_irreps; irrep_u++ ){
5435             const int d_u     = indices->getDMRGcumulative( irrep_u );
5436             const int num_u   = indices->getNDMRG( irrep_u );
5437             const int nocc_u  = indices->getNOCC( irrep_u );
5438             const int irrep_v = Irreps::directProd( Irreps::directProd( irrep, irrep_t ), irrep_u );
5439             const int d_v     = indices->getDMRGcumulative( irrep_v );
5440             const int num_v   = indices->getNDMRG( irrep_v );
5441             int jump_row = 0;
5442             for ( int irrep_x = 0; irrep_x < num_irreps; irrep_x++ ){
5443                const int d_x    = indices->getDMRGcumulative( irrep_x );
5444                const int num_x  = indices->getNDMRG( irrep_x );
5445                const int nocc_x = indices->getNOCC( irrep_x );
5446                for ( int irrep_y = 0; irrep_y < num_irreps; irrep_y++ ){
5447                   const int d_y     = indices->getDMRGcumulative( irrep_y );
5448                   const int num_y   = indices->getNDMRG( irrep_y );
5449                   const int nocc_y  = indices->getNOCC( irrep_y );
5450                   const int irrep_z = Irreps::directProd( Irreps::directProd( irrep, irrep_x ), irrep_y );
5451                   const int d_z     = indices->getDMRGcumulative( irrep_z );
5452                   const int num_z   = indices->getNDMRG( irrep_z );
5453 
5454                   // SAA: - Gamma_{ztuyxv}
5455                   // FAA: - f_dot_4dm[ ztuyxv ] + ( f_tt + f_uu + f_xx + f_yy ) SAA[ Ii ][ xyztuv ]
5456                   if ( OVLP ){
5457                      #pragma omp parallel for schedule(static)
5458                      for ( int v = 0; v < num_v; v++ ){
5459                         for ( int u = 0; u < num_u; u++ ){
5460                            for ( int t = 0; t < num_t; t++ ){
5461                               for ( int z = 0; z < num_z; z++ ){
5462                                  for ( int y = 0; y < num_y; y++ ){
5463                                     for ( int x = 0; x < num_x; x++ ){
5464                                        const int ptr = jump_row + x + num_x * ( y + num_y * z ) + SIZE * ( jump_col + t + num_t * ( u + num_u * v ) );
5465                                        SAA[ irrep ][ ptr ] = - three_rdm[ d_z + z + LAS * ( d_t + t + LAS * ( d_u + u + LAS * ( d_y + y + LAS * ( d_x + x + LAS * ( d_v + v ))))) ];
5466                                     }
5467                                  }
5468                               }
5469                            }
5470                         }
5471                      }
5472                   } else {
5473                      #pragma omp parallel for schedule(static)
5474                      for ( int v = 0; v < num_v; v++ ){
5475                         for ( int u = 0; u < num_u; u++ ){
5476                            const double f_uu = fock->get( irrep_u, nocc_u + u, nocc_u + u );
5477                            for ( int t = 0; t < num_t; t++ ){
5478                               const double f_tt = fock->get( irrep_t, nocc_t + t, nocc_t + t );
5479                               for ( int z = 0; z < num_z; z++ ){
5480                                  for ( int y = 0; y < num_y; y++ ){
5481                                     const double f_yy = fock->get( irrep_y, nocc_y + y, nocc_y + y );
5482                                     for ( int x = 0; x < num_x; x++ ){
5483                                        const double f_xx = fock->get( irrep_x, nocc_x + x, nocc_x + x );
5484                                        const int ptr = jump_row + x + num_x * ( y + num_y * z ) + SIZE * ( jump_col + t + num_t * ( u + num_u * v ) );
5485                                        FAA[ irrep ][ ptr ] = - f_dot_4dm[ d_z + z + LAS * ( d_t + t + LAS * ( d_u + u + LAS * ( d_y + y + LAS * ( d_x + x + LAS * ( d_v + v ))))) ]
5486                                                              + ( f_tt + f_uu + f_xx + f_yy ) * SAA[ irrep ][ ptr ];
5487                                     }
5488                                  }
5489                               }
5490                            }
5491                         }
5492                      }
5493                   }
5494 
5495                   // SCC: + Gamma_{zxuytv}
5496                   // FCC: + f_dot_4dm[ zxuytv ] + ( f_uu + f_yy ) SCC[ Ia ][ xyztuv ]
5497                   if ( OVLP ){
5498                      #pragma omp parallel for schedule(static)
5499                      for ( int v = 0; v < num_v; v++ ){
5500                         for ( int u = 0; u < num_u; u++ ){
5501                            for ( int t = 0; t < num_t; t++ ){
5502                               for ( int z = 0; z < num_z; z++ ){
5503                                  for ( int y = 0; y < num_y; y++ ){
5504                                     for ( int x = 0; x < num_x; x++ ){
5505                                        const int ptr = jump_row + x + num_x * ( y + num_y * z ) + SIZE * ( jump_col + t + num_t * ( u + num_u * v ) );
5506                                        SCC[ irrep ][ ptr ] = three_rdm[ d_z + z + LAS * ( d_x + x + LAS * ( d_u + u + LAS * ( d_y + y + LAS * ( d_t + t + LAS * ( d_v + v ))))) ];
5507                                     }
5508                                  }
5509                               }
5510                            }
5511                         }
5512                      }
5513                   } else {
5514                      #pragma omp parallel for schedule(static)
5515                      for ( int v = 0; v < num_v; v++ ){
5516                         for ( int u = 0; u < num_u; u++ ){
5517                            const double f_uu = fock->get( irrep_u, nocc_u + u, nocc_u + u );
5518                            for ( int t = 0; t < num_t; t++ ){
5519                               for ( int z = 0; z < num_z; z++ ){
5520                                  for ( int y = 0; y < num_y; y++ ){
5521                                     const double f_yy = fock->get( irrep_y, nocc_y + y, nocc_y + y );
5522                                     for ( int x = 0; x < num_x; x++ ){
5523                                        const int ptr = jump_row + x + num_x * ( y + num_y * z ) + SIZE * ( jump_col + t + num_t * ( u + num_u * v ) );
5524                                        FCC[ irrep ][ ptr ] = f_dot_4dm[ d_z + z + LAS * ( d_x + x + LAS * ( d_u + u + LAS * ( d_y + y + LAS * ( d_t + t + LAS * ( d_v + v ))))) ]
5525                                                            + ( f_yy + f_uu ) * SCC[ irrep ][ ptr ];
5526                                     }
5527                                  }
5528                               }
5529                            }
5530                         }
5531                      }
5532                   }
5533 
5534                   // SAA: + 2 delta_tx Gamma_{zuyv}
5535                   // FAA: + 2 delta_tx ( f_dot_3dm[ zuyv ] - f_tt Gamma_{zuyv} )
5536                   if ( irrep_t == irrep_x ){
5537                      if ( OVLP ){
5538                         #pragma omp parallel for schedule(static)
5539                         for ( int v = 0; v < num_v; v++ ){
5540                            for ( int u = 0; u < num_u; u++ ){
5541                               for ( int xt = 0; xt < num_t; xt++ ){
5542                                  for ( int z = 0; z < num_z; z++ ){
5543                                     for ( int y = 0; y < num_y; y++ ){
5544                                        SAA[ irrep ][ jump_row + xt + num_x * ( y + num_y * z ) + SIZE * ( jump_col + xt + num_t * ( u + num_u * v ) ) ]
5545                                           += 2 * two_rdm[ d_z + z + LAS * ( d_u + u + LAS * ( d_y + y + LAS * ( d_v + v ))) ];
5546                                     }
5547                                  }
5548                               }
5549                            }
5550                         }
5551                      } else {
5552                         #pragma omp parallel for schedule(static)
5553                         for ( int v = 0; v < num_v; v++ ){
5554                            for ( int u = 0; u < num_u; u++ ){
5555                               for ( int tx = 0; tx < num_t; tx++ ){
5556                                  const double f_tt = fock->get( irrep_t, nocc_t + tx, nocc_t + tx );
5557                                  for ( int z = 0; z < num_z; z++ ){
5558                                     for ( int y = 0; y < num_y; y++ ){
5559                                        FAA[ irrep ][ jump_row + tx + num_x * ( y + num_y * z ) + SIZE * ( jump_col + tx + num_t * ( u + num_u * v ) ) ]
5560                                           += 2 * ( f_dot_3dm[ d_z + z + LAS * ( d_u + u + LAS * ( d_y + y + LAS * ( d_v + v ))) ]
5561                                             - f_tt * two_rdm[ d_z + z + LAS * ( d_u + u + LAS * ( d_y + y + LAS * ( d_v + v ))) ] );
5562                                     }
5563                                  }
5564                               }
5565                            }
5566                         }
5567                      }
5568                   }
5569 
5570                   // SAA: - delta_uy Gamma_{ztvx}
5571                   // FAA: - delta_uy ( f_dot_3dm[ ztvx ] - f_uu Gamma_{ztvx} )
5572                   if ( irrep_u == irrep_y ){
5573                      if ( OVLP ){
5574                         #pragma omp parallel for schedule(static)
5575                         for ( int v = 0; v < num_v; v++ ){
5576                            for ( int uy = 0; uy < num_u; uy++ ){
5577                               for ( int t = 0; t < num_t; t++ ){
5578                                  for ( int z = 0; z < num_z; z++ ){
5579                                     for ( int x = 0; x < num_x; x++ ){
5580                                        SAA[ irrep ][ jump_row + x + num_x * ( uy + num_y * z ) + SIZE * ( jump_col + t + num_t * ( uy + num_u * v ) ) ]
5581                                           -= two_rdm[ d_t + t + LAS * ( d_z + z + LAS * ( d_x + x + LAS * ( d_v + v ))) ];
5582                                     }
5583                                  }
5584                               }
5585                            }
5586                         }
5587                      } else {
5588                         #pragma omp parallel for schedule(static)
5589                         for ( int v = 0; v < num_v; v++ ){
5590                            for ( int uy = 0; uy < num_u; uy++ ){
5591                               const double f_uu = fock->get( irrep_u, nocc_u + uy, nocc_u + uy );
5592                               for ( int t = 0; t < num_t; t++ ){
5593                                  for ( int z = 0; z < num_z; z++ ){
5594                                     for ( int x = 0; x < num_x; x++ ){
5595                                        FAA[ irrep ][ jump_row + x + num_x * ( uy + num_y * z ) + SIZE * ( jump_col + t + num_t * ( uy + num_u * v ) ) ]
5596                                           -= ( f_dot_3dm[ d_t + t + LAS * ( d_z + z + LAS * ( d_x + x + LAS * ( d_v + v ))) ]
5597                                         - f_uu * two_rdm[ d_t + t + LAS * ( d_z + z + LAS * ( d_x + x + LAS * ( d_v + v ))) ] );
5598                                     }
5599                                  }
5600                               }
5601                            }
5602                         }
5603                      }
5604                   }
5605 
5606                   // SAA: - delta_ty Gamma_{zuxv}
5607                   // FAA: - delta_ty ( f_dot_3dm[ zuxv ] - f_tt Gamma_{zuxv} )
5608                   if ( irrep_t == irrep_y ){
5609                      if ( OVLP ){
5610                         #pragma omp parallel for schedule(static)
5611                         for ( int v = 0; v < num_v; v++ ){
5612                            for ( int u = 0; u < num_u; u++ ){
5613                               for ( int ty = 0; ty < num_t; ty++ ){
5614                                  for ( int z = 0; z < num_z; z++ ){
5615                                     for ( int x = 0; x < num_x; x++ ){
5616                                        SAA[ irrep ][ jump_row + x + num_x * ( ty + num_y * z ) + SIZE * ( jump_col + ty + num_t * ( u + num_u * v ) ) ]
5617                                           -= two_rdm[ d_z + z + LAS * ( d_u + u + LAS * ( d_x + x + LAS * ( d_v + v ))) ];
5618                                     }
5619                                  }
5620                               }
5621                            }
5622                         }
5623                      } else {
5624                         #pragma omp parallel for schedule(static)
5625                         for ( int v = 0; v < num_v; v++ ){
5626                            for ( int u = 0; u < num_u; u++ ){
5627                               for ( int ty = 0; ty < num_t; ty++ ){
5628                                  const double f_tt = fock->get( irrep_t, nocc_t + ty, nocc_t + ty );
5629                                  for ( int z = 0; z < num_z; z++ ){
5630                                     for ( int x = 0; x < num_x; x++ ){
5631                                        FAA[ irrep ][ jump_row + x + num_x * ( ty + num_y * z ) + SIZE * ( jump_col + ty + num_t * ( u + num_u * v ) ) ]
5632                                           -= ( f_dot_3dm[ d_z + z + LAS * ( d_u + u + LAS * ( d_x + x + LAS * ( d_v + v ))) ]
5633                                         - f_tt * two_rdm[ d_z + z + LAS * ( d_u + u + LAS * ( d_x + x + LAS * ( d_v + v ))) ] );
5634                                     }
5635                                  }
5636                               }
5637                            }
5638                         }
5639                      }
5640                   }
5641 
5642                   // SAA: - delta_ux Gamma_{ztyv}
5643                   // FAA: - delta_ux ( f_dot_3dm[ ztyv ] - f_uu Gamma_{ztyv} )
5644                   if ( irrep_u == irrep_x ){
5645                      if ( OVLP ){
5646                         #pragma omp parallel for schedule(static)
5647                         for ( int v = 0; v < num_v; v++ ){
5648                            for ( int ux = 0; ux < num_u; ux++ ){
5649                               for ( int t = 0; t < num_t; t++ ){
5650                                  for ( int z = 0; z < num_z; z++ ){
5651                                     for ( int y = 0; y < num_y; y++ ){
5652                                        SAA[ irrep ][ jump_row + ux + num_x * ( y + num_y * z ) + SIZE * ( jump_col + t + num_t * ( ux + num_u * v ) ) ]
5653                                           -= two_rdm[ d_z + z + LAS * ( d_t + t + LAS * ( d_y + y + LAS * ( d_v + v ))) ];
5654                                     }
5655                                  }
5656                               }
5657                            }
5658                         }
5659                      } else {
5660                         #pragma omp parallel for schedule(static)
5661                         for ( int v = 0; v < num_v; v++ ){
5662                            for ( int ux = 0; ux < num_u; ux++ ){
5663                               const double f_uu = fock->get( irrep_u, nocc_u + ux, nocc_u + ux );
5664                               for ( int t = 0; t < num_t; t++ ){
5665                                  for ( int z = 0; z < num_z; z++ ){
5666                                     for ( int y = 0; y < num_y; y++ ){
5667                                        FAA[ irrep ][ jump_row + ux + num_x * ( y + num_y * z ) + SIZE * ( jump_col + t + num_t * ( ux + num_u * v ) ) ]
5668                                           -= ( f_dot_3dm[ d_z + z + LAS * ( d_t + t + LAS * ( d_y + y + LAS * ( d_v + v ))) ]
5669                                         - f_uu * two_rdm[ d_z + z + LAS * ( d_t + t + LAS * ( d_y + y + LAS * ( d_v + v ))) ] );
5670                                     }
5671                                  }
5672                               }
5673                            }
5674                         }
5675                      }
5676                   }
5677 
5678                   // SCC: + delta_uy Gamma_{xztv}
5679                   // FCC: + delta_uy ( f_dot_3dm[ xztv ] - f_uu Gamma_{xztv} )
5680                   if ( irrep_u == irrep_y ){
5681                      if ( OVLP ){
5682                         #pragma omp parallel for schedule(static)
5683                         for ( int v = 0; v < num_v; v++ ){
5684                            for ( int uy = 0; uy < num_u; uy++ ){
5685                               for ( int t = 0; t < num_t; t++ ){
5686                                  for ( int z = 0; z < num_z; z++ ){
5687                                     for ( int x = 0; x < num_x; x++ ){
5688                                        SCC[ irrep ][ jump_row + x + num_x * ( uy + num_y * z ) + SIZE * ( jump_col + t + num_t * ( uy + num_u * v ) ) ]
5689                                           += two_rdm[ d_x + x + LAS * ( d_z + z + LAS * ( d_t + t + LAS * ( d_v + v ))) ];
5690                                     }
5691                                  }
5692                               }
5693                            }
5694                         }
5695                      } else {
5696                         #pragma omp parallel for schedule(static)
5697                         for ( int v = 0; v < num_v; v++ ){
5698                            for ( int uy = 0; uy < num_y; uy++ ){
5699                               const double f_uu = fock->get( irrep_y, nocc_y + uy, nocc_y + uy );
5700                               for ( int t = 0; t < num_t; t++ ){
5701                                  for ( int z = 0; z < num_z; z++ ){
5702                                     for ( int x = 0; x < num_x; x++ ){
5703                                        FCC[ irrep ][ jump_row + x + num_x * ( uy + num_y * z ) + SIZE * ( jump_col + t + num_t * ( uy + num_y * v ) ) ]
5704                                           += ( f_dot_3dm[ d_x + x + LAS * ( d_z + z + LAS * ( d_t + t + LAS * ( d_v + v ))) ]
5705                                         - f_uu * two_rdm[ d_x + x + LAS * ( d_z + z + LAS * ( d_t + t + LAS * ( d_v + v ))) ] );
5706                                     }
5707                                  }
5708                               }
5709                            }
5710                         }
5711                      }
5712                   }
5713 
5714                   // SCC: + delta_xy Gamma_{zutv}
5715                   // FCC: + delta_xy ( f_dot_3dm[ zutv ] - f_xx Gamma_{zutv} )
5716                   if ( irrep_x == irrep_y ){
5717                      if ( OVLP ){
5718                         #pragma omp parallel for schedule(static)
5719                         for ( int v = 0; v < num_v; v++ ){
5720                            for ( int u = 0; u < num_u; u++ ){
5721                               for ( int t = 0; t < num_t; t++ ){
5722                                  for ( int z = 0; z < num_z; z++ ){
5723                                     for ( int xy = 0; xy < num_x; xy++ ){
5724                                        SCC[ irrep ][ jump_row + xy + num_x * ( xy + num_y * z ) + SIZE * ( jump_col + t + num_t * ( u + num_u * v ) ) ]
5725                                           += two_rdm[ d_z + z + LAS * ( d_u + u + LAS * ( d_t + t + LAS * ( d_v + v ))) ];
5726                                     }
5727                                  }
5728                               }
5729                            }
5730                         }
5731                      } else {
5732                         #pragma omp parallel for schedule(static)
5733                         for ( int v = 0; v < num_v; v++ ){
5734                            for ( int u = 0; u < num_u; u++ ){
5735                               for ( int t = 0; t < num_t; t++ ){
5736                                  for ( int z = 0; z < num_z; z++ ){
5737                                     for ( int xy = 0; xy < num_x; xy++ ){
5738                                        const double f_xx = fock->get( irrep_x, nocc_x + xy, nocc_x + xy );
5739                                        FCC[ irrep ][ jump_row + xy + num_x * ( xy + num_x * z ) + SIZE * ( jump_col + t + num_t * ( u + num_u * v ) ) ]
5740                                           += ( f_dot_3dm[ d_z + z + LAS * ( d_u + u + LAS * ( d_t + t + LAS * ( d_v + v ))) ]
5741                                         - f_xx * two_rdm[ d_z + z + LAS * ( d_u + u + LAS * ( d_t + t + LAS * ( d_v + v ))) ] );
5742                                     }
5743                                  }
5744                               }
5745                            }
5746                         }
5747                      }
5748                   }
5749 
5750                   // SCC: + delta_ut Gamma_{zxyv}
5751                   // FCC: + delta_ut ( f_dot_3dm[ zxyv ] - f_tt Gamma_{zxyv} )
5752                   if ( irrep_u == irrep_t ){
5753                      if ( OVLP ){
5754                         #pragma omp parallel for schedule(static)
5755                         for ( int v = 0; v < num_v; v++ ){
5756                            for ( int ut = 0; ut < num_u; ut++ ){
5757                               for ( int z = 0; z < num_z; z++ ){
5758                                  for ( int y = 0; y < num_y; y++ ){
5759                                     for ( int x = 0; x < num_x; x++ ){
5760                                        SCC[ irrep ][ jump_row + x + num_x * ( y + num_y * z ) + SIZE * ( jump_col + ut + num_t * ( ut + num_u * v ) ) ]
5761                                           += two_rdm[ d_z + z + LAS * ( d_x + x + LAS * ( d_y + y + LAS * ( d_v + v ))) ];
5762                                     }
5763                                  }
5764                               }
5765                            }
5766                         }
5767                      } else {
5768                         #pragma omp parallel for schedule(static)
5769                         for ( int v = 0; v < num_v; v++ ){
5770                            for ( int ut = 0; ut < num_u; ut++ ){
5771                               const double f_tt = fock->get( irrep_t, nocc_t + ut, nocc_t + ut );
5772                               for ( int z = 0; z < num_z; z++ ){
5773                                  for ( int y = 0; y < num_y; y++ ){
5774                                     for ( int x = 0; x < num_x; x++ ){
5775                                        FCC[ irrep ][ jump_row + x + num_x * ( y + num_y * z ) + SIZE * ( jump_col + ut + num_u * ( ut + num_u * v ) ) ]
5776                                           += ( f_dot_3dm[ d_z + z + LAS * ( d_x + x + LAS * ( d_y + y + LAS * ( d_v + v ))) ]
5777                                         - f_tt * two_rdm[ d_z + z + LAS * ( d_x + x + LAS * ( d_y + y + LAS * ( d_v + v ))) ] );
5778                                     }
5779                                  }
5780                               }
5781                            }
5782                         }
5783                      }
5784                   }
5785 
5786                   // SAA: + 2 delta_tx delta_uy Gamma_{zv}
5787                   // FAA: + 2 delta_tx delta_uy ( f_dot_2dm[ zv ] - ( f_tt + f_uu ) Gamma_{zv} )
5788                   if (( irrep_t == irrep_x ) && ( irrep_u == irrep_y ) && ( irrep_z == irrep_v )){
5789                      if ( OVLP ){
5790                         #pragma omp parallel for schedule(static)
5791                         for ( int v = 0; v < num_v; v++ ){
5792                            for ( int z = 0; z < num_z; z++ ){
5793                               for ( int uy = 0; uy < num_u; uy++ ){
5794                                  for ( int xt = 0; xt < num_t; xt++ ){
5795                                     SAA[ irrep ][ jump_row + xt + num_x * ( uy + num_y * z ) + SIZE * ( jump_col + xt + num_t * ( uy + num_u * v ) ) ]
5796                                        += 2 * one_rdm[ d_z + z + LAS * ( d_v + v ) ];
5797                                  }
5798                               }
5799                            }
5800                         }
5801                      } else {
5802                         #pragma omp parallel for schedule(static)
5803                         for ( int v = 0; v < num_v; v++ ){
5804                            for ( int z = 0; z < num_z; z++ ){
5805                               for ( int uy = 0; uy < num_u; uy++ ){
5806                                  const double f_uu = fock->get( irrep_u, nocc_u + uy, nocc_u + uy );
5807                                  for ( int xt = 0; xt < num_t; xt++ ){
5808                                     const double f_tt = fock->get( irrep_t, nocc_t + xt, nocc_t + xt );
5809                                     FAA[ irrep ][ jump_row + xt + num_x * ( uy + num_y * z ) + SIZE * ( jump_col + xt + num_t * ( uy + num_u * v ) ) ]
5810                                        += 2 * ( f_dot_2dm[ d_z + z + LAS * ( d_v + v ) ] - ( f_tt + f_uu ) * one_rdm[ d_z + z + LAS * ( d_v + v ) ] );
5811                                  }
5812                               }
5813                            }
5814                         }
5815                      }
5816                   }
5817 
5818                   // SAA: - delta_ux delta_ty Gamma_{zv}
5819                   // FAA: - delta_ux delta_ty ( f_dot_2dm[ zv ] - ( f_tt + f_uu ) Gamma_{zv} )
5820                   if (( irrep_u == irrep_x ) && ( irrep_t == irrep_y ) && ( irrep_z == irrep_v )){
5821                      if ( OVLP ){
5822                         #pragma omp parallel for schedule(static)
5823                         for ( int v = 0; v < num_v; v++ ){
5824                            for ( int z = 0; z < num_z; z++ ){
5825                               for ( int ux = 0; ux < num_u; ux++ ){
5826                                  for ( int ty = 0; ty < num_t; ty++ ){
5827                                     SAA[ irrep ][ jump_row + ux + num_x * ( ty + num_y * z ) + SIZE * ( jump_col + ty + num_t * ( ux + num_u * v ) ) ]
5828                                        -= one_rdm[ d_z + z + LAS * ( d_v + v ) ];
5829                                  }
5830                               }
5831                            }
5832                         }
5833                      } else {
5834                         #pragma omp parallel for schedule(static)
5835                         for ( int v = 0; v < num_v; v++ ){
5836                            for ( int z = 0; z < num_z; z++ ){
5837                               for ( int ux = 0; ux < num_u; ux++ ){
5838                                  const double f_uu = fock->get( irrep_u, nocc_u + ux, nocc_u + ux );
5839                                  for ( int ty = 0; ty < num_t; ty++ ){
5840                                     const double f_tt = fock->get( irrep_t, nocc_t + ty, nocc_t + ty );
5841                                     FAA[ irrep ][ jump_row + ux + num_x * ( ty + num_y * z ) + SIZE * ( jump_col + ty + num_t * ( ux + num_u * v ) ) ]
5842                                        -= ( f_dot_2dm[ d_z + z + LAS * ( d_v + v ) ] - ( f_tt + f_uu ) * one_rdm[ d_z + z + LAS * ( d_v + v ) ] );
5843                                  }
5844                               }
5845                            }
5846                         }
5847                      }
5848                   }
5849 
5850                   // SCC: + delta_ut delta_xy Gamma_{zv}
5851                   // FCC: + delta_ut delta_xy ( f_dot_2dm[ zv ] - ( f_tt + f_xx ) Gamma_{zv} )
5852                   if (( irrep_u == irrep_t ) && ( irrep_x == irrep_y ) && ( irrep_z == irrep_v )){
5853                      if ( OVLP ){
5854                         #pragma omp parallel for schedule(static)
5855                         for ( int v = 0; v < num_z; v++ ){
5856                            for ( int z = 0; z < num_z; z++ ){
5857                               for ( int tu = 0; tu < num_t; tu++ ){
5858                                  for ( int xy = 0; xy < num_x; xy++ ){
5859                                     SCC[ irrep ][ jump_row + xy + num_x * ( xy + num_x * z ) + SIZE * ( jump_col + tu + num_t * ( tu + num_t * v ) ) ]
5860                                        += one_rdm[ d_z + z + LAS * ( d_v + v ) ];
5861                                  }
5862                               }
5863                            }
5864                         }
5865                      } else {
5866                         #pragma omp parallel for schedule(static)
5867                         for ( int v = 0; v < num_v; v++ ){
5868                            for ( int z = 0; z < num_z; z++ ){
5869                               for ( int tu = 0; tu < num_t; tu++ ){
5870                                  const double f_tt = fock->get( irrep_t, nocc_t + tu, nocc_t + tu );
5871                                  for ( int xy = 0; xy < num_x; xy++ ){
5872                                     const double f_xx = fock->get( irrep_x, nocc_x + xy, nocc_x + xy );
5873                                     FCC[ irrep ][ jump_row + xy + num_x * ( xy + num_x * z ) + SIZE * ( jump_col + tu + num_t * ( tu + num_t * v ) ) ]
5874                                        += ( f_dot_2dm[ d_z + z + LAS * ( d_v + v ) ] - ( f_tt + f_xx ) * one_rdm[ d_z + z + LAS * ( d_v + v ) ] );
5875                                  }
5876                               }
5877                            }
5878                         }
5879                      }
5880                   }
5881                   jump_row += num_x * num_y * num_z;
5882                }
5883             }
5884             if (( OVLP == false ) && ( fabs( IPEA ) > 0.0 )){
5885                // A: E_ti E_uv | 0 >   --->   t: excitation into,   u: excitation into, v: excitation out of
5886                // C: E_at E_uv | 0 >   --->   t: excitation out of, u: excitation into, v: excitation out of
5887                #pragma omp parallel for schedule(static)
5888                for ( int v = 0; v < num_v; v++ ){
5889                   const double gamma_vv = one_rdm[ ( d_v + v ) * ( 1 + LAS ) ];
5890                   for ( int u = 0; u < num_u; u++ ){
5891                      const double gamma_uu = one_rdm[ ( d_u + u ) * ( 1 + LAS ) ];
5892                      for ( int t = 0; t < num_t; t++ ){
5893                         const double gamma_tt = one_rdm[ ( d_t + t ) * ( 1 + LAS ) ];
5894                         const double ipea_A_tuv = 0.5 * IPEA * ( 2.0 + gamma_tt + gamma_uu - gamma_vv );
5895                         const double ipea_C_tuv = 0.5 * IPEA * ( 4.0 - gamma_tt + gamma_uu - gamma_vv );
5896                         const int ptr = ( jump_col + t + num_t * ( u + num_u * v ) ) * ( 1 + SIZE );
5897                         FAA[ irrep ][ ptr ] += ipea_A_tuv * SAA[ irrep ][ ptr ];
5898                         FCC[ irrep ][ ptr ] += ipea_C_tuv * SCC[ irrep ][ ptr ];
5899                      }
5900                   }
5901                }
5902             }
5903             jump_col += num_t * num_u * num_v;
5904          }
5905       }
5906    }
5907 
5908 }
5909 
make_DD(const bool OVLP,const double IPEA)5910 void CheMPS2::CASPT2::make_DD( const bool OVLP, const double IPEA ){
5911 
5912    /*
5913       SDD: < D(bjxy) | 1 | D(aitu) > = delta_ab delta_ij ( SDD[ It x Iu ][ xytu ] )
5914       FDD: < D(bjxy) | F | D(aitu) > = delta_ab delta_ij ( FDD[ It x Iu ][ xytu ] + ( 2 sum_k f_kk + f_aa - f_ii ) SDD[ It x Iu ][ xytu ] )
5915 
5916             SD1D1[ It x Iu ][ xytu ] = ( + 2 * Gamma_{ytxu}
5917                                          + 2 * delta_tx Gamma_{yu}
5918                                        )
5919 
5920             FD1D1[ It x Iu ][ xytu ] = ( + 2 * f_dot_3dm[ ytxu ] + ( f_xx + f_tt ) SD1D1[ It x Iu ][ xytu ]
5921                                          + 2 * delta_tx ( f_dot_2dm[ yu ] - f_tt Gamma_{yu} )
5922                                        )
5923 
5924             SD2D2[ It x Iu ][ xytu ] = ( - Gamma_{ytux}
5925                                          + 2 * delta_tx Gamma_{yu}
5926                                        )
5927 
5928             FD2D2[ It x Iu ][ xytu ] = ( - f_dot_3dm[ ytux ] + ( f_xx + f_tt ) SD2D2[ It x Iu ][ xytu ]
5929                                          + 2 * delta_tx ( f_dot_2dm[ yu ] - f_tt Gamma_{yu} )
5930                                        )
5931 
5932             SD1D2[ It x Iu ][ xytu ] = ( - Gamma_{ytxu}
5933                                          - delta_tx Gamma_{yu}
5934                                        )
5935 
5936             FD1D2[ It x Iu ][ xytu ] = ( - f_dot_3dm[ ytxu ] + ( f_xx + f_tt ) SD1D2[ It x Iu ][ xytu ]
5937                                          - delta_tx ( f_dot_2dm[ yu ] - f_tt Gamma_{yu} )
5938                                        )
5939 
5940             SD2D1[ It x Iu ][ xytu ] = ( - Gamma_{ytxu}
5941                                          - delta_tx Gamma_{yu}
5942                                        )
5943 
5944             FD2D1[ It x Iu ][ xytu ] = ( - f_dot_3dm[ ytxu ] + ( f_xx + f_tt ) SD2D1[ It x Iu ][ xytu ]
5945                                          - delta_tx ( f_dot_2dm[ yu ] - f_tt Gamma_{yu} )
5946                                        )
5947    */
5948 
5949    if ( OVLP ){ SDD = new double*[ num_irreps ]; }
5950    else {       FDD = new double*[ num_irreps ]; }
5951 
5952    const int LAS = indices->getDMRGcumulative( num_irreps );
5953 
5954    for ( int irrep = 0; irrep < num_irreps; irrep++ ){
5955 
5956       const int SIZE   = size_D[ irrep ];
5957       const int D2JUMP = SIZE / 2;
5958       if ( OVLP ){ SDD[ irrep ] = new double[ SIZE * SIZE ]; }
5959       else {       FDD[ irrep ] = new double[ SIZE * SIZE ]; }
5960 
5961       int jump_col = 0;
5962       for ( int irrep_t = 0; irrep_t < num_irreps; irrep_t++ ){
5963          const int d_t     = indices->getDMRGcumulative( irrep_t );
5964          const int num_t   = indices->getNDMRG( irrep_t );
5965          const int nocc_t  = indices->getNOCC( irrep_t );
5966          const int irrep_u = Irreps::directProd( irrep, irrep_t );
5967          const int d_u     = indices->getDMRGcumulative( irrep_u );
5968          const int num_u   = indices->getNDMRG( irrep_u );
5969          int jump_row = 0;
5970          for ( int irrep_x = 0; irrep_x < num_irreps; irrep_x++ ){
5971             const int d_x     = indices->getDMRGcumulative( irrep_x );
5972             const int num_x   = indices->getNDMRG( irrep_x );
5973             const int nocc_x  = indices->getNOCC( irrep_x );
5974             const int irrep_y = Irreps::directProd( irrep, irrep_x );
5975             const int d_y     = indices->getDMRGcumulative( irrep_y );
5976             const int num_y   = indices->getNDMRG( irrep_y );
5977 
5978             if ( OVLP ){
5979                #pragma omp parallel for schedule(static)
5980                for ( int u = 0; u < num_u; u++ ){
5981                   for ( int t = 0; t < num_t; t++ ){
5982                      for ( int y = 0; y < num_y; y++ ){
5983                         for ( int x = 0; x < num_x; x++ ){
5984                            const double gamma_ytxu = two_rdm[ d_y + y + LAS * ( d_t + t + LAS * ( d_x + x + LAS * ( d_u + u ))) ];
5985                            const double gamma_ytux = two_rdm[ d_t + t + LAS * ( d_y + y + LAS * ( d_x + x + LAS * ( d_u + u ))) ];
5986                            const int ptr = jump_row + x + num_x * y + SIZE * ( jump_col + t + num_t * u );
5987                            SDD[ irrep ][ ptr                          ] = 2 * gamma_ytxu;
5988                            SDD[ irrep ][ ptr +          SIZE * D2JUMP ] =   - gamma_ytxu;
5989                            SDD[ irrep ][ ptr + D2JUMP                 ] =   - gamma_ytxu;
5990                            SDD[ irrep ][ ptr + D2JUMP + SIZE * D2JUMP ] =   - gamma_ytux;
5991                         }
5992                      }
5993                   }
5994                }
5995             } else {
5996                #pragma omp parallel for schedule(static)
5997                for ( int u = 0; u < num_u; u++ ){
5998                   for ( int t = 0; t < num_t; t++ ){
5999                      const double f_tt = fock->get( irrep_t, nocc_t + t, nocc_t + t );
6000                      for ( int y = 0; y < num_y; y++ ){
6001                         for ( int x = 0; x < num_x; x++ ){
6002                            const double f_xx = fock->get( irrep_x, nocc_x + x, nocc_x + x );
6003                            const double f_3dm_ytxu = f_dot_3dm[ d_y + y + LAS * ( d_t + t + LAS * ( d_x + x + LAS * ( d_u + u ))) ];
6004                            const double f_3dm_ytux = f_dot_3dm[ d_t + t + LAS * ( d_y + y + LAS * ( d_x + x + LAS * ( d_u + u ))) ];
6005                            const int ptr = jump_row + x + num_x * y + SIZE * ( jump_col + t + num_t * u );
6006                            FDD[ irrep ][ ptr                          ] = 2 * f_3dm_ytxu + ( f_xx + f_tt ) * SDD[ irrep ][ ptr                          ];
6007                            FDD[ irrep ][ ptr +          SIZE * D2JUMP ] =   - f_3dm_ytxu + ( f_xx + f_tt ) * SDD[ irrep ][ ptr +          SIZE * D2JUMP ];
6008                            FDD[ irrep ][ ptr + D2JUMP                 ] =   - f_3dm_ytxu + ( f_xx + f_tt ) * SDD[ irrep ][ ptr + D2JUMP                 ];
6009                            FDD[ irrep ][ ptr + D2JUMP + SIZE * D2JUMP ] =   - f_3dm_ytux + ( f_xx + f_tt ) * SDD[ irrep ][ ptr + D2JUMP + SIZE * D2JUMP ];
6010                         }
6011                      }
6012                   }
6013                }
6014             }
6015 
6016             if (( irrep_x == irrep_t ) && ( irrep_y == irrep_u )){
6017                if ( OVLP ){
6018                   #pragma omp parallel for schedule(static)
6019                   for ( int u = 0; u < num_y; u++ ){
6020                      for ( int xt = 0; xt < num_x; xt++ ){
6021                         for ( int y = 0; y < num_y; y++ ){
6022                            const double gamma_yu = one_rdm[ d_y + y + LAS * ( d_y + u ) ];
6023                            const int ptr = jump_row + xt + num_x * y + SIZE * ( jump_col + xt + num_x * u );
6024                            SDD[ irrep ][ ptr                          ] += 2 * gamma_yu;
6025                            SDD[ irrep ][ ptr +          SIZE * D2JUMP ] -=     gamma_yu;
6026                            SDD[ irrep ][ ptr + D2JUMP                 ] -=     gamma_yu;
6027                            SDD[ irrep ][ ptr + D2JUMP + SIZE * D2JUMP ] += 2 * gamma_yu;
6028                         }
6029                      }
6030                   }
6031                } else {
6032                   #pragma omp parallel for schedule(static)
6033                   for ( int u = 0; u < num_y; u++ ){
6034                      for ( int xt = 0; xt < num_x; xt++ ){
6035                         const double f_tt = fock->get( irrep_t, nocc_t + xt, nocc_t + xt );
6036                         for ( int y = 0; y < num_y; y++ ){
6037                            const double val_yu = ( f_dot_2dm[ d_y + y + LAS * ( d_y + u ) ]
6038                                             - f_tt * one_rdm[ d_y + y + LAS * ( d_y + u ) ] );
6039                            const int ptr = jump_row + xt + num_x * y + SIZE * ( jump_col + xt + num_x * u );
6040                            FDD[ irrep ][ ptr                          ] += 2 * val_yu;
6041                            FDD[ irrep ][ ptr +          SIZE * D2JUMP ] -=     val_yu;
6042                            FDD[ irrep ][ ptr + D2JUMP                 ] -=     val_yu;
6043                            FDD[ irrep ][ ptr + D2JUMP + SIZE * D2JUMP ] += 2 * val_yu;
6044                         }
6045                      }
6046                   }
6047                }
6048             }
6049             jump_row += num_x * num_y;
6050          }
6051          if (( OVLP == false ) && ( fabs( IPEA ) > 0.0 )){
6052             // D1: E_ai E_tu | 0 >   --->   t: excitation into, u excitation out of
6053             // D2: E_ti E_au | 0 >   --->   t: excitation into, u excitation out of
6054             #pragma omp parallel for schedule(static)
6055             for ( int u = 0; u < num_u; u++ ){
6056                const double gamma_uu = one_rdm[ ( d_u + u ) * ( 1 + LAS ) ];
6057                for ( int t = 0; t < num_t; t++ ){
6058                   const double gamma_tt = one_rdm[ ( d_t + t ) * ( 1 + LAS ) ];
6059                   const double ipea_tu = 0.5 * IPEA * ( 2.0 + gamma_tt - gamma_uu );
6060                   const int ptr1 = (          jump_col + t + num_t * u ) * ( 1 + SIZE );
6061                   const int ptr2 = ( D2JUMP + jump_col + t + num_t * u ) * ( 1 + SIZE );
6062                   FDD[ irrep ][ ptr1 ] += ipea_tu * SDD[ irrep ][ ptr1 ];
6063                   FDD[ irrep ][ ptr2 ] += ipea_tu * SDD[ irrep ][ ptr2 ];
6064                }
6065             }
6066          }
6067          jump_col += num_t * num_u;
6068       }
6069    }
6070 
6071 }
6072 
make_BB_FF_singlet(const bool OVLP,const double IPEA)6073 void CheMPS2::CASPT2::make_BB_FF_singlet( const bool OVLP, const double IPEA ){
6074 
6075    /*
6076       | SB_tiuj > = ( E_ti E_uj + E_tj E_ui ) / sqrt( 1 + delta_ij ) | 0 >  with  i <= j and t <= u
6077 
6078       SBB singlet: < SB_xkyl | 1 | SB_tiuj > = 2 delta_ik delta_jl ( SBB_singlet[ Itu ][ xytu ] )
6079       FBB singlet: < SB_xkyl | F | SB_tiuj > = 2 delta_ik delta_jl ( FBB_singlet[ Itu ][ xytu ] + ( 2 sum_n f_nn - f_ii - f_jj ) * SBB_singlet[ Itu ][ xytu ] )
6080 
6081             SBB_singlet[ Itu ][ xytu ] = ( + Gamma_{utyx}
6082                                            + Gamma_{utxy}
6083                                            + 2 ( delta_uy delta_tx + delta_ux delta_ty )
6084                                            - delta_uy Gamma_{tx}
6085                                            - delta_tx Gamma_{uy}
6086                                            - delta_ux Gamma_{ty}
6087                                            - delta_ty Gamma_{ux}
6088                                          )
6089 
6090             FBB_singlet[ Itu ][ xytu ] = ( + f_dot_3dm[ utyx ]
6091                                            + f_dot_3dm[ tuyx ]
6092                                            + ( f_tt + f_uu + f_xx + f_yy ) SBB_singlet[ xytu ]
6093                                            + 2 ( delta_uy delta_tx + delta_ux delta_ty ) ( f_dot_1dm - f_tt - f_uu )
6094                                            -   delta_uy ( f_dot_2dm[ tx ] - f_uu Gamma_{tx} )
6095                                            -   delta_tx ( f_dot_2dm[ uy ] - f_tt Gamma_{uy} )
6096                                            -   delta_ux ( f_dot_2dm[ ty ] - f_uu Gamma_{ty} )
6097                                            -   delta_ty ( f_dot_2dm[ ux ] - f_tt Gamma_{ux} )
6098                                          )
6099 
6100       | SF_atbu > = ( E_at E_bu + E_bt E_au ) / sqrt( 1 + delta_ab ) | 0 >  with  a <= b and t <= u
6101 
6102       SFF singlet: < SF_cxdy | 1 | SF_atbu > = 2 delta_ac delta_bd ( SFF_singlet[ Itu ][ xytu ] )
6103       FFF singlet: < SF_cxdy | F | SF_atbu > = 2 delta_ac delta_bd ( FFF_singlet[ Iab ][ xytu ] + ( 2 sum_n f_nn + f_aa + f_bb ) * SFF_singlet[ Iab ][ xytu ] )
6104 
6105             SFF_singlet[ Itu ][ xytu ] = ( + Gamma_{yxut}
6106                                            + Gamma_{yxtu}
6107                                          )
6108 
6109             FFF_singlet[ Itu ][ xytu ] = ( + f_dot_3dm[ yxut ]
6110                                            + f_dot_3dm[ yxtu ]
6111                                          )
6112    */
6113 
6114    if ( OVLP ){ SBB_singlet = new double*[ num_irreps ];
6115                 SFF_singlet = new double*[ num_irreps ]; }
6116    else {       FBB_singlet = new double*[ num_irreps ];
6117                 FFF_singlet = new double*[ num_irreps ]; }
6118 
6119    const int LAS = indices->getDMRGcumulative( num_irreps );
6120 
6121    { // First do irrep == Iia x Ijb == 0 -->  It == Iu and Ix == Iy
6122       assert( size_B_singlet[ 0 ] == size_F_singlet[ 0 ] ); // At construction
6123       const int SIZE = size_B_singlet[ 0 ];
6124       if ( OVLP ){ SBB_singlet[ 0 ] = new double[ SIZE * SIZE ];
6125                    SFF_singlet[ 0 ] = new double[ SIZE * SIZE ]; }
6126       else {       FBB_singlet[ 0 ] = new double[ SIZE * SIZE ];
6127                    FFF_singlet[ 0 ] = new double[ SIZE * SIZE ]; }
6128 
6129       int jump_col = 0;
6130       for ( int irrep_ut = 0; irrep_ut < num_irreps; irrep_ut++ ){
6131          const int d_ut    = indices->getDMRGcumulative( irrep_ut );
6132          const int num_ut  = indices->getNDMRG( irrep_ut );
6133          const int nocc_ut = indices->getNOCC( irrep_ut );
6134          int jump_row = 0;
6135          for ( int irrep_xy = 0; irrep_xy < num_irreps; irrep_xy++ ){
6136             const int d_xy    = indices->getDMRGcumulative( irrep_xy );
6137             const int num_xy  = indices->getNDMRG( irrep_xy );
6138             const int nocc_xy = indices->getNOCC( irrep_xy );
6139             const int shift   = jump_row + SIZE * jump_col;
6140 
6141             // SBB singlet: + Gamma_{utyx} + Gamma_{utxy}
6142             // SFF singlet: + Gamma_{yxut} + Gamma_{yxtu}
6143             // FBB singlet: + f_dot_3dm[ utyx ] + f_dot_3dm[ tuyx ] + ( f_tt + f_uu + f_xx + f_yy ) SBB_singlet[ xytu ]
6144             // FFF singlet: + f_dot_3dm[ yxut ] + f_dot_3dm[ yxtu ]
6145             if ( OVLP ){
6146                for ( int t = 0; t < num_ut; t++ ){
6147                   for ( int u = t; u < num_ut; u++ ){ // 0 <= t <= u < num_ut
6148                      for ( int x = 0; x < num_xy; x++ ){
6149                         for ( int y = x; y < num_xy; y++ ){ // 0 <= x <= y < num_xy
6150                            const double value = ( two_rdm[ d_xy + x + LAS * ( d_xy + y + LAS * ( d_ut + t + LAS * ( d_ut + u ))) ]
6151                                                 + two_rdm[ d_xy + x + LAS * ( d_xy + y + LAS * ( d_ut + u + LAS * ( d_ut + t ))) ] );
6152                            const int ptr = shift + x + ( y * ( y + 1 ) ) / 2 + SIZE * ( t + ( u * ( u + 1 ) ) / 2 );
6153                            SBB_singlet[ 0 ][ ptr ] = value;
6154                            SFF_singlet[ 0 ][ ptr ] = value;
6155                         }
6156                      }
6157                   }
6158                }
6159             } else {
6160                for ( int t = 0; t < num_ut; t++ ){
6161                   const double f_tt = fock->get( irrep_ut, nocc_ut + t, nocc_ut + t );
6162                   for ( int u = t; u < num_ut; u++ ){ // 0 <= t <= u < num_ut
6163                      const double f_uu = fock->get( irrep_ut, nocc_ut + u, nocc_ut + u );
6164                      for ( int x = 0; x < num_xy; x++ ){
6165                         const double f_xx = fock->get( irrep_xy, nocc_xy + x, nocc_xy + x );
6166                         for ( int y = x; y < num_xy; y++ ){ // 0 <= x <= y < num_xy
6167                            const double f_yy = fock->get( irrep_xy, nocc_xy + y, nocc_xy + y );
6168                            const int ptr = shift + x + ( y * ( y + 1 ) ) / 2 + SIZE * ( t + ( u * ( u + 1 ) ) / 2 );
6169                            const double fdotsum = ( f_dot_3dm[ d_xy + x + LAS * ( d_xy + y + LAS * ( d_ut + t + LAS * ( d_ut + u ))) ]
6170                                                   + f_dot_3dm[ d_xy + x + LAS * ( d_xy + y + LAS * ( d_ut + u + LAS * ( d_ut + t ))) ] );
6171                            FBB_singlet[ 0 ][ ptr ] = fdotsum + ( f_tt + f_uu + f_xx + f_yy ) * SBB_singlet[ 0 ][ ptr ];
6172                            FFF_singlet[ 0 ][ ptr ] = fdotsum;
6173                         }
6174                      }
6175                   }
6176                }
6177             }
6178 
6179             if ( irrep_ut == irrep_xy ){ // All four irreps are equal: num_ut = num_xy and d_ut = d_xy
6180 
6181                // SBB singlet: + 2 ( delta_uy delta_tx + delta_ux delta_ty )
6182                // FBB singlet: + 2 ( delta_uy delta_tx + delta_ux delta_ty ) ( f_dot_1dm - ( f_tt + f_uu ) )
6183                if ( OVLP ){
6184                   for ( int t = 0; t < num_ut; t++ ){
6185                      SBB_singlet[ 0 ][ shift + t + ( t * ( t + 1 ) ) / 2 + SIZE * ( t + ( t * ( t + 1 ) ) / 2 ) ] += 4.0;
6186                      for ( int u = t+1; u < num_ut; u++ ){
6187                         SBB_singlet[ 0 ][ shift + t + ( u * ( u + 1 ) ) / 2 + SIZE * ( t + ( u * ( u + 1 ) ) / 2 ) ] += 2.0;
6188                      }
6189                   }
6190                } else {
6191                   for ( int t = 0; t < num_ut; t++ ){
6192                      const double f_tt = fock->get( irrep_ut, nocc_ut + t, nocc_ut + t );
6193                      FBB_singlet[ 0 ][ shift + t + ( t * ( t + 1 ) ) / 2 + SIZE * ( t + ( t * ( t + 1 ) ) / 2 ) ] += 4 * ( f_dot_1dm - 2 * f_tt );
6194                      for ( int u = t+1; u < num_ut; u++ ){
6195                         const double f_uu = fock->get( irrep_ut, nocc_ut + u, nocc_ut + u );
6196                         FBB_singlet[ 0 ][ shift + t + ( u * ( u + 1 ) ) / 2 + SIZE * ( t + ( u * ( u + 1 ) ) / 2 ) ] += 2 * ( f_dot_1dm - f_tt - f_uu );
6197                      }
6198                   }
6199                }
6200 
6201                // SBB singlet: - delta_uy Gamma_{tx}
6202                // FBB singlet: - delta_uy ( f_dot_2dm[ tx ] - f_uu Gamma_{tx} )
6203                if ( OVLP ){
6204                   for ( int uy = 0; uy < num_ut; uy++ ){
6205                      for ( int t = 0; t <= uy; t++ ){ // 0 <= t <= uy < num_ut
6206                         for ( int x = 0; x <= uy; x++ ){ // 0 <= x <= uy < num_ut
6207                            const double gamma_tx = one_rdm[ d_ut + t + LAS * ( d_ut + x ) ];
6208                            SBB_singlet[ 0 ][ shift + x + ( uy * ( uy + 1 ) ) / 2 + SIZE * ( t + ( uy * ( uy + 1 ) ) / 2 ) ] -= gamma_tx;
6209                         }
6210                      }
6211                   }
6212                } else {
6213                   for ( int uy = 0; uy < num_ut; uy++ ){
6214                      const double f_uu = fock->get( irrep_ut, nocc_ut + uy, nocc_ut + uy );
6215                      for ( int t = 0; t <= uy; t++ ){ // 0 <= t <= uy < num_ut
6216                         for ( int x = 0; x <= uy; x++ ){ // 0 <= x <= uy < num_ut
6217                            const double val_tx = ( f_dot_2dm[ d_ut + t + LAS * ( d_ut + x ) ]
6218                                             - f_uu * one_rdm[ d_ut + t + LAS * ( d_ut + x ) ] );
6219                            FBB_singlet[ 0 ][ shift + x + ( uy * ( uy + 1 ) ) / 2 + SIZE * ( t + ( uy * ( uy + 1 ) ) / 2 ) ] -= val_tx;
6220                         }
6221                      }
6222                   }
6223                }
6224 
6225                // SBB singlet: - delta_tx Gamma_{uy}
6226                // FBB singlet: - delta_tx ( f_dot_2dm[ uy ] - f_tt Gamma_{uy} )
6227                if ( OVLP ){
6228                   for ( int tx = 0; tx < num_ut; tx++ ){
6229                      for ( int u = tx; u < num_ut; u++ ){ // 0 <= tx <= u < num_ut
6230                         for ( int y = tx; y < num_ut; y++ ){ // 0 <= tx <= y < num_xy = num_ut
6231                            const double gamma_uy = one_rdm[ d_ut + u + LAS * ( d_ut + y ) ];
6232                            SBB_singlet[ 0 ][ shift + tx + ( y * ( y + 1 ) ) / 2 + SIZE * ( tx + ( u * ( u + 1 ) ) / 2 ) ] -= gamma_uy;
6233                         }
6234                      }
6235                   }
6236                } else {
6237                   for ( int tx = 0; tx < num_ut; tx++ ){
6238                      const double f_tt = fock->get( irrep_ut, nocc_ut + tx, nocc_ut + tx );
6239                      for ( int u = tx; u < num_ut; u++ ){ // 0 <= tx <= u < num_ut
6240                         for ( int y = tx; y < num_ut; y++ ){ // 0 <= tx <= y < num_xy = num_ut
6241                            const double val_uy = ( f_dot_2dm[ d_ut + u + LAS * ( d_ut + y ) ]
6242                                             - f_tt * one_rdm[ d_ut + u + LAS * ( d_ut + y ) ] );
6243                            FBB_singlet[ 0 ][ shift + tx + ( y * ( y + 1 ) ) / 2 + SIZE * ( tx + ( u * ( u + 1 ) ) / 2 ) ] -= val_uy;
6244                         }
6245                      }
6246                   }
6247                }
6248 
6249                // SBB singlet: - delta_ux Gamma_{ty}
6250                // FBB singlet: - delta_ux ( f_dot_2dm[ ty ] - f_uu Gamma_{ty} )
6251                if ( OVLP ){
6252                   for ( int ux = 0; ux < num_ut; ux++ ){
6253                      for ( int t = 0; t <= ux; t++ ){ // 0 <= t <= ux < num_ut
6254                         for ( int y = ux; y < num_ut; y++ ){ // 0 <= ux <= y < num_xy = num_ut
6255                            const double gamma_ty = one_rdm[ d_ut + t + LAS * ( d_ut + y ) ];
6256                            SBB_singlet[ 0 ][ shift + ux + ( y * ( y + 1 ) ) / 2 + SIZE * ( t + ( ux * ( ux + 1 ) ) / 2 ) ] -= gamma_ty;
6257                         }
6258                      }
6259                   }
6260                } else {
6261                   for ( int ux = 0; ux < num_ut; ux++ ){
6262                      const double f_uu = fock->get( irrep_ut, nocc_ut + ux, nocc_ut + ux );
6263                      for ( int t = 0; t <= ux; t++ ){ // 0 <= t <= ux < num_ut
6264                         for ( int y = ux; y < num_ut; y++ ){ // 0 <= ux <= y < num_xy = num_ut
6265                            const double val_ty = ( f_dot_2dm[ d_ut + t + LAS * ( d_ut + y ) ]
6266                                             - f_uu * one_rdm[ d_ut + t + LAS * ( d_ut + y ) ] );
6267                            FBB_singlet[ 0 ][ shift + ux + ( y * ( y + 1 ) ) / 2 + SIZE * ( t + ( ux * ( ux + 1 ) ) / 2 ) ] -= val_ty;
6268                         }
6269                      }
6270                   }
6271                }
6272 
6273                // SBB singlet: - delta_ty Gamma_{ux}
6274                // FBB singlet: - delta_ty ( f_dot_2dm[ ux ] - f_tt Gamma_{ux} )
6275                if ( OVLP ){
6276                   for ( int ty = 0; ty < num_ut; ty++ ){
6277                      for ( int u = ty; u < num_ut; u++ ){ // 0 <= ty <= u < num_ut
6278                         for ( int x = 0; x <= ty; x++ ){ // 0 <= x <= ty < num_ut
6279                            const double gamma_ux = one_rdm[ d_ut + u + LAS * ( d_ut + x ) ];
6280                            SBB_singlet[ 0 ][ shift + x + ( ty * ( ty + 1 ) ) / 2 + SIZE * ( ty + ( u * ( u + 1 ) ) / 2 ) ] -= gamma_ux;
6281                         }
6282                      }
6283                   }
6284                } else {
6285                   for ( int ty = 0; ty < num_ut; ty++ ){
6286                      const double f_tt = fock->get( irrep_ut, nocc_ut + ty, nocc_ut + ty );
6287                      for ( int u = ty; u < num_ut; u++ ){ // 0 <= ty <= u < num_ut
6288                         for ( int x = 0; x <= ty; x++ ){ // 0 <= x <= ty < num_ut
6289                            const double val_ux = ( f_dot_2dm[ d_ut + u + LAS * ( d_ut + x ) ]
6290                                             - f_tt * one_rdm[ d_ut + u + LAS * ( d_ut + x ) ] );
6291                            FBB_singlet[ 0 ][ shift + x + ( ty * ( ty + 1 ) ) / 2 + SIZE * ( ty + ( u * ( u + 1 ) ) / 2 ) ] -= val_ux;
6292                         }
6293                      }
6294                   }
6295                }
6296             }
6297             jump_row += ( num_xy * ( num_xy + 1 ) ) / 2;
6298          }
6299          if (( OVLP == false ) && ( fabs( IPEA ) > 0.0 )){
6300             // B: E_ti E_uj | 0 >   --->   tu: excitation into
6301             // F: E_at E_bu | 0 >   --->   tu: excitation out of
6302             for ( int u = 0; u < num_ut; u++ ){
6303                const double gamma_uu = one_rdm[ ( d_ut + u ) * ( 1 + LAS ) ];
6304                for ( int t = 0; t <= u; t++ ){ // 0 <= t <= u < num_ut
6305                   const double gamma_tt = one_rdm[ ( d_ut + t ) * ( 1 + LAS ) ];
6306                   const double ipea_B_tu = 0.5 * IPEA * ( gamma_tt + gamma_uu );
6307                   const double ipea_F_tu = 0.5 * IPEA * ( 4.0 - gamma_tt - gamma_uu );
6308                   const int ptr = ( jump_col + t + ( u * ( u + 1 ) ) / 2 ) * ( 1 + SIZE );
6309                   FBB_singlet[ 0 ][ ptr ] += ipea_B_tu * SBB_singlet[ 0 ][ ptr ];
6310                   FFF_singlet[ 0 ][ ptr ] += ipea_F_tu * SFF_singlet[ 0 ][ ptr ];
6311                }
6312             }
6313          }
6314          jump_col += ( num_ut * ( num_ut + 1 ) ) / 2;
6315       }
6316    }
6317 
6318    for ( int irrep = 1; irrep < num_irreps; irrep++ ){ // Then do irrep == Iia x Ijb != 0 -->  It != Iu and Ix != Iy
6319 
6320       assert( size_B_singlet[ irrep ] == size_F_singlet[ irrep ] ); // At construction
6321       const int SIZE = size_B_singlet[ irrep ];
6322       if ( OVLP ){ SBB_singlet[ irrep ] = new double[ SIZE * SIZE ];
6323                    SFF_singlet[ irrep ] = new double[ SIZE * SIZE ]; }
6324       else {       FBB_singlet[ irrep ] = new double[ SIZE * SIZE ];
6325                    FFF_singlet[ irrep ] = new double[ SIZE * SIZE ]; }
6326 
6327       int jump_col = 0;
6328       for ( int irrep_t = 0; irrep_t < num_irreps; irrep_t++ ){
6329          const int irrep_u = Irreps::directProd( irrep, irrep_t );
6330          if ( irrep_t < irrep_u ){
6331             const int d_t    = indices->getDMRGcumulative( irrep_t );
6332             const int num_t  = indices->getNDMRG( irrep_t );
6333             const int nocc_t = indices->getNOCC( irrep_t );
6334             const int d_u    = indices->getDMRGcumulative( irrep_u );
6335             const int num_u  = indices->getNDMRG( irrep_u );
6336             const int nocc_u = indices->getNOCC( irrep_u );
6337             int jump_row = 0;
6338             for ( int irrep_x = 0; irrep_x < num_irreps; irrep_x++ ){
6339                const int irrep_y = Irreps::directProd( irrep, irrep_x );
6340                if ( irrep_x < irrep_y ){
6341                   const int d_x    = indices->getDMRGcumulative( irrep_x );
6342                   const int num_x  = indices->getNDMRG( irrep_x );
6343                   const int nocc_x = indices->getNOCC( irrep_x );
6344                   const int d_y    = indices->getDMRGcumulative( irrep_y );
6345                   const int num_y  = indices->getNDMRG( irrep_y );
6346                   const int nocc_y = indices->getNOCC( irrep_y );
6347                   const int shift  = jump_row + SIZE * jump_col;
6348 
6349                   // SBB singlet: + Gamma_{utyx} + Gamma_{utxy}
6350                   // SFF singlet: + Gamma_{yxut} + Gamma_{yxtu}
6351                   // FBB singlet: + f_dot_3dm[ utyx ] + f_dot_3dm[ tuyx ] + ( f_tt + f_uu + f_xx + f_yy ) SBB_singlet[ xytu ]
6352                   // FFF singlet: + f_dot_3dm[ yxut ] + f_dot_3dm[ yxtu ]
6353                   if ( OVLP ){
6354                      for ( int t = 0; t < num_t; t++ ){
6355                         for ( int u = 0; u < num_u; u++ ){
6356                            for ( int x = 0; x < num_x; x++ ){
6357                               for ( int y = 0; y < num_y; y++ ){
6358                                  const double value = ( two_rdm[ d_x + x + LAS * ( d_y + y + LAS * ( d_t + t + LAS * ( d_u + u ))) ]
6359                                                       + two_rdm[ d_x + x + LAS * ( d_y + y + LAS * ( d_u + u + LAS * ( d_t + t ))) ] );
6360                                  const int ptr = shift + x + num_x * y + SIZE * ( t + num_t * u );
6361                                  SBB_singlet[ irrep ][ ptr ] = value;
6362                                  SFF_singlet[ irrep ][ ptr ] = value;
6363                               }
6364                            }
6365                         }
6366                      }
6367                   } else {
6368                      for ( int t = 0; t < num_t; t++ ){
6369                         const double f_tt = fock->get( irrep_t, nocc_t + t, nocc_t + t );
6370                         for ( int u = 0; u < num_u; u++ ){
6371                            const double f_uu = fock->get( irrep_u, nocc_u + u, nocc_u + u );
6372                            for ( int x = 0; x < num_x; x++ ){
6373                               const double f_xx = fock->get( irrep_x, nocc_x + x, nocc_x + x );
6374                               for ( int y = 0; y < num_y; y++ ){
6375                                  const double f_yy = fock->get( irrep_y, nocc_y + y, nocc_y + y );
6376                                  const int ptr = shift + x + num_x * y + SIZE * ( t + num_t * u );
6377                                  const double fdotsum = ( f_dot_3dm[ d_x + x + LAS * ( d_y + y + LAS * ( d_t + t + LAS * ( d_u + u ))) ]
6378                                                         + f_dot_3dm[ d_x + x + LAS * ( d_y + y + LAS * ( d_u + u + LAS * ( d_t + t ))) ] );
6379                                  FBB_singlet[ irrep ][ ptr ] = fdotsum + ( f_xx + f_yy + f_tt + f_uu ) * SBB_singlet[ irrep ][ ptr ];
6380                                  FFF_singlet[ irrep ][ ptr ] = fdotsum;
6381                               }
6382                            }
6383                         }
6384                      }
6385                   }
6386 
6387                   if (( irrep_u == irrep_y ) && ( irrep_t == irrep_x )){ // num_t == num_x  and  num_u == num_y
6388 
6389                      // SBB singlet: + 2 delta_uy delta_tx
6390                      // FBB singlet: + 2 delta_uy delta_tx ( f_dot_1dm - f_tt - f_uu )
6391                      if ( OVLP ){
6392                         for ( int xt = 0; xt < num_x; xt++ ){
6393                            for ( int yu = 0; yu < num_y; yu++ ){
6394                               SBB_singlet[ irrep ][ shift + xt + num_x * yu + SIZE * ( xt + num_x * yu ) ] += 2.0;
6395                            }
6396                         }
6397                      } else {
6398                         for ( int xt = 0; xt < num_x; xt++ ){
6399                            const double f_tt = fock->get( irrep_x, nocc_x + xt, nocc_x + xt );
6400                            for ( int yu = 0; yu < num_y; yu++ ){
6401                               const double f_uu = fock->get( irrep_y, nocc_y + yu, nocc_y + yu );
6402                               FBB_singlet[ irrep ][ shift + xt + num_x * yu + SIZE * ( xt + num_x * yu ) ] += 2 * ( f_dot_1dm - f_tt - f_uu );
6403                            }
6404                         }
6405                      }
6406 
6407                      // SBB singlet: - delta_tx Gamma_{uy}
6408                      // FBB singlet: - delta_tx ( f_dot_2dm[ uy ] - f_tt Gamma_{uy} )
6409                      if ( OVLP ){
6410                         for ( int xt = 0; xt < num_x; xt++ ){
6411                            for ( int u = 0; u < num_y; u++ ){
6412                               for ( int y = 0; y < num_y; y++ ){
6413                                  const double gamma_uy = one_rdm[ d_u + u + LAS * ( d_u + y ) ];
6414                                  SBB_singlet[ irrep ][ shift + xt + num_x * y + SIZE * ( xt + num_x * u ) ] -= gamma_uy;
6415                               }
6416                            }
6417                         }
6418                      } else {
6419                         for ( int xt = 0; xt < num_x; xt++ ){
6420                            const double f_tt = fock->get( irrep_x, nocc_x + xt, nocc_x + xt );
6421                            for ( int u = 0; u < num_y; u++ ){
6422                               for ( int y = 0; y < num_y; y++ ){
6423                                  const double val_uy = ( f_dot_2dm[ d_u + u + LAS * ( d_u + y ) ]
6424                                                   - f_tt * one_rdm[ d_u + u + LAS * ( d_u + y ) ] );
6425                                  FBB_singlet[ irrep ][ shift + xt + num_x * y + SIZE * ( xt + num_x * u ) ] -= val_uy;
6426                               }
6427                            }
6428                         }
6429                      }
6430 
6431                      // SBB singlet: - delta_uy Gamma_{tx}
6432                      // FBB singlet: - delta_uy ( f_dot_2dm[ tx ] - f_uu Gamma_{tx} )
6433                      if ( OVLP ){
6434                         for ( int yu = 0; yu < num_y; yu++ ){
6435                            for ( int t = 0; t < num_x; t++ ){
6436                               for ( int x = 0; x < num_x; x++ ){
6437                                  const double gamma_tx = one_rdm[ d_t + t + LAS * ( d_t + x ) ];
6438                                  SBB_singlet[ irrep ][ shift + x + num_x * yu + SIZE * ( t + num_t * yu ) ] -= gamma_tx;
6439                               }
6440                            }
6441                         }
6442                      } else {
6443                         for ( int yu = 0; yu < num_y; yu++ ){
6444                            const double f_uu = fock->get( irrep_y, nocc_y + yu, nocc_y + yu );
6445                            for ( int t = 0; t < num_x; t++ ){
6446                               for ( int x = 0; x < num_x; x++ ){
6447                                  const double val_tx = ( f_dot_2dm[ d_t + t + LAS * ( d_t + x ) ]
6448                                                   - f_uu * one_rdm[ d_t + t + LAS * ( d_t + x ) ] );
6449                                  FBB_singlet[ irrep ][ shift + x + num_x * yu + SIZE * ( t + num_t * yu ) ] -= val_tx;
6450                               }
6451                            }
6452                         }
6453                      }
6454                   }
6455                   jump_row += num_x * num_y;
6456                }
6457             }
6458             if (( OVLP == false ) && ( fabs( IPEA ) > 0.0 )){
6459                // B: E_ti E_uj | 0 >   --->   tu: excitation into
6460                // F: E_at E_bu | 0 >   --->   tu: excitation out of
6461                for ( int u = 0; u < num_u; u++ ){
6462                   const double gamma_uu = one_rdm[ ( d_u + u ) * ( 1 + LAS ) ];
6463                   for ( int t = 0; t < num_t; t++ ){
6464                      const double gamma_tt = one_rdm[ ( d_t + t ) * ( 1 + LAS ) ];
6465                      const double ipea_B_tu = 0.5 * IPEA * ( gamma_tt + gamma_uu );
6466                      const double ipea_F_tu = 0.5 * IPEA * ( 4.0 - gamma_tt - gamma_uu );
6467                      const int ptr = ( jump_col + t + num_t * u ) * ( 1 + SIZE );
6468                      FBB_singlet[ irrep ][ ptr ] += ipea_B_tu * SBB_singlet[ irrep ][ ptr ];
6469                      FFF_singlet[ irrep ][ ptr ] += ipea_F_tu * SFF_singlet[ irrep ][ ptr ];
6470                   }
6471                }
6472             }
6473             jump_col += num_t * num_u;
6474          }
6475       }
6476    }
6477 
6478 }
6479 
make_BB_FF_triplet(const bool OVLP,const double IPEA)6480 void CheMPS2::CASPT2::make_BB_FF_triplet( const bool OVLP, const double IPEA ){
6481 
6482    /*
6483       | TB_tiuj > = ( E_ti E_uj - E_tj E_ui ) / sqrt( 1 + delta_ij ) | 0 >  with  i < j and t < u
6484 
6485       SBB triplet: < TB_xkyl | 1 | TB_tiuj > = 2 delta_ik delta_jl ( SBB_triplet[ Itu ][ xytu ] )
6486       FBB triplet: < TB_xkyl | F | TB_tiuj > = 2 delta_ik delta_jl ( FBB_triplet[ Itu ][ xytu ] + ( 2 sum_n f_nn - f_ii - f_jj ) * SBB_triplet[ Itu ][ xytu ] )
6487 
6488             SBB_triplet[ Itu ][ xytu ] = ( + Gamma_{utyx}
6489                                            - Gamma_{utxy}
6490                                            + 6 delta_uy delta_tx
6491                                            - 6 delta_ux delta_ty
6492                                            - 3 delta_uy Gamma_{tx}
6493                                            - 3 delta_tx Gamma_{uy}
6494                                            + 3 delta_ux Gamma_{ty}
6495                                            + 3 delta_ty Gamma_{ux}
6496                                          )
6497 
6498             FBB_triplet[ Itu ][ xytu ] = ( + f_dot_3dm[ utyx ]
6499                                            - f_dot_3dm[ tuyx ]
6500                                            + ( f_tt + f_uu + f_xx + f_yy ) SBB_triplet[ xytu ]
6501                                            + 6 delta_uy delta_tx ( f_dot_1dm - f_tt - f_uu )
6502                                            - 6 delta_ux delta_ty ( f_dot_1dm - f_tt - f_uu )
6503                                            - 3 delta_uy ( f_dot_2dm[ tx ] - f_uu Gamma_{tx} )
6504                                            - 3 delta_tx ( f_dot_2dm[ uy ] - f_tt Gamma_{uy} )
6505                                            + 3 delta_ux ( f_dot_2dm[ ty ] - f_uu Gamma_{ty} )
6506                                            + 3 delta_ty ( f_dot_2dm[ ux ] - f_tt Gamma_{ux} )
6507                                          )
6508 
6509       | TF_atbu > = ( E_at E_bu - E_bt E_au ) / sqrt( 1 + delta_ab ) | 0 >  with  a < b and t < u
6510 
6511       SFF triplet: < TF_cxdy | 1 | TF_atbu > = 2 delta_ac delta_bd ( SFF_triplet[ Itu ][ xytu ] )
6512       FFF triplet: < TF_cxdy | F | TF_atbu > = 2 delta_ac delta_bd ( FFF_triplet[ Itu ][ xytu ] + ( 2 sum_n f_nn + f_aa + f_bb ) * SFF_triplet[ Itu ][ xytu ] )
6513 
6514             SFF_triplet[ Itu ][ xytu ] = ( + Gamma_{yxut}
6515                                            - Gamma_{yxtu}
6516                                          )
6517 
6518             FFF_triplet[ Itu ][ xytu ] = ( + f_dot_3dm[ yxut ]
6519                                            - f_dot_3dm[ yxtu ]
6520                                          )
6521    */
6522 
6523    if ( OVLP ){ SBB_triplet = new double*[ num_irreps ];
6524                 SFF_triplet = new double*[ num_irreps ]; }
6525    else {       FBB_triplet = new double*[ num_irreps ];
6526                 FFF_triplet = new double*[ num_irreps ]; }
6527 
6528    const int LAS = indices->getDMRGcumulative( num_irreps );
6529 
6530    { // First do irrep == Iia x Ijb == 0 -->  It == Iu and Ix == Iy
6531       assert( size_B_triplet[ 0 ] == size_F_triplet[ 0 ] ); // At construction
6532       const int SIZE = size_B_triplet[ 0 ];
6533       if ( OVLP ){ SBB_triplet[ 0 ] = new double[ SIZE * SIZE ];
6534                    SFF_triplet[ 0 ] = new double[ SIZE * SIZE ]; }
6535       else {       FBB_triplet[ 0 ] = new double[ SIZE * SIZE ];
6536                    FFF_triplet[ 0 ] = new double[ SIZE * SIZE ]; }
6537 
6538       int jump_col = 0;
6539       for ( int irrep_ut = 0; irrep_ut < num_irreps; irrep_ut++ ){
6540          const int d_ut    = indices->getDMRGcumulative( irrep_ut );
6541          const int num_ut  = indices->getNDMRG( irrep_ut );
6542          const int nocc_ut = indices->getNOCC( irrep_ut );
6543          int jump_row = 0;
6544          for ( int irrep_xy = 0; irrep_xy < num_irreps; irrep_xy++ ){
6545             const int d_xy    = indices->getDMRGcumulative( irrep_xy );
6546             const int num_xy  = indices->getNDMRG( irrep_xy );
6547             const int nocc_xy = indices->getNOCC( irrep_xy );
6548             const int shift   = jump_row + SIZE * jump_col;
6549 
6550             // SBB triplet: + Gamma_{utyx} - Gamma_{utxy}
6551             // SFF triplet: + Gamma_{yxut} - Gamma_{yxtu}
6552             // FBB triplet: + f_dot_3dm[ utyx ] - f_dot_3dm[ tuyx ] + ( f_tt + f_uu + f_xx + f_yy ) SBB_triplet[ xytu ]
6553             // FFF triplet: + f_dot_3dm[ yxut ] - f_dot_3dm[ yxtu ]
6554             if ( OVLP ){
6555                for ( int t = 0; t < num_ut; t++ ){
6556                   for ( int u = t+1; u < num_ut; u++ ){ // 0 <= t < u < num_ut
6557                      for ( int x = 0; x < num_xy; x++ ){
6558                         for ( int y = x+1; y < num_xy; y++ ){ // 0 <= x < y < num_xy
6559                            const int ptr = shift + x + ( y * ( y - 1 ) ) / 2 + SIZE * ( t + ( u * ( u - 1 ) ) / 2 );
6560                            const double value = ( two_rdm[ d_xy + x + LAS * ( d_xy + y + LAS * ( d_ut + t + LAS * ( d_ut + u ))) ]
6561                                                 - two_rdm[ d_xy + x + LAS * ( d_xy + y + LAS * ( d_ut + u + LAS * ( d_ut + t ))) ] );
6562                            SBB_triplet[ 0 ][ ptr ] = value;
6563                            SFF_triplet[ 0 ][ ptr ] = value;
6564                         }
6565                      }
6566                   }
6567                }
6568             } else {
6569                for ( int t = 0; t < num_ut; t++ ){
6570                   const double f_tt = fock->get( irrep_ut, nocc_ut + t, nocc_ut + t );
6571                   for ( int u = t+1; u < num_ut; u++ ){ // 0 <= t < u < num_ut
6572                      const double f_uu = fock->get( irrep_ut, nocc_ut + u, nocc_ut + u );
6573                      for ( int x = 0; x < num_xy; x++ ){
6574                         const double f_xx = fock->get( irrep_xy, nocc_xy + x, nocc_xy + x );
6575                         for ( int y = x+1; y < num_xy; y++ ){ // 0 <= x < y < num_xy
6576                            const double f_yy = fock->get( irrep_xy, nocc_xy + y, nocc_xy + y );
6577                            const int ptr = shift + x + ( y * ( y - 1 ) ) / 2 + SIZE * ( t + ( u * ( u - 1 ) ) / 2 );
6578                            const double fdotdiff = ( f_dot_3dm[ d_xy + x + LAS * ( d_xy + y + LAS * ( d_ut + t + LAS * ( d_ut + u ))) ]
6579                                                    - f_dot_3dm[ d_xy + x + LAS * ( d_xy + y + LAS * ( d_ut + u + LAS * ( d_ut + t ))) ] );
6580                            FBB_triplet[ 0 ][ ptr ] = fdotdiff + ( f_tt + f_uu + f_xx + f_yy ) * SBB_triplet[ 0 ][ ptr ];
6581                            FFF_triplet[ 0 ][ ptr ] = fdotdiff;
6582                         }
6583                      }
6584                   }
6585                }
6586             }
6587 
6588             if ( irrep_ut == irrep_xy ){ // All four irreps are equal --> num_ut == num_xy
6589 
6590                // SBB triplet: + 6 ( delta_uy delta_tx - delta_ux delta_ty )
6591                // FBB triplet: + 6 ( delta_uy delta_tx - delta_ux delta_ty ) ( f_dot_1dm - f_tt - f_uu )
6592                if ( OVLP ){
6593                   for ( int tx = 0; tx < num_ut; tx++ ){
6594                      for ( int uy = tx+1; uy < num_ut; uy++ ){
6595                         SBB_triplet[ 0 ][ shift + tx + ( uy * ( uy - 1 ) ) / 2 + SIZE * ( tx + ( uy * ( uy - 1 ) ) / 2 ) ] += 6.0;
6596                      }
6597                   }
6598                } else {
6599                   for ( int tx = 0; tx < num_ut; tx++ ){
6600                      const double f_tt = fock->get( irrep_ut, nocc_ut + tx, nocc_ut + tx );
6601                      for ( int uy = tx+1; uy < num_ut; uy++ ){
6602                         const double f_uu = fock->get( irrep_ut, nocc_ut + uy, nocc_ut + uy );
6603                         FBB_triplet[ 0 ][ shift + tx + ( uy * ( uy - 1 ) ) / 2 + SIZE * ( tx + ( uy * ( uy - 1 ) ) / 2 ) ] += 6 * ( f_dot_1dm - f_tt - f_uu );
6604                      }
6605                   }
6606                }
6607 
6608                // SBB triplet: - 3 delta_uy Gamma_{tx}
6609                // FBB triplet: - 3 delta_uy ( f_dot_2dm[ tx ] - f_yu Gamma_{tx} )
6610                if ( OVLP ){
6611                   for ( int uy = 0; uy < num_ut; uy++ ){
6612                      for ( int t = 0; t < uy; t++ ){ // 0 <= t < uy < num_ut
6613                         for ( int x = 0; x < uy; x++ ){ // 0 <= x < uy < num_ut
6614                            const double gamma_tx = one_rdm[ d_ut + t + LAS * ( d_xy + x ) ];
6615                            SBB_triplet[ 0 ][ shift + x + ( uy * ( uy - 1 ) ) / 2 + SIZE * ( t + ( uy * ( uy - 1 ) ) / 2 ) ] -= 3 * gamma_tx;
6616                         }
6617                      }
6618                   }
6619                } else {
6620                   for ( int uy = 0; uy < num_ut; uy++ ){
6621                      const double f_uu = fock->get( irrep_ut, nocc_ut + uy, nocc_ut + uy );
6622                      for ( int t = 0; t < uy; t++ ){ // 0 <= t < uy < num_ut
6623                         for ( int x = 0; x < uy; x++ ){ // 0 <= x < uy < num_ut
6624                            const double val_tx = ( f_dot_2dm[ d_ut + t + LAS * ( d_ut + x ) ]
6625                                             - f_uu * one_rdm[ d_ut + t + LAS * ( d_ut + x ) ] );
6626                            FBB_triplet[ 0 ][ shift + x + ( uy * ( uy - 1 ) ) / 2 + SIZE * ( t + ( uy * ( uy - 1 ) ) / 2 ) ] -= 3 * val_tx;
6627                         }
6628                      }
6629                   }
6630                }
6631 
6632                // SBB triplet: - 3 delta_tx Gamma_{uy}
6633                // FBB triplet: - 3 delta_tx ( f_dot_2dm[ uy ] - f_tt Gamma_{uy} )
6634                if ( OVLP ){
6635                   for ( int tx = 0; tx < num_ut; tx++ ){
6636                      for ( int u = tx+1; u < num_ut; u++ ){ // 0 <= tx < u < num_ut
6637                         for ( int y = tx+1; y < num_ut; y++ ){ // 0 <= tx < y < num_xy = num_ut
6638                            const double gamma_uy = one_rdm[ d_ut + u + LAS * ( d_xy + y ) ];
6639                            SBB_triplet[ 0 ][ shift + tx + ( y * ( y - 1 ) ) / 2 + SIZE * ( tx + ( u * ( u - 1 ) ) / 2 ) ] -= 3 * gamma_uy;
6640                         }
6641                      }
6642                   }
6643                } else {
6644                   for ( int tx = 0; tx < num_ut; tx++ ){
6645                      const double f_tt = fock->get( irrep_ut, nocc_ut + tx, nocc_ut + tx );
6646                      for ( int u = tx+1; u < num_ut; u++ ){ // 0 <= tx < u < num_ut
6647                         for ( int y = tx+1; y < num_ut; y++ ){ // 0 <= tx < y < num_xy = num_ut
6648                            const double val_uy = ( f_dot_2dm[ d_ut + u + LAS * ( d_ut + y ) ]
6649                                             - f_tt * one_rdm[ d_ut + u + LAS * ( d_ut + y ) ] );
6650                            FBB_triplet[ 0 ][ shift + tx + ( y * ( y - 1 ) ) / 2 + SIZE * ( tx + ( u * ( u - 1 ) ) / 2 ) ] -= 3 * val_uy;
6651                         }
6652                      }
6653                   }
6654                }
6655 
6656                // SBB triplet: + 3 delta_ux Gamma_{ty}
6657                // FBB triplet: + 3 delta_ux ( f_dot_2dm[ ty ] - f_uu Gamma_{ty} )
6658                if ( OVLP ){
6659                   for ( int ux = 0; ux < num_ut; ux++ ){
6660                      for ( int t = 0; t < ux; t++ ){ // 0 <= t < ux < num_ut
6661                         for ( int y = ux+1; y < num_ut; y++ ){ // 0 <= ux < y < num_xy = num_ut
6662                            const double gamma_ty = one_rdm[ d_ut + t + LAS * ( d_xy + y ) ];
6663                            SBB_triplet[ 0 ][ shift + ux + ( y * ( y - 1 ) ) / 2 + SIZE * ( t + ( ux * ( ux - 1 ) ) / 2 ) ] += 3 * gamma_ty;
6664                         }
6665                      }
6666                   }
6667                } else {
6668                   for ( int ux = 0; ux < num_ut; ux++ ){
6669                      const double f_uu = fock->get( irrep_ut, nocc_ut + ux, nocc_ut + ux );
6670                      for ( int t = 0; t < ux; t++ ){ // 0 <= t < ux < num_ut
6671                         for ( int y = ux+1; y < num_ut; y++ ){ // 0 <= ux < y < num_xy = num_ut
6672                            const double val_ty = ( f_dot_2dm[ d_ut + t + LAS * ( d_ut + y ) ]
6673                                             - f_uu * one_rdm[ d_ut + t + LAS * ( d_ut + y ) ] );
6674                            FBB_triplet[ 0 ][ shift + ux + ( y * ( y - 1 ) ) / 2 + SIZE * ( t + ( ux * ( ux - 1 ) ) / 2 ) ] += 3 * val_ty;
6675                         }
6676                      }
6677                   }
6678                }
6679 
6680                // SBB triplet: + 3 delta_ty Gamma_{ux}
6681                // FBB triplet: + 3 delta_ty ( f_dot_2dm[ ux ] - f_tt Gamma_{ux} )
6682                if ( OVLP ){
6683                   for ( int ty = 0; ty < num_ut; ty++ ){
6684                      for ( int u = ty+1; u < num_ut; u++ ){ // 0 <= ty < u < num_ut
6685                         for ( int x = 0; x < ty; x++ ){ // 0 <= x < ty < num_ut
6686                            const double gamma_ux = one_rdm[ d_ut + u + LAS * ( d_xy + x ) ];
6687                            SBB_triplet[ 0 ][ shift + x + ( ty * ( ty - 1 ) ) / 2 + SIZE * ( ty + ( u * ( u - 1 ) ) / 2 ) ] += 3 * gamma_ux;
6688                         }
6689                      }
6690                   }
6691                } else {
6692                   for ( int ty = 0; ty < num_ut; ty++ ){
6693                      const double f_tt = fock->get( irrep_ut, nocc_ut + ty, nocc_ut + ty );
6694                      for ( int u = ty+1; u < num_ut; u++ ){ // 0 <= ty < u < num_ut
6695                         for ( int x = 0; x < ty; x++ ){ // 0 <= x < ty < num_ut
6696                            const double val_ux = ( f_dot_2dm[ d_ut + u + LAS * ( d_ut + x ) ]
6697                                             - f_tt * one_rdm[ d_ut + u + LAS * ( d_ut + x ) ] );
6698                            FBB_triplet[ 0 ][ shift + x + ( ty * ( ty - 1 ) ) / 2 + SIZE * ( ty + ( u * ( u - 1 ) ) / 2 ) ] += 3 * val_ux;
6699                         }
6700                      }
6701                   }
6702                }
6703             }
6704             jump_row += ( num_xy * ( num_xy - 1 ) ) / 2;
6705          }
6706          if (( OVLP == false ) && ( fabs( IPEA ) > 0.0 )){
6707             // B: E_ti E_uj | 0 >   --->   tu: excitation into
6708             // F: E_at E_bu | 0 >   --->   tu: excitation out of
6709             for ( int u = 0; u < num_ut; u++ ){
6710                const double gamma_uu = one_rdm[ ( d_ut + u ) * ( 1 + LAS ) ];
6711                for ( int t = 0; t < u; t++ ){ // 0 <= t < u < num_ut
6712                   const double gamma_tt = one_rdm[ ( d_ut + t ) * ( 1 + LAS ) ];
6713                   const double ipea_B_tu = 0.5 * IPEA * ( gamma_tt + gamma_uu );
6714                   const double ipea_F_tu = 0.5 * IPEA * ( 4.0 - gamma_tt - gamma_uu );
6715                   const int ptr = ( jump_col + t + ( u * ( u - 1 ) ) / 2 ) * ( 1 + SIZE );
6716                   FBB_triplet[ 0 ][ ptr ] += ipea_B_tu * SBB_triplet[ 0 ][ ptr ];
6717                   FFF_triplet[ 0 ][ ptr ] += ipea_F_tu * SFF_triplet[ 0 ][ ptr ];
6718                }
6719             }
6720          }
6721          jump_col += ( num_ut * ( num_ut - 1 ) ) / 2;
6722       }
6723    }
6724 
6725    for ( int irrep = 1; irrep < num_irreps; irrep++ ){ // Then do irrep == Iia x Ijb != 0 -->  It != Iu and Ix != Iy
6726 
6727       assert( size_B_triplet[ irrep ] == size_F_triplet[ irrep ] ); // At construction
6728       const int SIZE = size_B_triplet[ irrep ];
6729       if ( OVLP ){ SBB_triplet[ irrep ] = new double[ SIZE * SIZE ];
6730                    SFF_triplet[ irrep ] = new double[ SIZE * SIZE ]; }
6731       else {       FBB_triplet[ irrep ] = new double[ SIZE * SIZE ];
6732                    FFF_triplet[ irrep ] = new double[ SIZE * SIZE ]; }
6733 
6734       int jump_col = 0;
6735       for ( int irrep_t = 0; irrep_t < num_irreps; irrep_t++ ){
6736          const int irrep_u = Irreps::directProd( irrep, irrep_t );
6737          if ( irrep_t < irrep_u ){
6738             const int d_t    = indices->getDMRGcumulative( irrep_t );
6739             const int num_t  = indices->getNDMRG( irrep_t );
6740             const int nocc_t = indices->getNOCC( irrep_t );
6741             const int d_u    = indices->getDMRGcumulative( irrep_u );
6742             const int num_u  = indices->getNDMRG( irrep_u );
6743             const int nocc_u = indices->getNOCC( irrep_u );
6744             int jump_row = 0;
6745             for ( int irrep_x = 0; irrep_x < num_irreps; irrep_x++ ){
6746                const int irrep_y = Irreps::directProd( irrep, irrep_x );
6747                if ( irrep_x < irrep_y ){
6748                   const int d_x    = indices->getDMRGcumulative( irrep_x );
6749                   const int num_x  = indices->getNDMRG( irrep_x );
6750                   const int nocc_x = indices->getNOCC( irrep_x );
6751                   const int d_y    = indices->getDMRGcumulative( irrep_y );
6752                   const int num_y  = indices->getNDMRG( irrep_y );
6753                   const int nocc_y = indices->getNOCC( irrep_y );
6754                   const int shift  = jump_row + SIZE * jump_col;
6755 
6756                   // SBB triplet: + Gamma_{utyx} - Gamma_{utxy}
6757                   // SFF triplet: + Gamma_{yxut} - Gamma_{yxtu}
6758                   // FBB triplet: + f_dot_3dm[ utyx ] - f_dot_3dm[ tuyx ] + ( f_tt + f_uu + f_xx + f_yy ) SBB_triplet[ xytu ]
6759                   // FFF triplet: + f_dot_3dm[ yxut ] - f_dot_3dm[ yxtu ]
6760                   if ( OVLP ){
6761                      for ( int t = 0; t < num_t; t++ ){
6762                         for ( int u = 0; u < num_u; u++ ){
6763                            for ( int x = 0; x < num_x; x++ ){
6764                               for ( int y = 0; y < num_y; y++ ){
6765                                  const int ptr = shift + x + num_x * y + SIZE * ( t + num_t * u );
6766                                  const double value = ( two_rdm[ d_x + x + LAS * ( d_y + y + LAS * ( d_t + t + LAS * ( d_u + u ))) ]
6767                                                       - two_rdm[ d_x + x + LAS * ( d_y + y + LAS * ( d_u + u + LAS * ( d_t + t ))) ] );
6768                                  SBB_triplet[ irrep ][ ptr ] = value;
6769                                  SFF_triplet[ irrep ][ ptr ] = value;
6770                               }
6771                            }
6772                         }
6773                      }
6774                   } else {
6775                      for ( int t = 0; t < num_t; t++ ){
6776                         const double f_tt = fock->get( irrep_t, nocc_t + t, nocc_t + t );
6777                         for ( int u = 0; u < num_u; u++ ){
6778                            const double f_uu = fock->get( irrep_u, nocc_u + u, nocc_u + u );
6779                            for ( int x = 0; x < num_x; x++ ){
6780                               const double f_xx = fock->get( irrep_x, nocc_x + x, nocc_x + x );
6781                               for ( int y = 0; y < num_y; y++ ){
6782                                  const double f_yy = fock->get( irrep_y, nocc_y + y, nocc_y + y );
6783                                  const int ptr = shift + x + num_x * y + SIZE * ( t + num_t * u );
6784                                  const double fdotdiff = ( f_dot_3dm[ d_x + x + LAS * ( d_y + y + LAS * ( d_t + t + LAS * ( d_u + u ))) ]
6785                                                          - f_dot_3dm[ d_x + x + LAS * ( d_y + y + LAS * ( d_u + u + LAS * ( d_t + t ))) ] );
6786                                  FBB_triplet[ irrep ][ ptr ] = fdotdiff + ( f_tt + f_uu + f_xx + f_yy ) * SBB_triplet[ irrep ][ ptr ];
6787                                  FFF_triplet[ irrep ][ ptr ] = fdotdiff;
6788                               }
6789                            }
6790                         }
6791                      }
6792                   }
6793 
6794                   if (( irrep_u == irrep_y ) && ( irrep_t == irrep_x )){ // --> num_t == num_x   and   num_u == num_y
6795 
6796                      // SBB triplet: + 6 delta_uy delta_tx
6797                      // FBB triplet: + 6 delta_uy delta_tx ( f_dot_1dm - f_tt - f_uu )
6798                      if ( OVLP ){
6799                         for ( int xt = 0; xt < num_x; xt++ ){
6800                            for ( int yu = 0; yu < num_y; yu++ ){
6801                               SBB_triplet[ irrep ][ shift + xt + num_x * yu + SIZE * ( xt + num_x * yu ) ] += 6.0;
6802                            }
6803                         }
6804                      } else {
6805                         for ( int xt = 0; xt < num_x; xt++ ){
6806                            const double f_tt = fock->get( irrep_x, nocc_x + xt, nocc_x + xt );
6807                            for ( int yu = 0; yu < num_y; yu++ ){
6808                               const double f_uu = fock->get( irrep_y, nocc_y + yu, nocc_y + yu );
6809                               FBB_triplet[ irrep ][ shift + xt + num_x * yu + SIZE * ( xt + num_x * yu ) ] += 6 * ( f_dot_1dm - f_tt - f_uu );
6810                            }
6811                         }
6812                      }
6813 
6814                      // SBB triplet: - 3 delta_tx Gamma_{uy}
6815                      // FBB triplet: - 3 delta_tx ( f_dot_2dm[ uy ] - f_tt Gamma_{uy} )
6816                      if ( OVLP ){
6817                         for ( int xt = 0; xt < num_x; xt++ ){
6818                            for ( int u = 0; u < num_y; u++ ){
6819                               for ( int y = 0; y < num_y; y++ ){
6820                                  const double gamma_uy = one_rdm[ d_u + u + LAS * ( d_y + y ) ];
6821                                  SBB_triplet[ irrep ][ shift + xt + num_x * y + SIZE * ( xt + num_x * u ) ] -= 3 * gamma_uy;
6822                               }
6823                            }
6824                         }
6825                      } else {
6826                         for ( int xt = 0; xt < num_x; xt++ ){
6827                            const double f_tt = fock->get( irrep_x, nocc_x + xt, nocc_x + xt );
6828                            for ( int u = 0; u < num_y; u++ ){
6829                               for ( int y = 0; y < num_y; y++ ){
6830                                  const double val_uy = ( f_dot_2dm[ d_u + u + LAS * ( d_u + y ) ]
6831                                                   - f_tt * one_rdm[ d_u + u + LAS * ( d_u + y ) ] );
6832                                  FBB_triplet[ irrep ][ shift + xt + num_x * y + SIZE * ( xt + num_x * u ) ] -= 3 * val_uy;
6833                               }
6834                            }
6835                         }
6836                      }
6837 
6838                      // SBB triplet: - 3 delta_uy Gamma_{tx}
6839                      // FBB triplet: - 3 delta_uy ( f_dot_2dm[ tx ] - f_uu Gamma_{tx} )
6840                      if ( OVLP ){
6841                         for ( int yu = 0; yu < num_y; yu++ ){
6842                            for ( int t = 0; t < num_x; t++ ){
6843                               for ( int x = 0; x < num_x; x++ ){
6844                                  const double gamma_tx = one_rdm[ d_t + t + LAS * ( d_x + x ) ];
6845                                  SBB_triplet[ irrep ][ shift + x + num_x * yu + SIZE * ( t + num_x * yu ) ] -= 3 * gamma_tx;
6846                               }
6847                            }
6848                         }
6849                      } else {
6850                         for ( int yu = 0; yu < num_y; yu++ ){
6851                            const double f_uu = fock->get( irrep_y, nocc_y + yu, nocc_y + yu );
6852                            for ( int t = 0; t < num_x; t++ ){
6853                               for ( int x = 0; x < num_x; x++ ){
6854                                  const double val_tx = ( f_dot_2dm[ d_t + t + LAS * ( d_t + x ) ]
6855                                                   - f_uu * one_rdm[ d_t + t + LAS * ( d_t + x ) ] );
6856                                  FBB_triplet[ irrep ][ shift + x + num_x * yu + SIZE * ( t + num_x * yu ) ] -= 3 * val_tx;
6857                               }
6858                            }
6859                         }
6860                      }
6861                   }
6862                   jump_row += num_x * num_y;
6863                }
6864             }
6865             if (( OVLP == false ) && ( fabs( IPEA ) > 0.0 )){
6866                // B: E_ti E_uj | 0 >   --->   tu: excitation into
6867                // F: E_at E_bu | 0 >   --->   tu: excitation out of
6868                for ( int u = 0; u < num_u; u++ ){
6869                   const double gamma_uu = one_rdm[ ( d_u + u ) * ( 1 + LAS ) ];
6870                   for ( int t = 0; t < num_t; t++ ){
6871                      const double gamma_tt = one_rdm[ ( d_t + t ) * ( 1 + LAS ) ];
6872                      const double ipea_B_tu = 0.5 * IPEA * ( gamma_tt + gamma_uu );
6873                      const double ipea_F_tu = 0.5 * IPEA * ( 4.0 - gamma_tt - gamma_uu );
6874                      const int ptr = ( jump_col + t + num_t * u ) * ( 1 + SIZE );
6875                      FBB_triplet[ irrep ][ ptr ] += ipea_B_tu * SBB_triplet[ irrep ][ ptr ];
6876                      FFF_triplet[ irrep ][ ptr ] += ipea_F_tu * SFF_triplet[ irrep ][ ptr ];
6877                   }
6878                }
6879             }
6880             jump_col += num_t * num_u;
6881          }
6882       }
6883    }
6884 
6885 }
6886 
make_EE_GG(const bool OVLP,const double IPEA)6887 void CheMPS2::CASPT2::make_EE_GG( const bool OVLP, const double IPEA ){
6888 
6889    /*
6890       | SE_tiaj > = ( E_ti E_aj + E_tj E_ai ) / sqrt( 1 + delta_ij ) | 0 >  with  i <= j
6891       | TE_tiaj > = ( E_ti E_aj - E_tj E_ai ) / sqrt( 1 + delta_ij ) | 0 >  with  i <  j
6892 
6893          < SE_ukbl | 1 | SE_tiaj > = 2 delta_ab delta_ik delta_jl ( SEE[ It ][ ut ] )
6894          < SE_ukbl | F | SE_tiaj > = 2 delta_ab delta_ik delta_jl ( FEE[ It ][ ut ] + ( 2 sum_k f_kk + f_aa - f_ii - f_jj ) SEE[ It ][ ut ] )
6895          < TE_ukbl | 1 | TE_tiaj > = 6 delta_ab delta_ik delta_jl ( SEE[ It ][ ut ] )
6896          < TE_ukbl | F | TE_tiaj > = 6 delta_ab delta_ik delta_jl ( FEE[ It ][ ut ] + ( 2 sum_k f_kk + f_aa - f_ii - f_jj ) SEE[ It ][ ut ] )
6897 
6898             SEE[ It ][ ut ] = ( + 2 delta_tu
6899                                 - Gamma_tu
6900                               )
6901 
6902             FEE[ It ][ ut ] = ( + 2 * delta_ut * ( f_dot_1dm + f_tt )
6903                                 - f_dot_2dm[ It ][ tu ] - ( f_tt + f_uu ) Gamma_ut
6904                               )
6905 
6906       | SG_aibt > = ( E_ai E_bt + E_bi E_at ) / sqrt( 1 + delta_ab ) | 0 >  with  a <= b
6907       | TG_aibt > = ( E_ai E_bt - E_bi E_at ) / sqrt( 1 + delta_ab ) | 0 >  with  a <  b
6908 
6909          < SG_cjdu | 1 | SG_aibt > = 2 delta_ij delta_ac delta_bd ( SGG[ It ][ ut ] )
6910          < SG_cjdu | F | SG_aibt > = 2 delta_ij delta_ac delta_bd ( FGG[ It ][ ut ] + ( 2 sum_k f_kk + f_aa + f_bb - f_ii ) SGG[ It ][ ut ] )
6911          < TG_cjdu | 1 | TG_aibt > = 6 delta_ij delta_ac delta_bd ( SGG[ It ][ ut ] )
6912          < TG_cjdu | F | TG_aibt > = 6 delta_ij delta_ac delta_bd ( FGG[ It ][ ut ] + ( 2 sum_k f_kk + f_aa + f_bb - f_ii ) SGG[ It ][ ut ] )
6913 
6914             SGG[ It ][ ut ] = ( + Gamma_ut
6915                               )
6916 
6917             FGG[ It ][ ut ] = ( + f_dot_2dm[ It ][ ut ]
6918                               )
6919    */
6920 
6921    if ( OVLP ){ SEE = new double*[ num_irreps ];
6922                 SGG = new double*[ num_irreps ]; }
6923    else {       FEE = new double*[ num_irreps ];
6924                 FGG = new double*[ num_irreps ]; }
6925 
6926    const int LAS = indices->getDMRGcumulative( num_irreps );
6927 
6928    for ( int irrep_ut = 0; irrep_ut < num_irreps; irrep_ut++ ){
6929       const int d_ut  = indices->getDMRGcumulative( irrep_ut );
6930       const int SIZE  = indices->getNDMRG( irrep_ut );
6931       const int NOCC  = indices->getNOCC( irrep_ut );
6932       if ( OVLP ){
6933          SEE[ irrep_ut ] = new double[ SIZE * SIZE ];
6934          SGG[ irrep_ut ] = new double[ SIZE * SIZE ];
6935          for ( int t = 0; t < SIZE; t++ ){
6936             for ( int u = 0; u < SIZE; u++ ){
6937                const double gamma_ut = one_rdm[ d_ut + u + LAS * ( d_ut + t ) ];
6938                SEE[ irrep_ut ][ u + SIZE * t ] = - gamma_ut;
6939                SGG[ irrep_ut ][ u + SIZE * t ] =   gamma_ut;
6940             }
6941             SEE[ irrep_ut ][ t + SIZE * t ] += 2.0;
6942          }
6943       } else {
6944          FEE[ irrep_ut ] = new double[ SIZE * SIZE ];
6945          FGG[ irrep_ut ] = new double[ SIZE * SIZE ];
6946          for ( int t = 0; t < SIZE; t++ ){
6947             const double f_tt = fock->get( irrep_ut, NOCC + t, NOCC + t );
6948             for ( int u = 0; u < SIZE; u++ ){
6949                const double f_uu = fock->get( irrep_ut, NOCC + u, NOCC + u );
6950                const double gamma_ut =   one_rdm[ d_ut + u + LAS * ( d_ut + t ) ];
6951                const double fdot2_ut = f_dot_2dm[ d_ut + u + LAS * ( d_ut + t ) ];
6952                FEE[ irrep_ut ][ u + SIZE * t ] = - fdot2_ut - ( f_tt + f_uu ) * gamma_ut;
6953                FGG[ irrep_ut ][ u + SIZE * t ] =   fdot2_ut;
6954             }
6955             FEE[ irrep_ut ][ t + SIZE * t ] += 2 * ( f_dot_1dm + f_tt );
6956          }
6957          if ( fabs( IPEA ) > 0.0 ){
6958             // E: E_ti E_aj | 0 >   --->   t: excitation into
6959             // G: E_ai E_bt | 0 >   --->   t: excitation out of
6960             for ( int t = 0; t < SIZE; t++ ){
6961                const double gamma_tt = one_rdm[ ( d_ut + t ) * ( 1 + LAS ) ];
6962                const double ipea_E_t = 0.5 * IPEA * ( gamma_tt );
6963                const double ipea_G_t = 0.5 * IPEA * ( 2.0 - gamma_tt );
6964                const int ptr = t * ( 1 + SIZE );
6965                FEE[ irrep_ut ][ ptr ] += ipea_E_t * SEE[ irrep_ut ][ ptr ];
6966                FGG[ irrep_ut ][ ptr ] += ipea_G_t * SGG[ irrep_ut ][ ptr ];
6967             }
6968          }
6969       }
6970    }
6971 
6972 }
6973 
6974 
6975