CheMPS2
Davidson.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 <iostream>
23 #include <assert.h>
24 
25 #include "Davidson.h"
26 #include "Lapack.h"
27 
28 using std::cout;
29 using std::endl;
30 
31 CheMPS2::Davidson::Davidson( const int veclength, const int MAX_NUM_VEC, const int NUM_VEC_KEEP, const double RTOL, const double DIAG_CUTOFF, const bool debug_print, const char problem_type ){
32 
33  assert( ( problem_type == 'E' ) || ( problem_type == 'L' ) );
34 
35  this->debug_print = debug_print;
36  this->veclength = veclength;
37  this->problem_type = problem_type;
38  this->MAX_NUM_VEC = MAX_NUM_VEC;
39  this->NUM_VEC_KEEP = NUM_VEC_KEEP;
40  this->DIAG_CUTOFF = DIAG_CUTOFF;
41  this->RTOL = RTOL;
42 
43  state = 'I'; // <I>nitialized Davidson
44  nMultiplications = 0;
45 
46  // To store the vectors and the matrix x vectors
47  num_vec = 0;
48  vecs = new double*[ MAX_NUM_VEC ];
49  Hvecs = new double*[ MAX_NUM_VEC ];
50  num_allocated = 0;
51 
52  // The projected problem
53  mxM = new double[ MAX_NUM_VEC * MAX_NUM_VEC ];
54  mxM_eigs = new double[ MAX_NUM_VEC ];
55  mxM_vecs = new double[ MAX_NUM_VEC * MAX_NUM_VEC ];
56  mxM_lwork = 3 * MAX_NUM_VEC - 1;
57  mxM_work = new double[ mxM_lwork ];
58  mxM_rhs = (( problem_type == 'L' ) ? new double[ MAX_NUM_VEC ] : NULL );
59 
60  // Vector spaces
61  diag = new double[ veclength ];
62  t_vec = new double[ veclength ];
63  u_vec = new double[ veclength ];
64  work_vec = new double[ veclength ];
65  RHS = (( problem_type == 'L' ) ? new double[ veclength ] : NULL );
66 
67  // For the deflation
68  Reortho_Lowdin = NULL;
69  Reortho_Overlap_eigs = NULL;
70  Reortho_Overlap = NULL;
71  Reortho_Eigenvecs = NULL;
72 
73 }
74 
76 
77  for (int cnt = 0; cnt < num_allocated; cnt++){
78  delete [] vecs[cnt];
79  delete [] Hvecs[cnt];
80  }
81  delete [] vecs;
82  delete [] Hvecs;
83 
84  delete [] mxM;
85  delete [] mxM_eigs;
86  delete [] mxM_vecs;
87  delete [] mxM_work;
88  if ( mxM_rhs != NULL ){ delete [] mxM_rhs; }
89 
90  delete [] diag;
91  delete [] t_vec;
92  delete [] u_vec;
93  delete [] work_vec;
94  if ( RHS != NULL ){ delete [] RHS; }
95 
96  if ( Reortho_Lowdin != NULL ){ delete [] Reortho_Lowdin; }
97  if ( Reortho_Overlap_eigs != NULL ){ delete [] Reortho_Overlap_eigs; }
98  if ( Reortho_Overlap != NULL ){ delete [] Reortho_Overlap; }
99  if ( Reortho_Eigenvecs != NULL ){ delete [] Reortho_Eigenvecs; }
100 
101 }
102 
103 int CheMPS2::Davidson::GetNumMultiplications() const{ return nMultiplications; }
104 
105 char CheMPS2::Davidson::FetchInstruction( double ** pointers ){
106 
107  /*
108  Possible states:
109  - I : just initialized
110  - U : just before the big loop, the initial guess and the diagonal are set
111  - N : a new vector has just been added to the list and a matrix-vector multiplication has been performed
112  - F : the space has been deflated and a few matrix-vector multiplications are required
113  - C : convergence was reached
114 
115  Possible instructions:
116  - A : copy the initial guess to pointers[0], the diagonal to pointers[1], and if problem_type=='L' the right-hand side of the linear problem to pointers[2]
117  - B : perform pointers[1] = symmetric matrix times pointers[0]
118  - C : copy the converged solution from pointers[0] back; pointers[1][0] contains the converged energy if problem_type=='E' and the residual norm if problem_type=='L'
119  - D : there was an error
120  */
121 
122  if ( state == 'I' ){
123  pointers[ 0 ] = t_vec;
124  pointers[ 1 ] = diag;
125  if ( problem_type == 'L' ){ pointers[ 2 ] = RHS; }
126  state = 'U';
127  return 'A';
128  }
129 
130  if ( state == 'U' ){
131  SafetyCheckGuess();
132  AddNewVec();
133  pointers[ 0 ] = vecs[ num_vec ];
134  pointers[ 1 ] = Hvecs[ num_vec ];
135  nMultiplications++;
136  state = 'N';
137  return 'B';
138  }
139 
140  if ( state == 'N' ){
141  const double rnorm = DiagonalizeSmallMatrixAndCalcResidual();
142  // if ( debug_print ){ cout << "WARNING AT DAVIDSON : Current residual norm = " << rnorm << endl; }
143  if ( rnorm > RTOL ){ // Not yet converged
144  CalculateNewVec();
145  if ( num_vec == MAX_NUM_VEC ){
146  Deflation();
147  pointers[ 0 ] = vecs[ num_vec ];
148  pointers[ 1 ] = Hvecs[ num_vec ];
149  nMultiplications++;
150  num_vec++;
151  state = 'F';
152  return 'B';
153  }
154  AddNewVec();
155  pointers[ 0 ] = vecs[ num_vec ];
156  pointers[ 1 ] = Hvecs[ num_vec ];
157  nMultiplications++;
158  state = 'N';
159  return 'B';
160  } else { // Converged
161  state = 'C';
162  pointers[ 0 ] = u_vec;
163  pointers[ 1 ] = work_vec;
164  if ( problem_type == 'E' ){ work_vec[ 0 ] = mxM_eigs[ 0 ]; }
165  if ( problem_type == 'L' ){ work_vec[ 0 ] = rnorm; }
166  return 'C';
167  }
168  }
169 
170  if ( state == 'F' ){
171  if ( num_vec == NUM_VEC_KEEP ){
172  MxMafterDeflation();
173  AddNewVec();
174  pointers[ 0 ] = vecs[ num_vec ];
175  pointers[ 1 ] = Hvecs[ num_vec ];
176  nMultiplications++;
177  state = 'N';
178  return 'B';
179  } else {
180  pointers[ 0 ] = vecs[ num_vec ];
181  pointers[ 1 ] = Hvecs[ num_vec ];
182  nMultiplications++;
183  num_vec++;
184  state = 'F';
185  return 'B';
186  }
187  }
188 
189  return 'D';
190 
191 }
192 
193 double CheMPS2::Davidson::FrobeniusNorm( double * current_vector ){
194 
195  char frobenius = 'F';
196  int inc1 = 1;
197  const double twonorm = dlange_( &frobenius, &veclength, &inc1, current_vector, &veclength, NULL ); // Work is not referenced for Frobenius norm
198  return twonorm;
199 
200 }
201 
202 void CheMPS2::Davidson::SafetyCheckGuess(){
203 
204  const double twonorm = FrobeniusNorm( t_vec );
205  if ( twonorm == 0.0 ){
206  for ( int cnt = 0; cnt < veclength; cnt++ ){ t_vec[ cnt ] = ( (double) rand() ) / RAND_MAX; }
207  if ( debug_print ){
208  cout << "WARNING AT DAVIDSON : Initial guess was a zero-vector. Now it is overwritten with random numbers." << endl;
209  }
210  }
211 
212 }
213 
214 void CheMPS2::Davidson::AddNewVec(){
215 
216  int inc1 = 1;
217 
218  // Orthogonalize the new vector w.r.t. the old basis
219  for ( int cnt = 0; cnt < num_vec; cnt++ ){
220  double minus_overlap = - ddot_( &veclength, t_vec, &inc1, vecs[ cnt ], &inc1 );
221  daxpy_( &veclength, &minus_overlap, vecs[ cnt ], &inc1, t_vec, &inc1 );
222  }
223 
224  // Normalize the new vector
225  double alpha = 1.0 / FrobeniusNorm( t_vec );
226  dscal_( &veclength, &alpha, t_vec, &inc1 );
227 
228  // The new vector becomes part of vecs
229  if ( num_vec < num_allocated ){
230  double * temp = vecs[ num_vec ];
231  vecs[ num_vec ] = t_vec;
232  t_vec = temp;
233  } else {
234  vecs[ num_allocated ] = t_vec;
235  Hvecs[ num_allocated ] = new double[ veclength ];
236  t_vec = new double[ veclength ];
237  num_allocated++;
238  }
239 
240 }
241 
242 double CheMPS2::Davidson::DiagonalizeSmallMatrixAndCalcResidual(){
243 
244  int inc1 = 1;
245 
246  if ( problem_type == 'E' ){ // EIGENVALUE PROBLEM
247  // mxM contains V^T . A . V
248  for ( int cnt = 0; cnt < num_vec; cnt++ ){
249  mxM[ cnt + MAX_NUM_VEC * num_vec ] = ddot_( &veclength, vecs[ num_vec ], &inc1, Hvecs[ cnt ], &inc1 );
250  mxM[ num_vec + MAX_NUM_VEC * cnt ] = mxM[ cnt + MAX_NUM_VEC * num_vec ];
251  }
252  mxM[ num_vec + MAX_NUM_VEC * num_vec ] = ddot_( &veclength, vecs[ num_vec ], &inc1, Hvecs[ num_vec ], &inc1 );
253  } else { // LINEAR PROBLEM
254  // mxM contains V^T . A^T . A . V
255  for ( int cnt = 0; cnt < num_vec; cnt++ ){
256  mxM[ cnt + MAX_NUM_VEC * num_vec ] = ddot_( &veclength, Hvecs[ num_vec ], &inc1, Hvecs[ cnt ], &inc1 );
257  mxM[ num_vec + MAX_NUM_VEC * cnt ] = mxM[ cnt + MAX_NUM_VEC * num_vec ];
258  }
259  mxM[ num_vec + MAX_NUM_VEC * num_vec ] = ddot_( &veclength, Hvecs[ num_vec ], &inc1, Hvecs[ num_vec ], &inc1 );
260  // mxM_rhs contains V^T . A^T . RHS
261  mxM_rhs[ num_vec ] = ddot_( &veclength, Hvecs[ num_vec ], &inc1, RHS, &inc1 );
262  }
263 
264  // When t-vec was added to vecs, the number of vecs was actually increased by one. Now the number is incremented.
265  num_vec++;
266 
267  // Diagonalize mxM ( always )
268  char jobz = 'V';
269  char uplo = 'U';
270  int info;
271  for ( int cnt1 = 0; cnt1 < num_vec; cnt1++ ){
272  for ( int cnt2 = 0; cnt2 < num_vec; cnt2++ ){
273  mxM_vecs[ cnt1 + MAX_NUM_VEC * cnt2 ] = mxM[ cnt1 + MAX_NUM_VEC * cnt2 ];
274  }
275  }
276  dsyev_( &jobz, &uplo, &num_vec, mxM_vecs, &MAX_NUM_VEC, mxM_eigs, mxM_work, &mxM_lwork, &info ); // Ascending order of eigenvalues
277 
278  // If problem_type == linear, solve the linear problem in least-squares sense!
279  if ( problem_type == 'L' ){
280  double one = 1.0;
281  double set = 0.0;
282  char trans = 'T';
283  char notra = 'N';
284  dgemm_( &trans, &notra, &num_vec, &inc1, &num_vec, &one, mxM_vecs, &MAX_NUM_VEC, mxM_rhs, &MAX_NUM_VEC, &set, mxM_work, &MAX_NUM_VEC ); // mxM_work = mxM_vecs^T * mxM_rhs
285  for ( int cnt = 0; cnt < num_vec; cnt++ ){
286  double current_eigenvalue = mxM_eigs[ cnt ];
287  if ( fabs( current_eigenvalue ) < DIAG_CUTOFF ){
288  current_eigenvalue = DIAG_CUTOFF * (( current_eigenvalue < 0.0 ) ? -1 : 1 );
289  if ( debug_print ){
290  cout << "WARNING AT DAVIDSON : The eigenvalue " << mxM_eigs[ cnt ] << " to solve Ax = b has been overwritten with " << current_eigenvalue << "." << endl;
291  }
292  }
293  mxM_work[ cnt ] = mxM_work[ cnt ] / current_eigenvalue;
294  }
295  dgemm_( &notra, &notra, &num_vec, &inc1, &num_vec, &one, mxM_vecs, &MAX_NUM_VEC, mxM_work, &MAX_NUM_VEC, &set, mxM_work + MAX_NUM_VEC, &MAX_NUM_VEC );
296  for ( int cnt = 0; cnt < num_vec; cnt++ ){ mxM_vecs[ cnt ] = mxM_work[ MAX_NUM_VEC + cnt ]; } // mxM_vecs = U * eigs^{-1} * U^T * RHS
297  }
298 
299  // Calculate u and r. r is stored in t_vec, u in u_vec.
300  for ( int cnt = 0; cnt < veclength; cnt++ ){ t_vec[ cnt ] = 0.0; }
301  for ( int cnt = 0; cnt < veclength; cnt++ ){ u_vec[ cnt ] = 0.0; }
302  for ( int cnt = 0; cnt < num_vec; cnt++ ){
303  double alpha = mxM_vecs[ cnt ]; // Eigenvector with lowest eigenvalue, hence mxM_vecs[ cnt + MAX_NUM_VEC * 0 ]
304  daxpy_( &veclength, &alpha, Hvecs[ cnt ], &inc1, t_vec, &inc1 );
305  daxpy_( &veclength, &alpha, vecs[ cnt ], &inc1, u_vec, &inc1 );
306  }
307  if ( problem_type == 'E' ){ // t_vec = H * x - lambda * x
308  double alpha = - mxM_eigs[ 0 ];
309  daxpy_( &veclength, &alpha, u_vec, &inc1, t_vec, &inc1 );
310  } else { // t_vec = H * x - RHS
311  double alpha = -1.0;
312  daxpy_( &veclength, &alpha, RHS, &inc1, t_vec, &inc1 );
313  }
314 
315  // Calculate the norm of r
316  const double rnorm = FrobeniusNorm( t_vec );
317  //cout << " Davidson :: rnorm = " << rnorm << endl;
318  return rnorm;
319 
320 }
321 
322 void CheMPS2::Davidson::CalculateNewVec(){
323 
324  int inc1 = 1;
325  const double shift = (( problem_type == 'E' ) ? mxM_eigs[ 0 ] : 0.0 );
326 
327  // Calculate the new t_vec based on the residual of the lowest eigenvalue, to add to the vecs.
328  for ( int cnt = 0; cnt < veclength; cnt++ ){
329  const double difference = diag[ cnt ] - shift;
330  const double fabsdiff = fabs( difference );
331  if ( fabsdiff > DIAG_CUTOFF ){
332  work_vec[ cnt ] = u_vec[ cnt ] / difference; // work_vec = K^(-1) u_vec
333  } else {
334  work_vec[ cnt ] = u_vec[ cnt ] / DIAG_CUTOFF;
335  if ( debug_print ){ cout << "WARNING AT DAVIDSON : fabs( precon[" << cnt << "] ) = " << fabsdiff << endl; }
336  }
337  }
338  double alpha = - ddot_( &veclength, work_vec, &inc1, t_vec, &inc1 ) / ddot_( &veclength, work_vec, &inc1, u_vec, &inc1 ); // alpha = - (u^T K^(-1) r) / (u^T K^(-1) u)
339  daxpy_( &veclength, &alpha, u_vec, &inc1, t_vec, &inc1 ); // t_vec = r - (u^T K^(-1) r) / (u^T K^(-1) u) u
340  for ( int cnt = 0; cnt < veclength; cnt++ ){
341  const double difference = diag[ cnt ] - shift;
342  const double fabsdiff = fabs( difference );
343  if ( fabsdiff > DIAG_CUTOFF ){
344  t_vec[ cnt ] = - t_vec[ cnt ] / difference; // t_vec = - K^(-1) (r - (u^T K^(-1) r) / (u^T K^(-1) u) u)
345  } else {
346  t_vec[ cnt ] = - t_vec[ cnt ] / DIAG_CUTOFF;
347  }
348  }
349 
350 }
351 
352 void CheMPS2::Davidson::Deflation(){
353 
354  int inc1 = 1;
355 
356  // When the maximum number of vectors is reached: construct the one with lowest eigenvalue & restart
357  if ( NUM_VEC_KEEP <= 1 ){
358 
359  double alpha = 1.0 / FrobeniusNorm( u_vec );
360  dscal_( &veclength, &alpha, u_vec, &inc1 );
361  dcopy_( &veclength, u_vec, &inc1, vecs[ 0 ], &inc1 );
362 
363  } else {
364 
365  if ( problem_type == 'L' ){ SolveLinearSystemDeflation( NUM_VEC_KEEP ); }
366 
367  if ( Reortho_Eigenvecs == NULL ){ Reortho_Eigenvecs = new double[ veclength * NUM_VEC_KEEP ]; }
368  if ( Reortho_Overlap == NULL ){ Reortho_Overlap = new double[ NUM_VEC_KEEP * NUM_VEC_KEEP ]; }
369  if ( Reortho_Overlap_eigs == NULL ){ Reortho_Overlap_eigs = new double[ NUM_VEC_KEEP ]; }
370  if ( Reortho_Lowdin == NULL ){ Reortho_Lowdin = new double[ NUM_VEC_KEEP * NUM_VEC_KEEP ]; }
371 
372  // Construct the lowest NUM_VEC_KEEP eigenvectors
373  dcopy_( &veclength, u_vec, &inc1, Reortho_Eigenvecs, &inc1 );
374  for ( int cnt = 1; cnt < NUM_VEC_KEEP; cnt++ ){
375  for ( int irow = 0; irow < veclength; irow++ ){
376  Reortho_Eigenvecs[ irow + veclength * cnt ] = 0.0;
377  for ( int ivec = 0; ivec < MAX_NUM_VEC; ivec++ ){
378  Reortho_Eigenvecs[ irow + veclength * cnt ] += vecs[ ivec ][ irow ] * mxM_vecs[ ivec + MAX_NUM_VEC * cnt ];
379  }
380  }
381  }
382 
383  // Calculate the overlap matrix
384  char trans = 'T';
385  char notr = 'N';
386  double one = 1.0;
387  double zero = 0.0; //set
388  dgemm_( &trans, &notr, &NUM_VEC_KEEP, &NUM_VEC_KEEP, &veclength, &one, Reortho_Eigenvecs, &veclength, Reortho_Eigenvecs, &veclength, &zero, Reortho_Overlap, &NUM_VEC_KEEP );
389 
390  // Calculate the Lowdin tfo
391  char jobz = 'V';
392  char uplo = 'U';
393  int info;
394  dsyev_( &jobz, &uplo, &NUM_VEC_KEEP, Reortho_Overlap, &NUM_VEC_KEEP, Reortho_Overlap_eigs, mxM_work, &mxM_lwork, &info ); // Ascending order of eigs
395  for ( int icnt = 0; icnt < NUM_VEC_KEEP; icnt++ ){
396  Reortho_Overlap_eigs[ icnt ] = pow( Reortho_Overlap_eigs[ icnt ], -0.25 );
397  dscal_( &NUM_VEC_KEEP, Reortho_Overlap_eigs + icnt, Reortho_Overlap + NUM_VEC_KEEP * icnt, &inc1 );
398  }
399  dgemm_( &notr, &trans, &NUM_VEC_KEEP, &NUM_VEC_KEEP, &NUM_VEC_KEEP, &one, Reortho_Overlap, &NUM_VEC_KEEP, Reortho_Overlap, &NUM_VEC_KEEP, &zero, Reortho_Lowdin, &NUM_VEC_KEEP );
400 
401  // Reortho: Put the Lowdin tfo eigenvecs in vecs
402  for ( int ivec = 0; ivec < NUM_VEC_KEEP; ivec++ ){
403  for ( int loop = 0; loop < veclength; loop++ ){ vecs[ ivec ][ loop ] = 0.0; }
404  for ( int ivec2 = 0; ivec2 < NUM_VEC_KEEP; ivec2++ ){
405  daxpy_( &veclength, Reortho_Lowdin + ivec2 + NUM_VEC_KEEP * ivec, Reortho_Eigenvecs + veclength * ivec2, &inc1, vecs[ ivec ], &inc1 );
406  }
407  }
408  }
409 
410  num_vec = 0;
411 
412 }
413 
414 void CheMPS2::Davidson::SolveLinearSystemDeflation( const int NUM_SOLUTIONS ){
415 
416  assert( problem_type == 'L' );
417  assert( num_vec == MAX_NUM_VEC );
418  assert( NUM_SOLUTIONS <= MAX_NUM_VEC );
419  assert( 2 <= NUM_SOLUTIONS );
420 
421  double * work1 = new double[ MAX_NUM_VEC * MAX_NUM_VEC ]; // projector
422  double * work3 = new double[ MAX_NUM_VEC * MAX_NUM_VEC ]; // projector times mxM_rhs
423  double * work2 = new double[ MAX_NUM_VEC * NUM_SOLUTIONS ]; // solutions
424 
425  for ( int solution = 0; solution < NUM_SOLUTIONS; solution++ ){
426 
427  // work1 = ( 1 - sum_j v_j v_j^T ) = projector
428  for ( int cntr = 0; cntr < MAX_NUM_VEC * MAX_NUM_VEC; cntr++ ){ work1[ cntr ] = 0.0; }
429  for ( int diag = 0; diag < MAX_NUM_VEC; diag++ ){ work1[ diag * ( 1 + MAX_NUM_VEC ) ] = 1.0; }
430  for ( int prev = 0; prev < solution; prev++ ){
431  for ( int col = 0; col < MAX_NUM_VEC; col++ ){
432  for ( int row = 0; row < MAX_NUM_VEC; row++ ){
433  work1[ row + MAX_NUM_VEC * col ] -= work2[ row + MAX_NUM_VEC * prev ] * work2[ col + MAX_NUM_VEC * prev ];
434  }
435  }
436  }
437 
438  // work3 = ( 1 - sum_j v_j v_j^T ) * [ U^T * A^T * b ] = work1 * mxM_rhs
439  // mxM_vecs = ( 1 - sum_j v_j v_j^T ) * [ U^T * A^T * A * U ] * ( 1 - sum_j v_j v_j^T ) = work1 * mxM * work1
440  {
441  double one = 1.0;
442  double set = 0.0;
443  char notrans = 'N';
444  int inc1 = 1;
445  dgemm_( &notrans, &notrans, &MAX_NUM_VEC, &MAX_NUM_VEC, &MAX_NUM_VEC, &one, work1, &MAX_NUM_VEC, mxM, &MAX_NUM_VEC, &set, work3, &MAX_NUM_VEC );
446  dgemm_( &notrans, &notrans, &MAX_NUM_VEC, &MAX_NUM_VEC, &MAX_NUM_VEC, &one, work3, &MAX_NUM_VEC, work1, &MAX_NUM_VEC, &set, mxM_vecs, &MAX_NUM_VEC );
447  dgemm_( &notrans, &notrans, &MAX_NUM_VEC, &inc1, &MAX_NUM_VEC, &one, work1, &MAX_NUM_VEC, mxM_rhs, &MAX_NUM_VEC, &set, work3, &MAX_NUM_VEC );
448  }
449 
450  // Diagonalize mxM_vecs = V * lambda * V^T ===> Ascending order of eigenvalues & ( V, lambda ) = ( mxM_vecs, mxM_eigs )
451  {
452  char jobz = 'V';
453  char uplo = 'U';
454  int info;
455  dsyev_( &jobz, &uplo, &MAX_NUM_VEC, mxM_vecs, &MAX_NUM_VEC, mxM_eigs, mxM_work, &mxM_lwork, &info );
456  }
457 
458  // work2[ :, solution ] = [ ( 1 - sum_j v_j v_j^T ) * [ U^T * A^T * A * U ] * ( 1 - sum_j v_j v_j^T ) ]^{-1} * ( 1 - sum_j v_j v_j^T ) * [ U^T * A^T * b ] = V * lambda^{-1} * V^T * work3
459  {
460  double one = 1.0;
461  double set = 0.0;
462  char trans = 'T';
463  char notrans = 'N';
464  int inc1 = 1;
465  dgemm_( &trans, &notrans, &MAX_NUM_VEC, &inc1, &MAX_NUM_VEC, &one, mxM_vecs, &MAX_NUM_VEC, work3, &MAX_NUM_VEC, &set, mxM_work, &MAX_NUM_VEC );
466  for ( int diag = 0; diag < MAX_NUM_VEC; diag++ ){
467  if ( diag < solution ){
468  mxM_work[ diag ] = 0.0; // PSEUDOINVERSE
469  } else {
470  double current_eigenvalue = mxM_eigs[ diag ];
471  if ( fabs( current_eigenvalue ) < DIAG_CUTOFF ){
472  current_eigenvalue = DIAG_CUTOFF * (( current_eigenvalue < 0.0 ) ? -1 : 1 );
473  if ( debug_print ){
474  cout << "WARNING AT DAVIDSON : The eigenvalue " << mxM_eigs[ diag ] << " to solve Ax = b has been overwritten with " << current_eigenvalue << "." << endl;
475  }
476  }
477  mxM_work[ diag ] = mxM_work[ diag ] / current_eigenvalue;
478  }
479  }
480  dgemm_( &notrans, &notrans, &MAX_NUM_VEC, &inc1, &MAX_NUM_VEC, &one, mxM_vecs, &MAX_NUM_VEC, mxM_work, &MAX_NUM_VEC, &set, work2 + MAX_NUM_VEC * solution, &MAX_NUM_VEC );
481  }
482 
483  // Normalize work2[ :, solution ]
484  {
485  int inc1 = 1;
486  double * ptr = work2 + MAX_NUM_VEC * solution;
487  const double twonorm = sqrt( ddot_( &MAX_NUM_VEC, ptr, &inc1, ptr, &inc1 ) );
488  /*if ( debug_print ){
489  cout << "Davidson :: Deflation :: Norm of solution " << solution << " = " << twonorm << endl;
490  }*/
491  double factor = 1.0 / twonorm;
492  dscal_( &MAX_NUM_VEC, &factor, ptr, &inc1 );
493  }
494 
495  }
496 
497  // Copy over work2 to mxM_vecs
498  {
499  int inc1 = 1;
500  int size = MAX_NUM_VEC * NUM_SOLUTIONS;
501  dcopy_( &size, work2, &inc1, mxM_vecs, &inc1 );
502  }
503 
504  delete [] work1;
505  delete [] work2;
506  delete [] work3;
507 
508 }
509 
510 void CheMPS2::Davidson::MxMafterDeflation(){
511 
512  int inc1 = 1;
513 
514  if ( problem_type == 'E' ){ // EIGENVALUE PROBLEM
515  // mxM contains V^T . A . V
516  for ( int ivec = 0; ivec < NUM_VEC_KEEP; ivec++ ){
517  for ( int ivec2 = ivec; ivec2 < NUM_VEC_KEEP; ivec2++ ){
518  mxM[ ivec + MAX_NUM_VEC * ivec2 ] = ddot_( &veclength, vecs[ ivec ], &inc1, Hvecs[ ivec2 ], &inc1 );
519  mxM[ ivec2 + MAX_NUM_VEC * ivec ] = mxM[ ivec + MAX_NUM_VEC * ivec2 ];
520  }
521  }
522  } else { // LINEAR PROBLEM
523  // mxM contains V^T . A^T . A . V
524  for ( int ivec = 0; ivec < NUM_VEC_KEEP; ivec++ ){
525  for ( int ivec2 = ivec; ivec2 < NUM_VEC_KEEP; ivec2++ ){
526  mxM[ ivec + MAX_NUM_VEC * ivec2 ] = ddot_( &veclength, Hvecs[ ivec ], &inc1, Hvecs[ ivec2 ], &inc1 );
527  mxM[ ivec2 + MAX_NUM_VEC * ivec ] = mxM[ ivec + MAX_NUM_VEC * ivec2 ];
528  }
529  }
530  // mxM_rhs contains V^T . A^T . RHS
531  for ( int ivec = 0; ivec < NUM_VEC_KEEP; ivec++ ){
532  mxM_rhs[ ivec ] = ddot_( &veclength, Hvecs[ ivec ], &inc1, RHS, &inc1 );
533  }
534  }
535 
536 }
537 
538 
virtual ~Davidson()
Destructor.
Definition: Davidson.cpp:75
int GetNumMultiplications() const
Get the number of matrix vector multiplications which have been performed.
Definition: Davidson.cpp:103
Davidson(const int veclength, const int MAX_NUM_VEC, const int NUM_VEC_KEEP, const double RTOL, const double DIAG_CUTOFF, const bool debug_print, const char problem_type= 'E')
Constructor.
Definition: Davidson.cpp:31
char FetchInstruction(double **pointers)
The iterator to converge the ground state vector.
Definition: Davidson.cpp:105