CheMPS2
Sobject.cpp
1 /*
2  CheMPS2: a spin-adapted implementation of DMRG for ab initio quantum chemistry
3  Copyright (C) 2013-2016 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 <math.h>
21 #include <stdlib.h>
22 #include <algorithm>
23 #include <assert.h>
24 
25 #include "Sobject.h"
26 #include "TensorT.h"
27 #include "SyBookkeeper.h"
28 #include "Lapack.h"
29 #include "MPIchemps2.h"
30 #include "Wigner.h"
31 #include "Special.h"
32 
33 using std::min;
34 using std::max;
35 
36 CheMPS2::Sobject::Sobject( const int index, SyBookkeeper * denBK ){
37 
38  this->index = index;
39  this->denBK = denBK;
40  Ilocal1 = denBK->gIrrep( index );
41  Ilocal2 = denBK->gIrrep( index + 1 );
42 
43  nKappa = 0;
44 
45  for ( int NL = denBK->gNmin( index ); NL <= denBK->gNmax( index ); NL++ ){
46  for ( int TwoSL = denBK->gTwoSmin( index, NL ); TwoSL <= denBK->gTwoSmax( index, NL ); TwoSL += 2 ){
47  for ( int IL = 0; IL < denBK->getNumberOfIrreps(); IL++ ){
48  const int dimL = denBK->gCurrentDim( index, NL, TwoSL, IL );
49  if ( dimL > 0 ){
50  for ( int N1 = 0; N1 <= 2; N1++ ){
51  for ( int N2 = 0; N2 <= 2; N2++ ){
52  const int NR = NL + N1 + N2;
53  const int IM = (( N1 == 1 ) ? Irreps::directProd( IL, Ilocal1 ) : IL );
54  const int IR = (( N2 == 1 ) ? Irreps::directProd( IM, Ilocal2 ) : IM );
55  const int TwoJmin = ( N1 + N2 ) % 2;
56  const int TwoJmax = ((( N1 == 1 ) && ( N2 == 1 )) ? 2 : TwoJmin );
57  for ( int TwoJ = TwoJmin; TwoJ <= TwoJmax; TwoJ += 2 ){
58  for ( int TwoSR = TwoSL - TwoJ; TwoSR <= TwoSL + TwoJ; TwoSR += 2 ){
59  if ( TwoSR >= 0 ){
60  const int dimR = denBK->gCurrentDim( index + 2, NR, TwoSR, IR );
61  if ( dimR > 0 ){
62  nKappa++;
63  }
64  }
65  }
66  }
67  }
68  }
69  }
70  }
71  }
72  }
73 
74  sectorNL = new int[ nKappa ];
75  sectorTwoSL = new int[ nKappa ];
76  sectorIL = new int[ nKappa ];
77  sectorN1 = new int[ nKappa ];
78  sectorN2 = new int[ nKappa ];
79  sectorTwoJ = new int[ nKappa ];
80  sectorNR = new int[ nKappa ];
81  sectorTwoSR = new int[ nKappa ];
82  sectorIR = new int[ nKappa ];
83  kappa2index = new int[ nKappa + 1 ];
84  kappa2index[ 0 ] = 0;
85 
86  nKappa = 0;
87 
88  for ( int NL = denBK->gNmin( index ); NL <= denBK->gNmax( index ); NL++ ){
89  for ( int TwoSL = denBK->gTwoSmin( index, NL ); TwoSL <= denBK->gTwoSmax( index, NL ); TwoSL += 2 ){
90  for ( int IL = 0; IL < denBK->getNumberOfIrreps(); IL++ ){
91  const int dimL = denBK->gCurrentDim( index, NL, TwoSL, IL );
92  if ( dimL > 0 ){
93  for ( int N1 = 0; N1 <= 2; N1++ ){
94  for ( int N2 = 0; N2 <= 2; N2++ ){
95  const int NR = NL + N1 + N2;
96  const int IM = (( N1 == 1 ) ? Irreps::directProd( IL, Ilocal1 ) : IL );
97  const int IR = (( N2 == 1 ) ? Irreps::directProd( IM, Ilocal2 ) : IM );
98  const int TwoJmin = ( N1 + N2 ) % 2;
99  const int TwoJmax = ((( N1 == 1 ) && ( N2 == 1 )) ? 2 : TwoJmin );
100  for ( int TwoJ = TwoJmin; TwoJ <= TwoJmax; TwoJ += 2 ){
101  for ( int TwoSR = TwoSL - TwoJ; TwoSR <= TwoSL + TwoJ; TwoSR += 2 ){
102  if ( TwoSR >= 0 ){
103  const int dimR = denBK->gCurrentDim( index + 2, NR, TwoSR, IR );
104  if ( dimR > 0 ){
105  sectorNL [ nKappa ] = NL;
106  sectorTwoSL[ nKappa ] = TwoSL;
107  sectorIL [ nKappa ] = IL;
108  sectorN1 [ nKappa ] = N1;
109  sectorN2 [ nKappa ] = N2;
110  sectorTwoJ [ nKappa ] = TwoJ;
111  sectorNR [ nKappa ] = NR;
112  sectorTwoSR[ nKappa ] = TwoSR;
113  sectorIR [ nKappa ] = IR;
114  nKappa++;
115  kappa2index[ nKappa ] = kappa2index[ nKappa - 1 ] + dimL * dimR;
116  }
117  }
118  }
119  }
120  }
121  }
122  }
123  }
124  }
125  }
126 
127  storage = new double[ kappa2index[ nKappa ] ];
128 
129  reorder = new int[ nKappa ];
130  for ( int cnt = 0; cnt < nKappa; cnt++ ){ reorder[ cnt ] = cnt; }
131  bool sorted = false;
132  while ( sorted == false ){ //Bubble sort so that blocksize(reorder[i]) >= blocksize(reorder[i+1]), with blocksize(k) = kappa2index[k+1]-kappa2index[k]
133  sorted = true;
134  for ( int cnt = 0; cnt < nKappa - 1; cnt++ ){
135  const int index1 = reorder[ cnt ];
136  const int index2 = reorder[ cnt + 1 ];
137  const int size1 = kappa2index[ index1 + 1 ] - kappa2index[ index1 ];
138  const int size2 = kappa2index[ index2 + 1 ] - kappa2index[ index2 ];
139  if ( size1 < size2 ){
140  sorted = false;
141  reorder[ cnt ] = index2;
142  reorder[ cnt + 1 ] = index1;
143  }
144  }
145  }
146 
147 }
148 
150 
151  delete [] sectorNL;
152  delete [] sectorTwoSL;
153  delete [] sectorIL;
154  delete [] sectorN1;
155  delete [] sectorN2;
156  delete [] sectorTwoJ;
157  delete [] sectorNR;
158  delete [] sectorTwoSR;
159  delete [] sectorIR;
160  delete [] kappa2index;
161  delete [] storage;
162  delete [] reorder;
163 
164 }
165 
166 int CheMPS2::Sobject::gNKappa() const { return nKappa; }
167 
168 double * CheMPS2::Sobject::gStorage() { return storage; }
169 
170 int CheMPS2::Sobject::gReorder( const int ikappa ) const{ return reorder[ ikappa ]; }
171 
172 int CheMPS2::Sobject::gKappa( const int NL, const int TwoSL, const int IL, const int N1, const int N2, const int TwoJ, const int NR, const int TwoSR, const int IR ) const{
173 
174  for ( int ikappa = 0; ikappa < nKappa; ikappa++ ){
175  if (( sectorNL [ ikappa ] == NL ) &&
176  ( sectorTwoSL[ ikappa ] == TwoSL ) &&
177  ( sectorIL [ ikappa ] == IL ) &&
178  ( sectorN1 [ ikappa ] == N1 ) &&
179  ( sectorN2 [ ikappa ] == N2 ) &&
180  ( sectorTwoJ [ ikappa ] == TwoJ ) &&
181  ( sectorNR [ ikappa ] == NR ) &&
182  ( sectorTwoSR[ ikappa ] == TwoSR ) &&
183  ( sectorIR [ ikappa ] == IR )){ return ikappa; }
184  }
185 
186  return -1;
187 
188 }
189 
190 int CheMPS2::Sobject::gKappa2index( const int kappa ) const{ return kappa2index[ kappa ]; }
191 
192 double * CheMPS2::Sobject::gStorage( const int NL, const int TwoSL, const int IL, const int N1, const int N2, const int TwoJ, const int NR, const int TwoSR, const int IR ){
193 
194  const int kappa = gKappa( NL, TwoSL, IL, N1, N2, TwoJ, NR, TwoSR, IR );
195  if ( kappa == -1 ){ return NULL; }
196  return storage + kappa2index[ kappa ];
197 
198 }
199 
200 int CheMPS2::Sobject::gIndex() const { return index; }
201 
202 int CheMPS2::Sobject::gNL ( const int ikappa ) const{ return sectorNL [ ikappa ]; }
203 int CheMPS2::Sobject::gTwoSL( const int ikappa ) const{ return sectorTwoSL[ ikappa ]; }
204 int CheMPS2::Sobject::gIL ( const int ikappa ) const{ return sectorIL [ ikappa ]; }
205 int CheMPS2::Sobject::gN1 ( const int ikappa ) const{ return sectorN1 [ ikappa ]; }
206 int CheMPS2::Sobject::gN2 ( const int ikappa ) const{ return sectorN2 [ ikappa ]; }
207 int CheMPS2::Sobject::gTwoJ ( const int ikappa ) const{ return sectorTwoJ [ ikappa ]; }
208 int CheMPS2::Sobject::gNR ( const int ikappa ) const{ return sectorNR [ ikappa ]; }
209 int CheMPS2::Sobject::gTwoSR( const int ikappa ) const{ return sectorTwoSR[ ikappa ]; }
210 int CheMPS2::Sobject::gIR ( const int ikappa ) const{ return sectorIR [ ikappa ]; }
211 
212 void CheMPS2::Sobject::Join( TensorT * Tleft, TensorT * Tright ){
213 
214  #pragma omp parallel for schedule(dynamic)
215  for ( int ikappa = 0; ikappa < nKappa; ikappa++ ){
216 
217  const int NL = sectorNL [ ikappa ];
218  const int TwoSL = sectorTwoSL[ ikappa ];
219  const int IL = sectorIL [ ikappa ];
220 
221  const int NR = sectorNR [ ikappa ];
222  const int TwoSR = sectorTwoSR[ ikappa ];
223  const int IR = sectorIR [ ikappa ];
224 
225  const int TwoJ = sectorTwoJ [ ikappa ];
226  const int N1 = sectorN1 [ ikappa ];
227  const int N2 = sectorN2 [ ikappa ];
228  const int TwoS1 = (( N1 == 1 ) ? 1 : 0 );
229  const int TwoS2 = (( N2 == 1 ) ? 1 : 0 );
230  const int fase = Special::phase( TwoSL + TwoSR + TwoS1 + TwoS2 );
231 
232  // Clear block
233  int dimL = denBK->gCurrentDim( index, NL, TwoSL, IL ); // dimL > 0, checked at creation
234  int dimR = denBK->gCurrentDim( index + 2, NR, TwoSR, IR ); // dimR > 0, checked at creation
235  double * block_s = storage + kappa2index[ ikappa ];
236  for ( int cnt = 0; cnt < dimL * dimR; cnt++ ){ block_s[ cnt ] = 0.0; }
237 
238  // Central symmetry sectors
239  const int NM = NL + N1;
240  const int IM = (( TwoS1 == 1 ) ? Irreps::directProd( IL, Ilocal1 ) : IL );
241  const int TwoJMlower = max( abs( TwoSL - TwoS1 ), abs( TwoSR - TwoS2 ) );
242  const int TwoJMupper = min( ( TwoSL + TwoS1 ), ( TwoSR + TwoS2 ) );
243  for ( int TwoJM = TwoJMlower; TwoJM <= TwoJMupper; TwoJM += 2 ){
244  int dimM = denBK->gCurrentDim( index + 1, NM, TwoJM, IM );
245  if ( dimM > 0 ){
246  double * block_left = Tleft ->gStorage( NL, TwoSL, IL, NM, TwoJM, IM );
247  double * block_right = Tright->gStorage( NM, TwoJM, IM, NR, TwoSR, IR );
248  double prefactor = fase
249  * sqrt( 1.0 * ( TwoJ + 1 ) * ( TwoJM + 1 ) )
250  * Wigner::wigner6j( TwoSL, TwoSR, TwoJ, TwoS2, TwoS1, TwoJM );
251  char notrans = 'N';
252  double add = 1.0;
253  dgemm_( &notrans, &notrans, &dimL, &dimR, &dimM, &prefactor, block_left, &dimL, block_right, &dimM, &add, block_s, &dimL );
254  }
255  }
256  }
257 
258 }
259 
260 double CheMPS2::Sobject::Split( TensorT * Tleft, TensorT * Tright, const int virtualdimensionD, const bool movingright, const bool change ){
261 
262  #ifdef CHEMPS2_MPI_COMPILATION
263  const bool am_i_master = ( MPIchemps2::mpi_rank() == MPI_CHEMPS2_MASTER );
264  #endif
265 
266  // Get the number of central sectors
267  int nCenterSectors = 0;
268  for ( int NM = denBK->gNmin( index + 1 ); NM <= denBK->gNmax( index + 1 ); NM++ ){
269  for ( int TwoSM = denBK->gTwoSmin( index + 1, NM ); TwoSM <= denBK->gTwoSmax( index + 1, NM ); TwoSM += 2 ){
270  for ( int IM = 0; IM < denBK->getNumberOfIrreps(); IM++ ){
271  const int dimM = denBK->gFCIdim( index + 1, NM, TwoSM, IM ); //FCIdim !! Whether possible hence.
272  if ( dimM > 0 ){
273  nCenterSectors++;
274  }
275  }
276  }
277  }
278 
279  // Get the labels of the central sectors
280  int * SplitSectNM = new int[ nCenterSectors ];
281  int * SplitSectTwoJM = new int[ nCenterSectors ];
282  int * SplitSectIM = new int[ nCenterSectors ];
283  nCenterSectors = 0;
284  for ( int NM = denBK->gNmin( index + 1 ); NM <= denBK->gNmax( index + 1 ); NM++ ){
285  for ( int TwoSM = denBK->gTwoSmin( index + 1, NM ); TwoSM <= denBK->gTwoSmax( index + 1, NM ); TwoSM += 2 ){
286  for ( int IM = 0; IM < denBK->getNumberOfIrreps(); IM++ ){
287  const int dimM = denBK->gFCIdim( index + 1, NM, TwoSM, IM ); //FCIdim !! Whether possible hence.
288  if ( dimM > 0 ){
289  SplitSectNM [ nCenterSectors ] = NM;
290  SplitSectTwoJM[ nCenterSectors ] = TwoSM;
291  SplitSectIM [ nCenterSectors ] = IM;
292  nCenterSectors++;
293  }
294  }
295  }
296  }
297 
298  // Only MPI_CHEMPS2_MASTER performs SVD --> Allocate memory
299  double ** Lambdas = NULL;
300  double ** Us = NULL;
301  double ** VTs = NULL;
302  int * CenterDims = NULL;
303  int * DimLtotal = NULL;
304  int * DimRtotal = NULL;
305 
306  #ifdef CHEMPS2_MPI_COMPILATION
307  if ( am_i_master ){
308  #endif
309 
310  Lambdas = new double*[ nCenterSectors ];
311  Us = new double*[ nCenterSectors ];
312  VTs = new double*[ nCenterSectors ];
313  CenterDims = new int[ nCenterSectors ];
314  DimLtotal = new int[ nCenterSectors ];
315  DimRtotal = new int[ nCenterSectors ];
316 
317  //PARALLEL
318  #pragma omp parallel for schedule(dynamic)
319  for ( int iCenter = 0; iCenter < nCenterSectors; iCenter++ ){
320 
321  //Determine left and right dimensions contributing to the center block iCenter
322  DimLtotal[ iCenter ] = 0;
323  for ( int NL = SplitSectNM[ iCenter ] - 2; NL <= SplitSectNM[ iCenter ]; NL++ ){
324  const int TwoS1 = (( NL + 1 == SplitSectNM[ iCenter ] ) ? 1 : 0 );
325  for ( int TwoSL = SplitSectTwoJM[ iCenter ] - TwoS1; TwoSL <= SplitSectTwoJM[ iCenter ] + TwoS1; TwoSL += 2 ){
326  if ( TwoSL >= 0 ){
327  const int IL = (( TwoS1 == 1 ) ? Irreps::directProd( Ilocal1, SplitSectIM[ iCenter ] ) : SplitSectIM[ iCenter ] );
328  const int dimL = denBK->gCurrentDim( index, NL, TwoSL, IL );
329  if ( dimL > 0 ){
330  DimLtotal[ iCenter ] += dimL;
331  }
332  }
333  }
334  }
335  DimRtotal[ iCenter ] = 0;
336  for ( int NR = SplitSectNM[ iCenter ]; NR <= SplitSectNM[ iCenter ] + 2; NR++ ){
337  const int TwoS2 = (( NR == SplitSectNM[ iCenter ] + 1 ) ? 1 : 0 );
338  for ( int TwoSR = SplitSectTwoJM[ iCenter ] - TwoS2; TwoSR <= SplitSectTwoJM[ iCenter ] + TwoS2; TwoSR += 2 ){
339  if ( TwoSR >= 0 ){
340  const int IR = (( TwoS2 == 1 ) ? Irreps::directProd( Ilocal2, SplitSectIM[ iCenter ] ) : SplitSectIM[ iCenter ] );
341  const int dimR = denBK->gCurrentDim( index + 2, NR, TwoSR, IR );
342  if ( dimR > 0 ){
343  DimRtotal[ iCenter ] += dimR;
344  }
345  }
346  }
347  }
348  CenterDims[ iCenter ] = min( DimLtotal[ iCenter ], DimRtotal[ iCenter ] ); // CenterDims contains the min. amount
349 
350  //Allocate memory to copy the different parts of the S-object. Use prefactor sqrt((2jR+1)/(2jM+1) * (2jM+1) * (2j+1)) W6J (-1)^(jL+jR+s1+s2) and sum over j.
351  if ( CenterDims[ iCenter ] > 0 ){
352 
353  // Only if CenterDims[ iCenter ] exists should you allocate the following three arrays
354  Lambdas[ iCenter ] = new double[ CenterDims[ iCenter ] ];
355  Us[ iCenter ] = new double[ CenterDims[ iCenter ] * DimLtotal[ iCenter ] ];
356  VTs[ iCenter ] = new double[ CenterDims[ iCenter ] * DimRtotal[ iCenter ] ];
357 
358  const int memsize = DimLtotal[ iCenter ] * DimRtotal[ iCenter ];
359  double * mem = new double[ memsize ];
360  for ( int cnt = 0; cnt < memsize; cnt++ ){ mem[ cnt ] = 0.0; }
361 
362  int dimLtotal2 = 0;
363  for ( int NL = SplitSectNM[ iCenter ] - 2; NL <= SplitSectNM[ iCenter ]; NL++ ){
364  const int TwoS1 = (( NL + 1 == SplitSectNM[ iCenter ] ) ? 1 : 0 );
365  for ( int TwoSL = SplitSectTwoJM[ iCenter ] - TwoS1; TwoSL <= SplitSectTwoJM[ iCenter ] + TwoS1; TwoSL += 2 ){
366  if ( TwoSL >= 0 ){
367  const int IL = (( TwoS1 == 1 ) ? Irreps::directProd( Ilocal1, SplitSectIM[ iCenter ] ) : SplitSectIM[ iCenter ] );
368  const int dimL = denBK->gCurrentDim( index, NL, TwoSL, IL );
369  if ( dimL > 0 ){
370  int dimRtotal2 = 0;
371  for ( int NR = SplitSectNM[ iCenter ]; NR <= SplitSectNM[ iCenter ] + 2; NR++ ){
372  const int TwoS2 = (( NR == SplitSectNM[ iCenter ] + 1 ) ? 1 : 0 );
373  for ( int TwoSR = SplitSectTwoJM[ iCenter ] - TwoS2; TwoSR <= SplitSectTwoJM[ iCenter ] + TwoS2; TwoSR += 2 ){
374  if ( TwoSR >= 0 ){
375  const int IR = (( TwoS2 == 1 ) ? Irreps::directProd( Ilocal2, SplitSectIM[ iCenter ] ) : SplitSectIM[ iCenter ] );
376  const int dimR = denBK->gCurrentDim( index + 2, NR, TwoSR, IR );
377  if ( dimR > 0 ){
378  // Loop over contributing TwoJ's
379  const int fase = Special::phase( TwoSL + TwoSR + TwoS1 + TwoS2 );
380  const int TwoJmin = max( abs( TwoSR - TwoSL ), abs( TwoS2 - TwoS1 ) );
381  const int TwoJmax = min( TwoS1 + TwoS2, TwoSL + TwoSR );
382  for ( int TwoJ = TwoJmin; TwoJ <= TwoJmax; TwoJ += 2 ){
383  // Calc prefactor
384  const double prefactor = fase
385  * sqrt( 1.0 * ( TwoJ + 1 ) * ( TwoSR + 1 ) )
386  * Wigner::wigner6j( TwoSL, TwoSR, TwoJ, TwoS2, TwoS1, SplitSectTwoJM[ iCenter ] );
387 
388  // Add them to mem --> += because several TwoJ
389  double * Block = gStorage( NL, TwoSL, IL, SplitSectNM[ iCenter ] - NL, NR - SplitSectNM[ iCenter ], TwoJ, NR, TwoSR, IR );
390  for ( int l = 0; l < dimL; l++ ){
391  for ( int r = 0; r < dimR; r++ ){
392  mem[ dimLtotal2 + l + DimLtotal[ iCenter ] * ( dimRtotal2 + r ) ] += prefactor * Block[ l + dimL * r ];
393  }
394  }
395  }
396  dimRtotal2 += dimR;
397  }
398  }
399  }
400  }
401  dimLtotal2 += dimL;
402  }
403  }
404  }
405  }
406 
407  // Now mem contains sqrt((2jR+1)/(2jM+1)) * (TT)^{jM nM IM) --> SVD per central symmetry
408  char jobz = 'S'; // M x min(M,N) in U and min(M,N) x N in VT
409  int lwork = 3 * CenterDims[ iCenter ] + max( max( DimLtotal[ iCenter ], DimRtotal[ iCenter ] ), 4 * CenterDims[ iCenter ] * ( CenterDims[ iCenter ] + 1 ) );
410  double * work = new double[ lwork ];
411  int * iwork = new int[ 8 * CenterDims[ iCenter ] ];
412  int info;
413 
414  // dgesdd is not thread-safe in every implementation ( intel MKL is safe, Atlas is not safe )
415  #pragma omp critical
416  dgesdd_( &jobz, DimLtotal + iCenter, DimRtotal + iCenter, mem, DimLtotal + iCenter,
417  Lambdas[ iCenter ], Us[ iCenter ], DimLtotal + iCenter, VTs[ iCenter ], CenterDims + iCenter, work, &lwork, iwork, &info );
418 
419  delete [] work;
420  delete [] iwork;
421  delete [] mem;
422  }
423  }
424 
425  #ifdef CHEMPS2_MPI_COMPILATION
426  }
427  #endif
428 
429  double discardedWeight = 0.0; // Only if change==true; will the discardedWeight be meaningful and different from zero.
430  int updateSectors = 0;
431  int * NewDims = NULL;
432 
433  // If change: determine new virtual dimensions.
434  #ifdef CHEMPS2_MPI_COMPILATION
435  if (( change ) && ( am_i_master )){
436  #else
437  if ( change ){
438  #endif
439 
440  NewDims = new int[ nCenterSectors ];
441  // First determine the total number of singular values
442  int totalDimSVD = 0;
443  for ( int iCenter = 0; iCenter < nCenterSectors; iCenter++ ){
444  NewDims[ iCenter ] = CenterDims[ iCenter ];
445  totalDimSVD += NewDims[ iCenter ];
446  }
447 
448  // If larger then the required virtualdimensionD, new virtual dimensions will be set in NewDims.
449  if ( totalDimSVD > virtualdimensionD ){
450  // Copy them all in 1 array
451  double * values = new double[ totalDimSVD ];
452  totalDimSVD = 0;
453  int inc = 1;
454  for ( int iCenter = 0; iCenter < nCenterSectors; iCenter++ ){
455  if ( NewDims[ iCenter ] > 0 ){
456  dcopy_( NewDims + iCenter, Lambdas[ iCenter ], &inc, values + totalDimSVD, &inc );
457  totalDimSVD += NewDims[ iCenter ];
458  }
459  }
460 
461  // Sort them in decreasing order
462  char ID = 'D';
463  int info;
464  dlasrt_( &ID, &totalDimSVD, values, &info ); // Quicksort
465 
466  // The D+1'th value becomes the lower bound Schmidt value. Every value smaller than or equal to the D+1'th value is thrown out (hence Dactual <= Ddesired).
467  const double lowerBound = values[ virtualdimensionD ];
468  for ( int iCenter = 0; iCenter < nCenterSectors; iCenter++ ){
469  for ( int cnt = 0; cnt < NewDims[ iCenter ]; cnt++ ){
470  if ( Lambdas[ iCenter ][ cnt ] <= lowerBound ){ NewDims[ iCenter ] = cnt; }
471  }
472  }
473 
474  // Discarded weight
475  double totalSum = 0.0;
476  double discardedSum = 0.0;
477  for ( int iCenter = 0; iCenter < nCenterSectors; iCenter++ ){
478  for ( int iLocal = 0; iLocal < CenterDims[ iCenter ]; iLocal++ ){
479  double temp = ( SplitSectTwoJM[ iCenter ] + 1 ) * Lambdas[ iCenter ][ iLocal ] * Lambdas[ iCenter ][ iLocal ];
480  totalSum += temp;
481  if ( Lambdas[ iCenter ][ iLocal ] <= lowerBound ){ discardedSum += temp; }
482  }
483  }
484  discardedWeight = discardedSum / totalSum;
485 
486  // Clean-up
487  delete [] values;
488  }
489 
490  // Check if there is a sector which differs
491  updateSectors = 0;
492  for ( int iCenter = 0; iCenter < nCenterSectors; iCenter++ ){
493  const int MPSdim = denBK->gCurrentDim( index + 1, SplitSectNM[ iCenter ], SplitSectTwoJM[ iCenter ], SplitSectIM[ iCenter ] );
494  if ( NewDims[ iCenter ] != MPSdim ){ updateSectors = 1; }
495  }
496 
497  }
498 
499  #ifdef CHEMPS2_MPI_COMPILATION
500  MPIchemps2::broadcast_array_int( &updateSectors, 1, MPI_CHEMPS2_MASTER );
501  MPIchemps2::broadcast_array_double( &discardedWeight, 1, MPI_CHEMPS2_MASTER );
502  #endif
503 
504  if ( updateSectors == 1 ){
505 
506  #ifdef CHEMPS2_MPI_COMPILATION
507  if ( NewDims == NULL ){ NewDims = new int[ nCenterSectors ]; }
508  MPIchemps2::broadcast_array_int( NewDims, nCenterSectors, MPI_CHEMPS2_MASTER );
509  #endif
510 
511  for ( int iCenter = 0; iCenter < nCenterSectors; iCenter++ ){
512  denBK->SetDim( index + 1, SplitSectNM[ iCenter ], SplitSectTwoJM[ iCenter ], SplitSectIM[ iCenter ], NewDims[ iCenter ] );
513  }
514  Tleft ->Reset();
515  Tright->Reset();
516 
517  }
518 
519  if ( NewDims != NULL ){ delete [] NewDims; }
520 
521  #ifdef CHEMPS2_MPI_COMPILATION
522  if ( am_i_master ){
523  #endif
524 
525  // Copy first dimM per central symmetry sector to the relevant parts
526  #pragma omp parallel for schedule(dynamic)
527  for ( int iCenter = 0; iCenter < nCenterSectors; iCenter++ ){
528  const int dimM = denBK->gCurrentDim( index + 1, SplitSectNM[ iCenter ], SplitSectTwoJM[ iCenter ], SplitSectIM[ iCenter ] );
529  if ( dimM > 0 ){
530  // U-part: copy
531  int dimLtotal2 = 0;
532  for ( int NL = SplitSectNM[ iCenter ] - 2; NL <= SplitSectNM[ iCenter ]; NL++ ){
533  const int TwoS1 = (( NL + 1 == SplitSectNM[ iCenter ] ) ? 1 : 0 );
534  for ( int TwoSL = SplitSectTwoJM[ iCenter ] - TwoS1; TwoSL <= SplitSectTwoJM[ iCenter ] + TwoS1; TwoSL += 2 ){
535  if ( TwoSL >= 0 ){
536  const int IL = (( TwoS1 == 1 ) ? Irreps::directProd( Ilocal1, SplitSectIM[ iCenter ] ) : SplitSectIM[ iCenter ] );
537  const int dimL = denBK->gCurrentDim( index, NL, TwoSL, IL );
538  if ( dimL > 0 ){
539  double * TleftBlock = Tleft->gStorage( NL, TwoSL, IL, SplitSectNM[ iCenter ], SplitSectTwoJM[ iCenter ], SplitSectIM[ iCenter ] );
540  const int dimension_limit_right = min( dimM, CenterDims[ iCenter ] );
541  for ( int r = 0; r < dimension_limit_right; r++ ){
542  const double factor = (( movingright ) ? 1.0 : Lambdas[ iCenter ][ r ] );
543  for ( int l = 0; l < dimL; l++ ){
544  TleftBlock[ l + dimL * r ] = factor * Us[ iCenter ][ dimLtotal2 + l + DimLtotal[ iCenter ] * r ];
545  }
546  }
547  for ( int r = dimension_limit_right; r < dimM; r++ ){
548  for ( int l = 0; l < dimL; l++ ){
549  TleftBlock[ l + dimL * r ] = 0.0;
550  }
551  }
552  dimLtotal2 += dimL;
553  }
554  }
555  }
556  }
557 
558  // VT-part: copy
559  int dimRtotal2 = 0;
560  for ( int NR = SplitSectNM[ iCenter ]; NR <= SplitSectNM[ iCenter ] + 2; NR++ ){
561  const int TwoS2 = (( NR == SplitSectNM[ iCenter ] + 1 ) ? 1 : 0 );
562  for ( int TwoSR = SplitSectTwoJM[ iCenter ] - TwoS2; TwoSR <= SplitSectTwoJM[ iCenter ] + TwoS2; TwoSR += 2 ){
563  if ( TwoSR >= 0 ){
564  const int IR = (( TwoS2 == 1 ) ? Irreps::directProd( Ilocal2, SplitSectIM[ iCenter ] ) : SplitSectIM[ iCenter ] );
565  const int dimR = denBK->gCurrentDim( index + 2, NR, TwoSR, IR );
566  if ( dimR > 0 ){
567  double * TrightBlock = Tright->gStorage( SplitSectNM[ iCenter ], SplitSectTwoJM[ iCenter ], SplitSectIM[ iCenter ], NR, TwoSR, IR );
568  const int dimension_limit_left = min( dimM, CenterDims[ iCenter ] );
569  const double factor_base = sqrt( ( SplitSectTwoJM[ iCenter ] + 1.0 ) / ( TwoSR + 1 ) );
570  for ( int l = 0; l < dimension_limit_left; l++ ){
571  const double factor = factor_base * (( movingright ) ? Lambdas[ iCenter ][ l ] : 1.0 );
572  for ( int r = 0; r < dimR; r++ ){
573  TrightBlock[ l + dimM * r ] = factor * VTs[ iCenter ][ l + CenterDims[ iCenter ] * ( dimRtotal2 + r ) ];
574  }
575  }
576  for ( int r = 0; r < dimR; r++ ){
577  for ( int l = dimension_limit_left; l < dimM; l++ ){
578  TrightBlock[ l + dimM * r ] = 0.0;
579  }
580  }
581  dimRtotal2 += dimR;
582  }
583  }
584  }
585  }
586  }
587  }
588 
589  #ifdef CHEMPS2_MPI_COMPILATION
590  }
591  MPIchemps2::broadcast_tensor( Tleft, MPI_CHEMPS2_MASTER );
592  MPIchemps2::broadcast_tensor( Tright, MPI_CHEMPS2_MASTER );
593  #endif
594 
595  // Clean up
596  delete [] SplitSectNM;
597  delete [] SplitSectTwoJM;
598  delete [] SplitSectIM;
599  #ifdef CHEMPS2_MPI_COMPILATION
600  if ( am_i_master )
601  #endif
602  {
603  for ( int iCenter = 0; iCenter < nCenterSectors; iCenter++ ){
604  if ( CenterDims[ iCenter ] > 0 ){
605  delete [] Us[ iCenter ];
606  delete [] Lambdas[ iCenter ];
607  delete [] VTs[ iCenter ];
608  }
609  }
610  delete [] Us;
611  delete [] Lambdas;
612  delete [] VTs;
613  delete [] CenterDims;
614  delete [] DimLtotal;
615  delete [] DimRtotal;
616  }
617 
618  return discardedWeight;
619 
620 }
621 
623 
624  #pragma omp parallel for schedule(dynamic)
625  for ( int ikappa = 0; ikappa < nKappa; ikappa++ ){
626 
627  int dim = kappa2index[ ikappa + 1 ] - kappa2index[ ikappa ];
628  double alpha = sqrt( sectorTwoSR[ ikappa ] + 1.0 );
629  int inc = 1;
630  dscal_( &dim, &alpha, storage + kappa2index[ ikappa ], &inc );
631 
632  }
633 
634 }
635 
637 
638  #pragma omp parallel for schedule(dynamic)
639  for ( int ikappa = 0; ikappa < nKappa; ikappa++ ){
640 
641  int dim = kappa2index[ ikappa + 1 ] - kappa2index[ ikappa ];
642  double alpha = 1.0 / sqrt( sectorTwoSR[ ikappa ] + 1.0 );
643  int inc = 1;
644  dscal_( &dim, &alpha, storage + kappa2index[ ikappa ], &inc );
645 
646  }
647 
648 }
649 
650 void CheMPS2::Sobject::addNoise( const double NoiseLevel ){
651 
652  for ( int cnt = 0; cnt < gKappa2index( gNKappa() ); cnt++ ){
653  const double RN = ( ( double ) rand() ) / RAND_MAX - 0.5;
654  gStorage()[ cnt ] += RN * NoiseLevel;
655  }
656 
657 }
658 
659 
int gTwoSL(const int ikappa) const
Get the left spin symmetry of block ikappa.
Definition: Sobject.cpp:203
Sobject(const int index, SyBookkeeper *denBK)
Constructor.
Definition: Sobject.cpp:36
static double wigner6j(const int two_ja, const int two_jb, const int two_jc, const int two_jd, const int two_je, const int two_jf)
Wigner-6j symbol (gsl api)
Definition: Wigner.cpp:294
int gNmin(const int boundary) const
Get the min. possible particle number for a certain boundary.
int gFCIdim(const int boundary, const int N, const int TwoS, const int irrep) const
Get the FCI virtual dimensions ( bound by SYBK_dimensionCutoff )
int gCurrentDim(const int boundary, const int N, const int TwoS, const int irrep) const
Get the current virtual dimensions.
int gKappa2index(const int kappa) const
Get the storage jump corresponding to a certain symmetry block.
Definition: Sobject.cpp:190
double Split(TensorT *Tleft, TensorT *Tright, const int virtualdimensionD, const bool movingright, const bool change)
SVD an S-object into 2 TensorT&#39;s.
Definition: Sobject.cpp:260
static int mpi_rank()
Get the rank of this MPI process.
Definition: MPIchemps2.h:131
void symm2prog()
Convert the storage from symmetric Hamiltonian convention to diagram convention.
Definition: Sobject.cpp:636
static int phase(const int power_times_two)
Phase function.
Definition: Special.h:36
int gTwoJ(const int ikappa) const
Get the central spin symmetry of block ikappa.
Definition: Sobject.cpp:207
int gTwoSmin(const int boundary, const int N) const
Get the minimum possible spin value for a certain boundary and particle number.
int gTwoSmax(const int boundary, const int N) const
Get the maximum possible spin value for a certain boundary and particle number.
int gNKappa() const
Get the number of symmetry blocks.
Definition: Sobject.cpp:166
int gIL(const int ikappa) const
Get the left irrep symmetry of block ikappa.
Definition: Sobject.cpp:204
double * gStorage()
Get the pointer to the storage.
Definition: Sobject.cpp:168
double * gStorage()
Get the pointer to the storage.
Definition: TensorT.cpp:134
int gKappa(const int NL, const int TwoSL, const int IL, const int N1, const int N2, const int TwoJ, const int NR, const int TwoSR, const int IR) const
Get the index corresponding to a certain symmetry block.
Definition: Sobject.cpp:172
void Join(TensorT *Tleft, TensorT *Tright)
Join two sites to form a composite S-object.
Definition: Sobject.cpp:212
static int directProd(const int Irrep1, const int Irrep2)
Get the direct product of the irreps with numbers Irrep1 and Irrep2: a bitwise XOR for psi4&#39;s convent...
Definition: Irreps.h:123
void prog2symm()
Convert the storage from diagram convention to symmetric Hamiltonian convention.
Definition: Sobject.cpp:622
virtual ~Sobject()
Destructor.
Definition: Sobject.cpp:149
int getNumberOfIrreps() const
Get the total number of irreps.
int gNL(const int ikappa) const
Get the left particle number symmetry of block ikappa.
Definition: Sobject.cpp:202
int gIR(const int ikappa) const
Get the right irrep symmetry of block ikappa.
Definition: Sobject.cpp:210
int gTwoSR(const int ikappa) const
Get the right spin symmetry of block ikappa.
Definition: Sobject.cpp:209
int gIndex() const
Get the location index.
Definition: Sobject.cpp:200
int gReorder(const int ikappa) const
Get the blocks from large to small: blocksize(reorder[i])>=blocksize(reorder[i+1]) ...
Definition: Sobject.cpp:170
int gNmax(const int boundary) const
Get the max. possible particle number for a certain boundary.
int gNR(const int ikappa) const
Get the right particle number symmetry of block ikappa.
Definition: Sobject.cpp:208
void SetDim(const int boundary, const int N, const int TwoS, const int irrep, const int value)
Get the current virtual dimensions.
int gIrrep(const int orbital) const
Get an orbital irrep.
void addNoise(const double NoiseLevel)
Add noise to the current S-object.
Definition: Sobject.cpp:650
int gN1(const int ikappa) const
Get the left local particle number symmetry of block ikappa.
Definition: Sobject.cpp:205
void Reset()
Reset the TensorT (if virtual dimensions are changed)
Definition: TensorT.cpp:125
int gN2(const int ikappa) const
Get the right local particle number symmetry of block ikappa.
Definition: Sobject.cpp:206