CheMPS2
DMRGSCFunitary.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 
25 #include "Lapack.h"
26 #include "DMRGSCFunitary.h"
27 
28 using std::cout;
29 using std::endl;
30 
32 
33  this->identity();
34 
35  //Find the unique indices for OCC-ACT, OCC-VIRT, and ACT-VIRT rotations
36  x_linearlength = 0;
37  jumper = new int*[ num_irreps ];
38  for ( int irrep = 0; irrep < num_irreps; irrep++ ){
39  jumper[ irrep ] = new int[ 3 ];
40  const int NOCC = iHandler->getNOCC( irrep );
41  const int NACT = iHandler->getNDMRG( irrep );
42  const int NVIR = iHandler->getNVIRT( irrep );
43  jumper[ irrep ][ 0 ] = x_linearlength;
44  x_linearlength += NOCC * NACT;
45  jumper[ irrep ][ 1 ] = x_linearlength;
46  x_linearlength += NACT * NVIR;
47  jumper[ irrep ][ 2 ] = x_linearlength;
48  x_linearlength += NOCC * NVIR;
49  }
50 
51 }
52 
54 
55  for ( int irrep = 0; irrep < num_irreps; irrep++ ){ delete [] jumper[ irrep ]; }
56  delete [] jumper;
57 
58 }
59 
60 int CheMPS2::DMRGSCFunitary::getNumVariablesX() const{ return x_linearlength; }
61 
62 void CheMPS2::DMRGSCFunitary::buildSkewSymmX( const int irrep, double * result, double * Xelem, const bool compact ) const{
63 
64  const int linsize = iHandler->getNORB( irrep );
65  const int NOCC = iHandler->getNOCC( irrep );
66  const int NACT = iHandler->getNDMRG( irrep );
67  const int NVIR = iHandler->getNVIRT( irrep );
68  for ( int cnt = 0; cnt < linsize * linsize; cnt++ ){ result[ cnt ] = 0.0; }
69 
70  if ( compact ){
71 
72  for ( int occ = 0; occ < NOCC; occ++ ){
73  for ( int act = 0; act < NACT; act++ ){
74  const int xsolindex = jumper[ irrep ][ 0 ] + act + NACT * occ;
75  const int index1 = NOCC + act;
76  result[ index1 + linsize * occ ] = Xelem[ xsolindex ];
77  result[ occ + linsize * index1 ] = - Xelem[ xsolindex ];
78  }
79  }
80  for ( int act = 0; act < NACT; act++ ){
81  for ( int vir = 0; vir < NVIR; vir++ ){
82  const int xsolindex = jumper[ irrep ][ 1 ] + vir + NVIR * act;
83  const int index1 = NOCC + NACT + vir;
84  const int index2 = NOCC + act;
85  result[ index1 + linsize * index2 ] = Xelem[ xsolindex ];
86  result[ index2 + linsize * index1 ] = - Xelem[ xsolindex ];
87  }
88  }
89  for ( int occ = 0; occ < NOCC; occ++ ){
90  for ( int vir = 0; vir < NVIR; vir++ ){
91  const int xsolindex = jumper[ irrep ][ 2 ] + vir + NVIR * occ;
92  const int index1 = NOCC + NACT + vir;
93  result[ index1 + linsize * occ ] = Xelem[ xsolindex ];
94  result[ occ + linsize * index1 ] = - Xelem[ xsolindex ];
95  }
96  }
97 
98  } else { //NOT compact
99 
100  int jump = 0;
101  for ( int cnt = 0; cnt < irrep; cnt++ ){
102  const int NORBx = iHandler->getNORB( cnt );
103  jump += ( NORBx * ( NORBx - 1 ) ) / 2;
104  }
105 
106  for ( int row = 0; row < linsize; row++ ){
107  for ( int col = row+1; col < linsize; col++ ){
108  result[ row + linsize * col ] = Xelem[ jump + row + ( col * ( col - 1 ) ) / 2 ];
109  result[ col + linsize * row ] = - Xelem[ jump + row + ( col * ( col - 1 ) ) / 2 ];
110  }
111  }
112 
113  }
114 
115 }
116 
117 void CheMPS2::DMRGSCFunitary::updateUnitary( double * temp1, double * temp2, double * vector, const bool multiply, const bool compact ){
118 
119  for ( int irrep = 0; irrep < num_irreps; irrep++ ){
120 
121  int linsize = iHandler->getNORB(irrep);
122  int size = linsize * linsize;
123 
124  if ( linsize > 1 ){ //linsize is op z'n minst 2 dus temp1, temp1+size, temp1+2*size,temp1+3*size zijn zeker ok
125 
126  double * xblock = temp1; //linsize*linsize
127  double * Bmat = temp1 + size; //linsize*linsize
128  double * work1 = temp1 + 2 * size; //linsize*linsize
129  double * work2 = temp1 + 3 * size; //linsize*linsize
130 
131  double * workLARGE = temp2; //4*size
132  int lworkLARGE = 4 * size; //4*size = 4*linsize*linsize > 3*linsize-1
133 
134  // Construct the antisymmetric x-matrix
135  buildSkewSymmX( irrep, xblock, vector, compact );
136 
137  //Bmat <= xblock * xblock
138  char notrans = 'N';
139  double one = 1.0;
140  double set = 0.0;
141  dgemm_( &notrans, &notrans, &linsize, &linsize, &linsize, &one, xblock, &linsize, xblock, &linsize, &set, Bmat, &linsize );
142 
143  //Bmat * work1 * Bmat^T <= xblock * xblock
144  char uplo = 'U';
145  char jobz = 'V';
146  int info;
147  dsyev_( &jobz, &uplo, &linsize, Bmat, &linsize, work1, workLARGE, &lworkLARGE, &info );
148 
149  //work2 <= Bmat^T * xblock * Bmat
150  dgemm_( &notrans, &notrans, &linsize, &linsize, &linsize, &one, xblock, &linsize, Bmat, &linsize, &set, work1, &linsize );
151  char trans = 'T';
152  dgemm_( &trans, &notrans, &linsize, &linsize, &linsize, &one, Bmat, &linsize, work1, &linsize, &set, work2, &linsize );
153 
154  if (CheMPS2::DMRGSCF_debugPrint){
155  cout << " DMRGSCFunitary::updateUnitary : Lambdas of irrep block " << irrep << " : " << endl;
156  for (int cnt=0; cnt<linsize/2; cnt++){
157  cout << " block = [ " << work2[2*cnt + linsize*2*cnt] << " , " << work2[2*cnt + linsize*(2*cnt+1)] << " ] " << endl;
158  cout << " [ " << work2[2*cnt+1 + linsize*2*cnt] << " , " << work2[2*cnt+1 + linsize*(2*cnt+1)] << " ] " << endl;
159  }
160  }
161 
162  //work1 <= values of the antisymmetric 2x2 blocks
163  for ( int cnt = 0; cnt < linsize/2; cnt++ ){
164  work1[cnt] = 0.5 * ( work2[ 2 * cnt + linsize * ( 2 * cnt + 1 ) ]
165  - work2[ 2 * cnt + 1 + linsize * ( 2 * cnt ) ] );
166  }
167 
168  if ( CheMPS2::DMRGSCF_debugPrint ){
169  for ( int cnt = 0; cnt < linsize/2; cnt++ ){
170  work2[ 2 * cnt + linsize * ( 2 * cnt + 1 ) ] -= work1[ cnt ];
171  work2[ 2 * cnt + 1 + linsize * ( 2 * cnt ) ] += work1[ cnt ];
172  }
173  double rms_diff = 0.0;
174  for ( int cnt = 0; cnt < size; cnt++ ){ rms_diff += work2[ cnt ] * work2[ cnt ]; }
175  rms_diff = sqrt( rms_diff );
176  cout << " DMRGSCFunitary::updateUnitary : RMSdeviation of irrep block " << irrep << " (should be 0.0) = " << rms_diff << endl;
177  }
178 
179  //work2 <= exp(Bmat^T * xblock * Bmat)
180  for ( int cnt = 0; cnt < size; cnt++ ){ work2[ cnt ] = 0.0; }
181  for ( int cnt = 0; cnt < linsize/2; cnt++ ){
182  double cosine = cos( work1[ cnt ] );
183  double sine = sin( work1[ cnt ] );
184  work2[ 2 * cnt + linsize * ( 2 * cnt ) ] = cosine;
185  work2[ 2 * cnt + 1 + linsize * ( 2 * cnt + 1) ] = cosine;
186  work2[ 2 * cnt + linsize * ( 2 * cnt + 1) ] = sine;
187  work2[ 2 * cnt + 1 + linsize * ( 2 * cnt ) ] = -sine;
188  }
189  for ( int cnt = 2*(linsize/2); cnt < linsize; cnt++ ){ work2[ cnt * (linsize + 1) ] = 1.0; }
190 
191  //work2 <= Bmat * exp(Bmat^T * xblock * Bmat) * Bmat^T = exp(xblock)
192  dgemm_( &notrans, &notrans, &linsize, &linsize, &linsize, &one, Bmat, &linsize, work2, &linsize, &set, work1, &linsize );
193  dgemm_( &notrans, &trans, &linsize, &linsize, &linsize, &one, work1, &linsize, Bmat, &linsize, &set, work2, &linsize );
194 
195  int inc1 = 1;
196  if ( multiply ){ //U <-- exp(x) * U
197  dgemm_( &notrans, &notrans, &linsize, &linsize, &linsize, &one, work2, &linsize, entries[ irrep ], &linsize, &set, work1, &linsize );
198  dcopy_( &size, work1, &inc1, entries[ irrep ], &inc1 );
199  } else { //U <-- exp(x)
200  dcopy_( &size, work2, &inc1, entries[ irrep ], &inc1 );
201  }
202 
203  }
204  }
205 
206  if (CheMPS2::DMRGSCF_debugPrint){ CheckDeviationFromUnitary( temp2 ); }
207 
208 }
209 
210 void CheMPS2::DMRGSCFunitary::rotateActiveSpaceVectors( double * eigenvecs, double * work ){
211 
212  int passed = 0;
213  int tot_dmrg = iHandler->getDMRGcumulative( num_irreps );
214  for ( int irrep = 0; irrep < num_irreps; irrep++ ){
215 
216  int NDMRG = iHandler->getNDMRG( irrep );
217  int NOCC = iHandler->getNOCC( irrep );
218  int NORB = iHandler->getNORB( irrep );
219 
220  if ( NDMRG > 1 ){
221 
222  double * rotation = eigenvecs + passed * ( tot_dmrg + 1 );
223 
224  char trans = 'T';
225  char notrans = 'N';
226  double one = 1.0;
227  double set = 0.0;
228  dgemm_( &trans, &notrans, &NDMRG, &NORB, &NDMRG, &one, rotation, &tot_dmrg, entries[ irrep ] + NOCC, &NORB, &set, work, &NDMRG );
229 
230  for ( int row = 0; row < NDMRG; row++ ){
231  for ( int col = 0; col < NORB; col++ ){
232  entries[ irrep ][ NOCC + row + NORB * col ] = work[ row + NDMRG * col ];
233  }
234  }
235 
236  }
237  passed += NDMRG;
238  }
239 
240  if (CheMPS2::DMRGSCF_debugPrint){ CheckDeviationFromUnitary( work ); }
241 
242 }
243 
245 
246  char trans = 'T';
247  char notrans = 'N';
248  double one = 1.0;
249  double set = 0.0;
250 
251  for ( int irrep = 0; irrep < num_irreps; irrep++ ){
252 
253  int linsize = iHandler->getNORB( irrep );
254  if ( linsize > 0 ){
255  dgemm_( &trans, &notrans, &linsize, &linsize, &linsize, &one, entries[ irrep ], &linsize, entries[ irrep ], &linsize, &set, work, &linsize );
256  for ( int diag = 0; diag < linsize; diag++ ){
257  work[ diag * ( 1 + linsize ) ] -= 1.0;
258  }
259  double rms_diff = 0.0;
260  for ( int cnt = 0; cnt < linsize * linsize; cnt++ ){
261  rms_diff += work[ cnt ] * work[ cnt ];
262  }
263  rms_diff = sqrt( rms_diff );
264  cout << " DMRGSCFunitary::CheckDeviationFromUnitary : 2-norm of U[" << irrep << "]^T * U[" << irrep << "] - I = " << rms_diff << endl;
265  }
266  }
267 
268 }
269 
270 double CheMPS2::DMRGSCFunitary::get_determinant( const int irrep, double * work1, double * work2, double * work_eig, int lwork_eig ) const{
271 
272  int NORB = iHandler->getNORB( irrep );
273 
274  // S = U + U^T --> work1
275  for ( int row = 0; row < NORB; row++ ){
276  for ( int col = 0; col < NORB; col++ ){
277  work1[ row + NORB * col ] = ( entries[ irrep ][ row + NORB * col ]
278  + entries[ irrep ][ col + NORB * row ] );
279  }
280  }
281 
282  // Get the eigenvectors of S = U + U^T: eigvals in work2, eigvecs in work1
283  char jobz = 'V';
284  char uplo = 'U';
285  int info;
286  dsyev_( &jobz, &uplo, &NORB, work1, &NORB, work2, work_eig, &lwork_eig, &info );
287 
288  // Transform the original orthogonal matrix with the real orthonormal eigenbasis of S --> the result V^T U V is stored in work2
289  char trans = 'T';
290  char notrans = 'N';
291  double one = 1.0;
292  double set = 0.0;
293  dgemm_( &trans, &notrans, &NORB, &NORB, &NORB, &one, work1, &NORB, entries[ irrep ], &NORB, &set, work_eig, &NORB );
294  dgemm_( &notrans, &notrans, &NORB, &NORB, &NORB, &one, work_eig, &NORB, work1, &NORB, &set, work2, &NORB );
295 
296  /* Work2 should contain
297  * 1x1 blocks containing [-1]
298  * 2x2 blocks containing [[cos(u) sin(u)][-sin(u) cos(u)]]
299  * 1x1 blocks containing [+1]
300  and is hence TRIDIAGONAL. Its determinant is easy to compute:
301  */
302  double f_low = 1.0;
303  double f_high = work2[ 0 ];
304  for ( int diag = 1; diag < NORB; diag++ ){
305  double temp = work2[ diag * ( 1 + NORB ) ] * f_high - work2[ diag + NORB * ( diag - 1 ) ] * work2[ diag - 1 + NORB * diag ] * f_low;
306  f_low = f_high;
307  f_high = temp;
308  }
309  const double determinant = f_high;
310 
311  if ( CheMPS2::DMRGSCF_debugPrint ){
312  cout << " DMRGSCFunitary::get_determinant : det( U[" << irrep << "] ) = " << determinant << endl;
313  }
314 
315  return determinant;
316 
317 }
318 
319 void CheMPS2::DMRGSCFunitary::getLog( double * vector, double * temp1, double * temp2 ) const{
320 
321  int jump = 0;
322 
323  for ( int irrep = 0; irrep < num_irreps; irrep++ ){
324 
325  int linsize = iHandler->getNORB( irrep );
326  int size = linsize * linsize;
327 
328  if ( linsize > 1 ){
329  /* linsize >= 2; hence temp1 is at least of size 4*linsize*linsize
330  if linsize <= 1; there corresponds no block in xmatrix to it */
331 
332  double * work1 = temp1; //linsize * linsize
333  double * work2 = temp1 + size; //linsize * linsize
334  double * work3 = temp1 + 2 * size; //linsize * linsize
335  double * work4 = temp1 + 3 * size; //linsize * linsize
336 
337  double * workLARGE = temp2;
338  int lworkLARGE = 4 * size;
339 
340  // work1 contains eigvec( U + U^T )
341  // work3 contains eigvec( U + U^T )^T * U * eigvec( U + U^T ) ( TRIDIAGONAL matrix )
342  const double determinant = get_determinant( irrep, work1, work3, workLARGE, lworkLARGE );
343  assert( determinant > 0.0 );
344 
345  //Fill work2 with ln(V^T U V) = ln(work3).
346  for ( int cnt = 0; cnt < size; cnt++ ){ work2[ cnt ] = 0.0; }
347  for ( int cnt = 0; cnt < linsize/2; cnt++ ){ //2x2 blocks
348  const double cosinus = 0.5 * ( work3[ 2 * cnt + linsize * ( 2 * cnt ) ] + work3[ 2 * cnt + 1 + linsize * ( 2 * cnt + 1 ) ] );
349  const double sinus = 0.5 * ( work3[ 2 * cnt + linsize * ( 2 * cnt + 1 ) ] - work3[ 2 * cnt + 1 + linsize * ( 2 * cnt ) ] );
350  const double theta = atan2( sinus, cosinus );
351  if (CheMPS2::DMRGSCF_debugPrint){
352  cout << " DMRGSCFunitary::getLog : block = [ " << work3[ 2 * cnt + linsize * ( 2 * cnt ) ] << " , "
353  << work3[ 2 * cnt + linsize * ( 2 * cnt + 1 ) ] << " ]" << endl;
354  cout << " [ " << work3[ 2 * cnt + 1 + linsize * ( 2 * cnt ) ] << " , "
355  << work3[ 2 * cnt + 1 + linsize * ( 2 * cnt + 1 ) ] << " ] and corresponds to theta = " << theta << endl;
356  }
357  work3[ 2 * cnt + linsize * ( 2 * cnt ) ] -= cosinus;
358  work3[ 2 * cnt + 1 + linsize * ( 2 * cnt + 1 ) ] -= cosinus;
359  work3[ 2 * cnt + linsize * ( 2 * cnt + 1 ) ] -= sinus;
360  work3[ 2 * cnt + 1 + linsize * ( 2 * cnt ) ] += sinus;
361  work2[ 2 * cnt + linsize * ( 2 * cnt + 1 ) ] = theta;
362  work2[ 2 * cnt + 1 + linsize * ( 2 * cnt ) ] = -theta;
363  } //The rest are 1x1 blocks, corresponding to eigenvalue 1 --> ln(1) = 0 --> work2 does not need to be updated.
364  for ( int cnt = 2*(linsize/2); cnt < linsize; cnt++ ){ work3[ cnt * ( linsize + 1 ) ] -= 1; }
365 
366  //Calculate the 2-norm of updated work3 (should be 0.0)
367  if (CheMPS2::DMRGSCF_debugPrint){
368  double RMSdeviation = 0.0;
369  for ( int cnt = 0; cnt < size; cnt++ ){ RMSdeviation += work3[ cnt ] * work3[ cnt ]; }
370  RMSdeviation = sqrt( RMSdeviation );
371  cout << " DMRGSCFunitary::getLog : 2-norm of [ V^T*U*V - assumed blocks ] for irrep " << irrep << " (should be 0.0) = " << RMSdeviation << endl;
372  }
373 
374  //Calculate V * ln(blocks) * V^T --> work4
375  char trans = 'T';
376  char notrans = 'N';
377  double one = 1.0;
378  double set = 0.0;
379  dgemm_(&notrans, &notrans, &linsize, &linsize, &linsize, &one, work1, &linsize, work2, &linsize, &set, workLARGE, &linsize ); //V*ln(blocks)
380  dgemm_(&notrans, &trans, &linsize, &linsize, &linsize, &one, workLARGE, &linsize, work1, &linsize, &set, work4, &linsize ); //V*ln(blocks)*V^T
381 
382  //Fill the vector with the just calculated ln(U) = work4
383  for ( int row = 0; row < linsize; row++ ){
384  for ( int col = row + 1; col < linsize; col++ ){
385  vector[ jump + row + ( col * ( col - 1 ) ) / 2 ] = 0.5 * ( work4[ row + linsize * col ] - work4[ col + linsize * row ] );
386  }
387  }
388 
389  jump += ( linsize * ( linsize - 1 ) ) / 2;
390 
391  }
392  }
393 
394  if ( true ){
396  tempU.updateUnitary( temp1, temp2, vector, false, false ); //multiply = compact = false
397  const double rms_diff = rms_deviation( &tempU );
398  cout << " DMRGSCFunitary::getLog : 2-norm of [ U - exp(ln(U)) ] (should be 0.0) = " << rms_diff << endl;
399  }
400 
401 }
402 
403 void CheMPS2::DMRGSCFunitary::makeSureAllBlocksDetOne( double * temp1, double * temp2 ){
404 
405  for ( int irrep = 0; irrep < num_irreps; irrep++ ){
406 
407  int linsize = iHandler->getNORB( irrep );
408  int size = linsize * linsize;
409 
410  if ( linsize > 1 ){
411  /* linsize >= 2; hence temp1 is at least of size 4*linsize*linsize
412  if linsize <= 1; there corresponds no block in xmatrix to it */
413 
414  double * work1 = temp1; //linsize * linsize
415  double * work2 = temp1 + size; //linsize * linsize
416 
417  double * workLARGE = temp2;
418  int lworkLARGE = 4 * size;
419 
420  const double determinant = get_determinant( irrep, work1, work2, workLARGE, lworkLARGE );
421  if ( determinant < 0.0 ){ // determinant = -1
422  // U <-- diag(-1, 1, 1, 1, ...) * U : First row of U changes sign!
423  for ( int cnt = 0; cnt < linsize; cnt++ ){ entries[ irrep ][ 0 + linsize * cnt ] *= -1; }
424  }
425 
426  }
427  }
428 
429 }
430 
431 void CheMPS2::DMRGSCFunitary::saveU( const string filename ) const{
432 
434 
435 }
436 
437 void CheMPS2::DMRGSCFunitary::loadU( const string filename ){
438 
440 
441 }
442 
443 
int getNORB(const int irrep) const
Get the number of orbitals for an irrep.
void rotateActiveSpaceVectors(double *eigenvecs, double *work)
Rotate the unitary matrix.
const DMRGSCFindices * iHandler
The information on the occupied, active, and virtual spaces.
Definition: DMRGSCFmatrix.h:94
DMRGSCFunitary(const DMRGSCFindices *iHandler)
Constructor.
void loadU(const string filename=DMRGSCF_unitary_storage_name)
Load the unitary from disk.
static void read(const string filename, const int n_irreps, double **storage)
Read the DMRGSCFmatrix from disk.
int getDMRGcumulative(const int irrep) const
Get the cumulative number of active orbitals for an irrep.
void getLog(double *vector, double *temp1, double *temp2) const
Obtain the logarithm of the unitary matrix.
void saveU(const string filename=DMRGSCF_unitary_storage_name) const
Save the unitary to disk.
virtual ~DMRGSCFunitary()
Destructor.
double ** entries
The matrix entries.
Definition: DMRGSCFmatrix.h:97
static void write(const string filename, const DMRGSCFindices *idx, double **storage)
Write a DMRGSCFmatrix to disk.
double rms_deviation(const DMRGSCFmatrix *other) const
Get the RMS deviation with another DMRGSCFmatrix.
void CheckDeviationFromUnitary(double *work) const
Calculate the two-norm of U^T*U - I.
void makeSureAllBlocksDetOne(double *temp1, double *temp2)
Orbitals are defined up to a phase factor. Make sure that the logarithm of each block of the unitary ...
int getNumVariablesX() const
Get the number of variables in the x-parametrization of the unitary update.
int num_irreps
The number of irreps.
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 getNDMRG(const int irrep) const
Get the number of active orbitals for an irrep.
int getNOCC(const int irrep) const
Get the number of occupied orbitals for an irrep.
int getNVIRT(const int irrep) const
Get the number of virtual orbitals for an irrep.