CheMPS2
EdmistonRuedenberg.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 <stdlib.h>
21 #include <assert.h>
22 #include <iostream>
23 #include <math.h>
24 #include <algorithm>
25 
26 #include "EdmistonRuedenberg.h"
27 #include "DMRGSCFrotations.h"
28 #include "Lapack.h"
29 
30 using std::cout;
31 using std::endl;
32 using std::max;
33 
34 CheMPS2::EdmistonRuedenberg::EdmistonRuedenberg( const FourIndex * Vmat, const int group, const int printLevelIn ){
35 
36  VMAT_ORIG = Vmat;
37  printLevel = printLevelIn;
38  SymmInfo.setGroup( group );
39 
40  int * Isizes = new int[ SymmInfo.getNumberOfIrreps() ];
41  int * Zeroes = new int[ SymmInfo.getNumberOfIrreps() ];
42  int L = 0;
43  for ( int irrep = 0; irrep < SymmInfo.getNumberOfIrreps(); irrep++ ){
44  Isizes[ irrep ] = VMAT_ORIG->get_irrep_size( irrep );
45  Zeroes[ irrep ] = 0;
46  L += Isizes[ irrep ];
47  }
48 
49  iHandler = new DMRGSCFindices( L, group, Zeroes, Isizes, Zeroes ); //Supposes all orbitals are active
50  unitary = new DMRGSCFunitary( iHandler );
51  VmatRotated = new FourIndex( group, Isizes );
52 
53  delete [] Zeroes;
54  delete [] Isizes;
55 
56 }
57 
59 
60  delete unitary; //unitary needs iHandler in its destructor!
61  delete VmatRotated;
62  delete iHandler;
63 
64 }
65 
66 double CheMPS2::EdmistonRuedenberg::Optimize(double * temp1, double * temp2, const bool startFromRandomUnitary, const double gradThreshold, const int maxIter){
67 
68  //Clear the unitary
69  unitary->identity();
70 
71  //Setting up the variables for the gradient
72  double gradNorm = 1.0;
73  const int numVariables = iHandler->getROTparamsize();
74  double * gradient = new double[numVariables];
75  for (int cnt=0; cnt<numVariables; cnt++){ gradient[cnt] = 0.0; }
76 
77  //Randomize the unitary if asked
78  if (startFromRandomUnitary){
79  for (int cnt=0; cnt<numVariables; cnt++){ gradient[cnt] = ((double) rand())/RAND_MAX - 0.5; }
80  unitary->updateUnitary(temp1, temp2, gradient, false, false); //multiply = compact = false
81  for (int cnt=0; cnt<numVariables; cnt++){ gradient[cnt] = 0.0; }
82  }
83 
84  const int mem_size = iHandler->getL() * iHandler->getL() * iHandler->getL() * iHandler->getL();
85  DMRGSCFrotations::rotate( VMAT_ORIG, VmatRotated, NULL, 'F', 'F', 'F', 'F', iHandler, unitary, temp1, temp2, mem_size, "edmistonruedenberg" );
86 
87  //Setting up the variables for the cost function
88  double Icost = costFunction();
89  if (printLevel>0){ cout << " EdmistonRuedenberg::Optimize : Cost function at start = " << Icost << endl; }
90  double Icost_previous = 0.0;
91 
92  int nIterations = 0;
93 
94  while ((gradNorm > gradThreshold) && (nIterations < maxIter)){
95 
96  nIterations++;
97 
98  //Update the unitary
99  unitary->updateUnitary(temp1, temp2, gradient, true, false); //multiply = true; compact = false
100 
101  //Rotate the Vmat
102  Icost_previous = Icost;
103  DMRGSCFrotations::rotate( VMAT_ORIG, VmatRotated, NULL, 'F', 'F', 'F', 'F', iHandler, unitary, temp1, temp2, mem_size, "edmistonruedenberg" );
104  Icost = costFunction();
105 
106  /* What if the cost function has dimished? Then make the rotation step a bit smaller!
107  Rotate back ever closer to unitary_effective = I; untill Icost_previous < Icost
108  1 - 0.5 - 0.5^2 - 0.5^3 - 0.5^4 - ... = 2 - ( 1 + 0.5 + 0.5^2 + ... ) = 2 - 1/(1-0.5) = 0 */
109  if (Icost_previous > Icost){
110  if (printLevel>1){ cout << " WARNING : Icost = " << Icost << " < Icost_previous = " << Icost_previous << endl; }
111  for (int cnt=0; cnt<numVariables; cnt++){ gradient[cnt] = -gradient[cnt]; } //Switch the sign of the update : rotate back!
112  int nIterationsBACK = 0;
113  while ((Icost_previous > Icost) && (nIterationsBACK < CheMPS2::EDMISTONRUED_maxIterBackTfo)){
114  nIterationsBACK++;
115  for (int cnt=0; cnt<numVariables; cnt++){ gradient[cnt] *= 0.5; }
116  unitary->updateUnitary(temp1, temp2, gradient, true, false); //multiply = true; compact = false
117  DMRGSCFrotations::rotate( VMAT_ORIG, VmatRotated, NULL, 'F', 'F', 'F', 'F', iHandler, unitary, temp1, temp2, mem_size, "edmistonruedenberg" );
118  Icost = costFunction();
119  }
120  if (printLevel>1){ cout << " WARNING : Rotated back a bit. Now Icost = " << Icost << endl; }
121  }
122 
123  //Calculate the gradient and Hessian and update
124  gradNorm = augmentedHessianNewtonRaphson(gradient, temp1, temp2);
125 
126  if (printLevel>1){
127  cout << " After iteration " << nIterations << " : Icost = " << Icost << " has gradNorm = " << gradNorm << endl;
128  }
129 
130  }
131 
132  delete [] gradient;
133 
134  if (printLevel>0){
135  cout << " Cost function at stop = " << Icost << endl;
136  cout << " Gradient norm = " << gradNorm << " after " << nIterations << " iterations." << endl;
137  }
138 
139  return Icost;
140 
141 }
142 
144 
145 double CheMPS2::EdmistonRuedenberg::costFunction() const{
146 
147  double Cost = 0.0;
148  for (int irrep=0; irrep<SymmInfo.getNumberOfIrreps(); irrep++){
149  for (int orb=0; orb<iHandler->getNORB(irrep); orb++){
150  Cost += VmatRotated->get(irrep,irrep,irrep,irrep,orb,orb,orb,orb);
151  }
152  }
153  return Cost;
154 
155 }
156 
157 double CheMPS2::EdmistonRuedenberg::calcGradientValue(const int irrep, const int p, const int q) const{
158 
159  return 4 * ( VmatRotated->get(irrep, irrep, irrep, irrep, p, p, p, q) - VmatRotated->get(irrep, irrep, irrep, irrep, q, q, q, p) );
160 
161 }
162 
163 double CheMPS2::EdmistonRuedenberg::calcHessianValue(const int irrep, const int p, const int q, const int r, const int s) const{
164 
165  double value = 0.0;
166  if (p==r){
167  value += ( 8*VmatRotated->get(irrep, irrep, irrep, irrep, p, p, q, s)
168  + 4*VmatRotated->get(irrep, irrep, irrep, irrep, p, q, p, s)
169  - 2*VmatRotated->get(irrep, irrep, irrep, irrep, q, q, q, s)
170  - 2*VmatRotated->get(irrep, irrep, irrep, irrep, s, s, s, q) );
171  }
172  if (q==s){
173  value += ( 8*VmatRotated->get(irrep, irrep, irrep, irrep, q, q, p, r)
174  + 4*VmatRotated->get(irrep, irrep, irrep, irrep, q, p, q, r)
175  - 2*VmatRotated->get(irrep, irrep, irrep, irrep, p, p, p, r)
176  - 2*VmatRotated->get(irrep, irrep, irrep, irrep, r, r, r, p) );
177  }
178  if (p==s){
179  value -= ( 8*VmatRotated->get(irrep, irrep, irrep, irrep, p, p, q, r)
180  + 4*VmatRotated->get(irrep, irrep, irrep, irrep, p, q, p, r)
181  - 2*VmatRotated->get(irrep, irrep, irrep, irrep, q, q, q, r)
182  - 2*VmatRotated->get(irrep, irrep, irrep, irrep, r, r, r, q) );
183  }
184  if (q==r){
185  value -= ( 8*VmatRotated->get(irrep, irrep, irrep, irrep, q, q, p, s)
186  + 4*VmatRotated->get(irrep, irrep, irrep, irrep, q, p, q, s)
187  - 2*VmatRotated->get(irrep, irrep, irrep, irrep, p, p, p, s)
188  - 2*VmatRotated->get(irrep, irrep, irrep, irrep, s, s, s, p) );
189  }
190  return value;
191 
192 }
193 
194 double CheMPS2::EdmistonRuedenberg::augmentedHessianNewtonRaphson(double * gradient, double * temp1, double * temp2) const{
195 
196  int jump = 0;
197  double gradNorm = 0.0;
198 
199  for (int irrep=0; irrep<SymmInfo.getNumberOfIrreps(); irrep++){
200 
201  int linsize = iHandler->getNORB(irrep);
202 
203  if (linsize>1){
204 
205  int hessianlinearsize = linsize*(linsize-1)/2 + 1;
206 
207  /* linsize linsize^4 hesslinsize hesslinsize^2 3*hesslinsize-1 4*hesslinsize-1 (hesslinsize+1)*hesslinsize
208  2 16 2 4 5 7 6
209  3 81 4 16 11 15 20
210  4 256 7 49 20 27 56
211  5 625 11 121 32 43 132
212  6 1296 16 256 47 63 272
213 
214  Conclusion : - temp1 is large enough for hessian (hessianlinearsize^2)
215  - temp2 is large enough for workspace (max(hesslinsize^2, 3*hesslinsize-1)) AND eigenvecs (hesslinsize)
216  */
217 
218  double * hessian = temp1;
219  double * eigen = temp2;
220  double * work = temp2 + hessianlinearsize;
221  int lwork = linsize*linsize*linsize*linsize - hessianlinearsize;
222 
223  /* Maximizing the cost function O(x) = O(0) + g^T * x + x^T * H * x / 2
224  == Minimizing the cost function -O(x) = -O(0) - g^T * x - x^T * H * x / 2 */
225  for (int row=0; row<linsize; row++){
226  for (int col=row+1; col<linsize; col++){
227  const double waarde = calcGradientValue(irrep, row, col);
228  gradient[jump + row + col*(col-1)/2] = waarde;
229  gradNorm += waarde * waarde;
230  for (int row2=0; row2<linsize; row2++){
231  for (int col2=row2+1; col2<linsize; col2++){
232  const double value = calcHessianValue(irrep, row, col, row2, col2);
233  hessian[row + col*(col-1)/2 + hessianlinearsize * (row2 + col2*(col2-1)/2)] = - value; //MINUS SIGN FOR OTHER DIRECTION
234  }
235  }
236  }
237  }
238 
239  for (int cnt=0; cnt<hessianlinearsize-1; cnt++){
240  hessian[hessianlinearsize-1 + hessianlinearsize * cnt] = - gradient[jump + cnt]; //MINUS SIGN FOR OTHER DIRECTION
241  hessian[cnt + hessianlinearsize * (hessianlinearsize-1)] = - gradient[jump + cnt]; //MINUS SIGN FOR OTHER DIRECTION
242  }
243  hessian[hessianlinearsize*hessianlinearsize-1] = 0.0;
244 
245  //Find its lowest eigenvalue and vector
246  char jobz = 'V';
247  char uplo = 'U';
248  int info;
249  dsyev_(&jobz,&uplo,&hessianlinearsize,hessian,&hessianlinearsize,eigen,work,&lwork,&info);
250 
251  double scalar = 1.0/hessian[hessianlinearsize-1];
252  int inc = 1;
253  dscal_(&hessianlinearsize,&scalar,hessian,&inc);
254 
255  for (int cnt=0; cnt<hessianlinearsize-1; cnt++){ gradient[jump + cnt] = hessian[cnt]; }
256 
257  jump += hessianlinearsize-1;
258 
259  }
260 
261  }
262 
263  gradNorm = sqrt(gradNorm);
264 
265  return gradNorm;
266 
267 }
268 
269 void CheMPS2::EdmistonRuedenberg::Fiedler(const int irrep, int * reorder, double * laplacian, double * temp2){
270 
271  //For information on the Fiedler vector: see http://web.eecs.utk.edu/~mberry/order/node9.html
272 
273  int linsize = iHandler->getNORB(irrep);
274  assert( linsize>=2 );
275 
276  //Preamble: linsize>=2 at this point
277  double * work = temp2; //temp2 at least of size 4*linsize*linsize
278  int lwork = 3*linsize*linsize; //3*linsize*linsize > 3*linsize-1 for linsize>=2
279  double * eigs = temp2 + lwork; //has size linsize*linsize > linsize for linsize>=2
280 
281  //Calculate the eigenspectrum of the Laplacian
282  char jobz = 'V';
283  char uplo = 'U';
284  int info;
285  dsyev_(&jobz, &uplo, &linsize, laplacian, &linsize, eigs, work, &lwork, &info);
286  if (printLevel>1){
287  cout << " EdmistonRuedenberg::Fiedler : Smallest eigs(Laplacian[" << irrep << "]) = [ " << eigs[0] << " , " << eigs[1] << " ]." << endl;
288  }
289 
290  //The second eigenvector is the Fiedler vector
291  double * FiedlerVector = laplacian + linsize;
292  double * FiedlerVectorCopy = laplacian;
293  for (int cnt=0; cnt<linsize; cnt++){ FiedlerVectorCopy[cnt] = FiedlerVector[cnt]; }
294 
295  //Find the Fiedler ordering
296  const double upperBound = 2.0; //Sum_i FiedlerVector[i]^2 = 1 --> fabs(FiedlerVector[i]) <= 1.0
297  for (int value=0; value<linsize; value++){
298  int index=0;
299  for (int cnt=1; cnt<linsize; cnt++){
300  if (FiedlerVectorCopy[cnt] < FiedlerVectorCopy[index]){ index = cnt; }
301  } //Index now corresponds to the smallest value in FiedlerVectorCopy
302  reorder[value] = index;
303  FiedlerVectorCopy[index] = upperBound; //Remove from the interesting entries in FiedlerVectorCopy
304  }
305 
306  if (printLevel>1){
307  bool isOK = true;
308  for (int cnt=0; cnt<linsize-1; cnt++){
309  if (FiedlerVector[reorder[cnt]] > FiedlerVector[reorder[cnt+1]]){ isOK = false; }
310  }
311  assert( isOK );
312  cout << " Reorder[" << irrep << "] = [ ";
313  for (int cnt=0; cnt<linsize-1; cnt++){ cout << reorder[cnt] << " , "; }
314  cout << reorder[linsize-1] << " ]." << endl;
315  }
316 
317  //Reorder the vectors in the unitary
318  double * blockU = unitary->getBlock(irrep);
319  for (int row=0; row<linsize; row++){
320  for (int col=0; col<linsize; col++){
321  work[row + linsize * col] = blockU[reorder[row] + linsize * col];
322  }
323  }
324 
325  int size = linsize*linsize;
326  int inc = 1;
327  dcopy_(&size, work, &inc, blockU, &inc);
328 
329  if (printLevel>1){
330  char trans = 'T';
331  char notrans = 'N';
332  double alpha = 1.0;
333  double beta = 0.0;
334  dgemm_(&trans, &notrans, &linsize, &linsize, &linsize, &alpha, blockU, &linsize, blockU, &linsize, &beta, work, &linsize);
335  double sum = 0.0;
336  for (int row=0; row<linsize; row++){
337  sum += (work[row+linsize*row]-1.0)*(work[row+linsize*row]-1.0);
338  for (int col=row+1; col<linsize; col++){
339  sum += work[row+linsize*col]*work[row+linsize*col] + work[col+linsize*row]*work[col+linsize*row];
340  }
341  }
342  sum = sqrt(sum);
343  cout << " 2-norm of Unitary[" << irrep << "]^T * Unitary[" << irrep << "] - I = " << sum << endl;
344  }
345 
346 }
347 
348 double CheMPS2::EdmistonRuedenberg::FiedlerExchangeCost() const{
349 
350  double Cost = 0.0;
351 
352  for (int irrep=0; irrep<SymmInfo.getNumberOfIrreps(); irrep++){
353  const int linsize = iHandler->getNORB(irrep);
354  if (linsize>1){
355  for (int row=0; row<linsize; row++){
356  for (int col=row+1; col<linsize; col++){
357  Cost += 2 * VmatRotated->get(irrep,irrep,irrep,irrep,row,col,col,row) * (col-row) * (col-row);
358  }
359  }
360  }
361  }
362 
363  return Cost;
364 
365 }
366 
367 void CheMPS2::EdmistonRuedenberg::FiedlerExchange(const int maxlinsize, double * temp1, double * temp2){
368 
369  //For information on the Fiedler vector: see http://web.eecs.utk.edu/~mberry/order/node9.html
370 
371  const int mem_size = iHandler->getL() * iHandler->getL() * iHandler->getL() * iHandler->getL();
372  DMRGSCFrotations::rotate( VMAT_ORIG, VmatRotated, NULL, 'F', 'F', 'F', 'F', iHandler, unitary, temp1, temp2, mem_size, "edmistonruedenberg" );
373 
374  if (printLevel>0){ cout << " EdmistonRuedenberg::FiedlerExchange : Cost function at start = " << FiedlerExchangeCost() << endl; }
375 
376  int * reorder = new int[maxlinsize];
377 
378  for (int irrep=0; irrep<SymmInfo.getNumberOfIrreps(); irrep++){
379 
380  const int linsize = iHandler->getNORB(irrep);
381  if (linsize>1){
382 
383  //temp1 and temp2 both at least of size 4*linsize*linsize
384  double * laplacian = temp1;
385 
386  //Fill the weighted graph Laplacian
387  for (int row=0; row<linsize; row++){
388  laplacian[row*(1+linsize)] = 0.0;
389  for (int col=row+1; col<linsize; col++){
390  laplacian[row + linsize*col] = - VmatRotated->get(irrep,irrep,irrep,irrep,row,col,col,row); //minus the exchange matrix
391  laplacian[col + linsize*row] = laplacian[row + linsize*col]; //Symmetric matrix
392  laplacian[row + linsize*row] -= laplacian[row + linsize*col]; //On diagonal : Minus the sum of the other terms on that row
393  }
394  for (int col=0; col<row; col++){
395  laplacian[row + linsize*row] -= laplacian[row + linsize*col]; //On diagonal : Minus the sum of the other terms on that row
396  }
397  }
398 
399  Fiedler(irrep, reorder, laplacian, temp2);
400 
401  }
402  }
403 
404  delete [] reorder;
405 
406  DMRGSCFrotations::rotate( VMAT_ORIG, VmatRotated, NULL, 'F', 'F', 'F', 'F', iHandler, unitary, temp1, temp2, mem_size, "edmistonruedenberg" );
407 
408  if (printLevel>0){ cout << " EdmistonRuedenberg::FiedlerExchange : Cost function at end = " << FiedlerExchangeCost() << endl; }
409 
410 }
411 
412 double CheMPS2::EdmistonRuedenberg::FiedlerGlobalCost( const DMRGSCFindices * idx, const FourIndex * VMAT_LOCAL, int * dmrg2ham ){
413 
414  double cost = 0.0;
415 
416  for ( int dmrg_row = 0; dmrg_row < idx->getL(); dmrg_row++ ){
417  for ( int dmrg_col = 0; dmrg_col < idx->getL(); dmrg_col++ ){
418  const int ham_row = dmrg2ham[ dmrg_row ];
419  const int ham_col = dmrg2ham[ dmrg_col ];
420  const int irrep_row = idx->getOrbitalIrrep( ham_row );
421  const int irrep_col = idx->getOrbitalIrrep( ham_col );
422  const int rel_row = ham_row - idx->getOrigNOCCstart( irrep_row );
423  const int rel_col = ham_col - idx->getOrigNOCCstart( irrep_col );
424  cost += VMAT_LOCAL->get( irrep_row, irrep_col, irrep_col, irrep_row, rel_row, rel_col, rel_col, rel_row ) * ( dmrg_row - dmrg_col ) * ( dmrg_row - dmrg_col );
425  }
426  }
427 
428  return cost;
429 
430 }
431 
432 void CheMPS2::EdmistonRuedenberg::FiedlerGlobal( int * dmrg2ham ) const{
433 
434  // For information on the Fiedler vector: see http://web.eecs.utk.edu/~mberry/order/node9.html
435 
436  for ( int orb = 0; orb < iHandler->getL(); orb++ ){ dmrg2ham[ orb ] = orb; }
437  if ( printLevel > 0 ){ cout << " EdmistonRuedenberg::FiedlerGlobal : Cost function at start = " << FiedlerGlobalCost( iHandler, VMAT_ORIG, dmrg2ham ) << endl; }
438 
439  // Build the Laplacian
440  double * laplacian = new double[ iHandler->getL() * iHandler->getL() ];
441  for ( int ham_row = 0; ham_row < iHandler->getL(); ham_row++ ){
442  double sum_over_column = 0.0;
443  for ( int ham_col = 0; ham_col < iHandler->getL(); ham_col++ ){
444  if ( ham_row != ham_col ){
445  const int irrep_row = iHandler->getOrbitalIrrep( ham_row );
446  const int irrep_col = iHandler->getOrbitalIrrep( ham_col );
447  const int rel_row = ham_row - iHandler->getOrigNOCCstart( irrep_row );
448  const int rel_col = ham_col - iHandler->getOrigNOCCstart( irrep_col );
449  const double value = fabs( VMAT_ORIG->get( irrep_row, irrep_col, irrep_col, irrep_row, rel_row, rel_col, rel_col, rel_row ) );
450  laplacian[ ham_row + iHandler->getL() * ham_col ] = - value;
451  sum_over_column += value;
452  } else {
453  laplacian[ ham_row + iHandler->getL() * ham_col ] = 0.0;
454  }
455  }
456  laplacian[ ham_row + iHandler->getL() * ham_row ] = sum_over_column;
457  }
458 
459  // Calculate the eigenspectrum of the Laplacian
460  int lwork = 3 * iHandler->getL() * iHandler->getL();
461  double * work = new double[ lwork ];
462  double * eigs = new double[ iHandler->getL() ];
463  char jobz = 'V';
464  char uplo = 'U';
465  int linsize = iHandler->getL();
466  int info;
467  dsyev_( &jobz, &uplo, &linsize, laplacian, &linsize, eigs, work, &lwork, &info );
468  delete [] work;
469  delete [] eigs;
470 
471  // Fill dmrg2ham
472  double * fiedler_vec = laplacian + linsize;
473  for ( int dummy = 0; dummy < linsize; dummy++ ){
474  int index = 0;
475  for ( int orb = 1; orb < linsize; orb++ ){
476  if ( fiedler_vec[ orb ] < fiedler_vec[ index ] ){ index = orb; }
477  }
478  dmrg2ham[ dummy ] = index;
479  fiedler_vec[ index ] = 2.0; // Eigenvectors are normalized to 1.0, so certainly OK
480  }
481 
482  delete [] laplacian;
483 
484  if ( printLevel > 0 ){
485  cout << " EdmistonRuedenberg::FiedlerGlobal : Cost function at end = " << FiedlerGlobalCost( iHandler, VMAT_ORIG, dmrg2ham ) << endl;
486  cout << " EdmistonRuedenberg::FiedlerGlobal : Reordering = [ ";
487  for ( int orb = 0; orb < linsize - 1; orb++ ){ cout << dmrg2ham[ orb ] << ", "; }
488  cout << dmrg2ham[ linsize - 1 ] << " ]." << endl;
489  }
490 
491 }
492 
493 
int getNORB(const int irrep) const
Get the number of orbitals for an irrep.
double * getBlock(const int irrep)
Get a matrix block.
bool setGroup(const int nGroup)
Set the group.
Definition: Irreps.cpp:50
static void rotate(const FourIndex *ORIG_VMAT, FourIndex *NEW_VMAT, DMRGSCFintegrals *ROT_TEI, const char space1, const char space2, const char space3, const char space4, DMRGSCFindices *idx, DMRGSCFunitary *umat, double *mem1, double *mem2, const int mem_size, const string filename)
Fill the rotated two-body matrix elements for the space. If the blocks become too large...
int getROTparamsize() const
Get the orbital rotation parameter space size.
int getNumberOfIrreps() const
Get the number of irreps for the currently activated group.
Definition: Irreps.cpp:94
virtual ~EdmistonRuedenberg()
Destructor.
double get(const int irrep_i, const int irrep_j, const int irrep_k, const int irrep_l, const int i, const int j, const int k, const int l) const
Get an element.
Definition: FourIndex.cpp:175
EdmistonRuedenberg(const FourIndex *Vmat, const int group, const int printLevelIn=1)
Constructor.
void FiedlerGlobal(int *dmrg2ham) const
Permute the orbitals so that the bandwidth of the exchange matrix is (approximately) minimized...
int getL() const
Get the number of orbitals.
int get_irrep_size(const int irrep) const
Get a given irrep size.
Definition: FourIndex.cpp:294
DMRGSCFunitary * getUnitary()
Get the pointer to the unitary to use in DMRGSCF.
void FiedlerExchange(const int maxlinsize, double *temp1, double *temp2)
Permute the orbitals so that the bandwidth of the exchange matrix is (approximately) minimized (per i...
double Optimize(double *temp1, double *temp2, const bool startFromRandomUnitary, const double gradThreshold=EDMISTONRUED_gradThreshold, const int maxIter=EDMISTONRUED_maxIter)
Maximize the Edmiston-Ruedenberg cost function.
void identity()
Make this matrix the identity matrix.
void updateUnitary(double *workmem1, double *workmem2, double *vector, const bool multiply, const bool compact)
Update the unitary transformation.
int getOrbitalIrrep(const int index) const
Get the irrep corresponding to a global orbital index.
int getOrigNOCCstart(const int irrep) const
Get in the original Hamiltonian index the start orbital for the occupied orbitals with a certain irre...