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 ){
33 assert( ( problem_type ==
'E' ) || ( problem_type ==
'L' ) );
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;
48 vecs =
new double*[ MAX_NUM_VEC ];
49 Hvecs =
new double*[ MAX_NUM_VEC ];
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 );
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 );
68 Reortho_Lowdin = NULL;
69 Reortho_Overlap_eigs = NULL;
70 Reortho_Overlap = NULL;
71 Reortho_Eigenvecs = NULL;
77 for (
int cnt = 0; cnt < num_allocated; cnt++){
88 if ( mxM_rhs != NULL ){
delete [] mxM_rhs; }
94 if ( RHS != NULL ){
delete [] RHS; }
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; }
123 pointers[ 0 ] = t_vec;
124 pointers[ 1 ] = diag;
125 if ( problem_type ==
'L' ){ pointers[ 2 ] = RHS; }
133 pointers[ 0 ] = vecs[ num_vec ];
134 pointers[ 1 ] = Hvecs[ num_vec ];
141 const double rnorm = DiagonalizeSmallMatrixAndCalcResidual();
145 if ( num_vec == MAX_NUM_VEC ){
147 pointers[ 0 ] = vecs[ num_vec ];
148 pointers[ 1 ] = Hvecs[ num_vec ];
155 pointers[ 0 ] = vecs[ num_vec ];
156 pointers[ 1 ] = Hvecs[ num_vec ];
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; }
171 if ( num_vec == NUM_VEC_KEEP ){
174 pointers[ 0 ] = vecs[ num_vec ];
175 pointers[ 1 ] = Hvecs[ num_vec ];
180 pointers[ 0 ] = vecs[ num_vec ];
181 pointers[ 1 ] = Hvecs[ num_vec ];
193 double CheMPS2::Davidson::FrobeniusNorm(
double * current_vector ){
195 char frobenius =
'F';
197 const double twonorm = dlange_( &frobenius, &veclength, &inc1, current_vector, &veclength, NULL );
202 void CheMPS2::Davidson::SafetyCheckGuess(){
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; }
208 cout <<
"WARNING AT DAVIDSON : Initial guess was a zero-vector. Now it is overwritten with random numbers." << endl;
214 void CheMPS2::Davidson::AddNewVec(){
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 );
225 double alpha = 1.0 / FrobeniusNorm( t_vec );
226 dscal_( &veclength, &alpha, t_vec, &inc1 );
229 if ( num_vec < num_allocated ){
230 double * temp = vecs[ num_vec ];
231 vecs[ num_vec ] = t_vec;
234 vecs[ num_allocated ] = t_vec;
235 Hvecs[ num_allocated ] =
new double[ veclength ];
236 t_vec =
new double[ veclength ];
242 double CheMPS2::Davidson::DiagonalizeSmallMatrixAndCalcResidual(){
246 if ( problem_type ==
'E' ){
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 ];
252 mxM[ num_vec + MAX_NUM_VEC * num_vec ] = ddot_( &veclength, vecs[ num_vec ], &inc1, Hvecs[ num_vec ], &inc1 );
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 ];
259 mxM[ num_vec + MAX_NUM_VEC * num_vec ] = ddot_( &veclength, Hvecs[ num_vec ], &inc1, Hvecs[ num_vec ], &inc1 );
261 mxM_rhs[ num_vec ] = ddot_( &veclength, Hvecs[ num_vec ], &inc1, RHS, &inc1 );
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 ];
276 dsyev_( &jobz, &uplo, &num_vec, mxM_vecs, &MAX_NUM_VEC, mxM_eigs, mxM_work, &mxM_lwork, &info );
279 if ( problem_type ==
'L' ){
284 dgemm_( &trans, ¬ra, &num_vec, &inc1, &num_vec, &one, mxM_vecs, &MAX_NUM_VEC, mxM_rhs, &MAX_NUM_VEC, &
set, mxM_work, &MAX_NUM_VEC );
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 );
290 cout <<
"WARNING AT DAVIDSON : The eigenvalue " << mxM_eigs[ cnt ] <<
" to solve Ax = b has been overwritten with " << current_eigenvalue <<
"." << endl;
293 mxM_work[ cnt ] = mxM_work[ cnt ] / current_eigenvalue;
295 dgemm_( ¬ra, ¬ra, &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 ]; }
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 ];
304 daxpy_( &veclength, &alpha, Hvecs[ cnt ], &inc1, t_vec, &inc1 );
305 daxpy_( &veclength, &alpha, vecs[ cnt ], &inc1, u_vec, &inc1 );
307 if ( problem_type ==
'E' ){
308 double alpha = - mxM_eigs[ 0 ];
309 daxpy_( &veclength, &alpha, u_vec, &inc1, t_vec, &inc1 );
312 daxpy_( &veclength, &alpha, RHS, &inc1, t_vec, &inc1 );
316 const double rnorm = FrobeniusNorm( t_vec );
322 void CheMPS2::Davidson::CalculateNewVec(){
325 const double shift = (( problem_type ==
'E' ) ? mxM_eigs[ 0 ] : 0.0 );
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;
334 work_vec[ cnt ] = u_vec[ cnt ] / DIAG_CUTOFF;
335 if ( debug_print ){ cout <<
"WARNING AT DAVIDSON : fabs( precon[" << cnt <<
"] ) = " << fabsdiff << endl; }
338 double alpha = - ddot_( &veclength, work_vec, &inc1, t_vec, &inc1 ) / ddot_( &veclength, work_vec, &inc1, u_vec, &inc1 );
339 daxpy_( &veclength, &alpha, u_vec, &inc1, t_vec, &inc1 );
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;
346 t_vec[ cnt ] = - t_vec[ cnt ] / DIAG_CUTOFF;
352 void CheMPS2::Davidson::Deflation(){
357 if ( NUM_VEC_KEEP <= 1 ){
359 double alpha = 1.0 / FrobeniusNorm( u_vec );
360 dscal_( &veclength, &alpha, u_vec, &inc1 );
361 dcopy_( &veclength, u_vec, &inc1, vecs[ 0 ], &inc1 );
365 if ( problem_type ==
'L' ){ SolveLinearSystemDeflation( NUM_VEC_KEEP ); }
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 ]; }
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 ];
388 dgemm_( &trans, ¬r, &NUM_VEC_KEEP, &NUM_VEC_KEEP, &veclength, &one, Reortho_Eigenvecs, &veclength, Reortho_Eigenvecs, &veclength, &zero, Reortho_Overlap, &NUM_VEC_KEEP );
394 dsyev_( &jobz, &uplo, &NUM_VEC_KEEP, Reortho_Overlap, &NUM_VEC_KEEP, Reortho_Overlap_eigs, mxM_work, &mxM_lwork, &info );
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 );
399 dgemm_( ¬r, &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 );
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 );
414 void CheMPS2::Davidson::SolveLinearSystemDeflation(
const int NUM_SOLUTIONS ){
416 assert( problem_type ==
'L' );
417 assert( num_vec == MAX_NUM_VEC );
418 assert( NUM_SOLUTIONS <= MAX_NUM_VEC );
419 assert( 2 <= NUM_SOLUTIONS );
421 double * work1 =
new double[ MAX_NUM_VEC * MAX_NUM_VEC ];
422 double * work3 =
new double[ MAX_NUM_VEC * MAX_NUM_VEC ];
423 double * work2 =
new double[ MAX_NUM_VEC * NUM_SOLUTIONS ];
425 for (
int solution = 0; solution < NUM_SOLUTIONS; solution++ ){
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 ];
445 dgemm_( ¬rans, ¬rans, &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_( ¬rans, ¬rans, &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_( ¬rans, ¬rans, &MAX_NUM_VEC, &inc1, &MAX_NUM_VEC, &one, work1, &MAX_NUM_VEC, mxM_rhs, &MAX_NUM_VEC, &
set, work3, &MAX_NUM_VEC );
455 dsyev_( &jobz, &uplo, &MAX_NUM_VEC, mxM_vecs, &MAX_NUM_VEC, mxM_eigs, mxM_work, &mxM_lwork, &info );
465 dgemm_( &trans, ¬rans, &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;
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 );
474 cout <<
"WARNING AT DAVIDSON : The eigenvalue " << mxM_eigs[ diag ] <<
" to solve Ax = b has been overwritten with " << current_eigenvalue <<
"." << endl;
477 mxM_work[ diag ] = mxM_work[ diag ] / current_eigenvalue;
480 dgemm_( ¬rans, ¬rans, &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 );
486 double * ptr = work2 + MAX_NUM_VEC * solution;
487 const double twonorm = sqrt( ddot_( &MAX_NUM_VEC, ptr, &inc1, ptr, &inc1 ) );
491 double factor = 1.0 / twonorm;
492 dscal_( &MAX_NUM_VEC, &factor, ptr, &inc1 );
500 int size = MAX_NUM_VEC * NUM_SOLUTIONS;
501 dcopy_( &size, work2, &inc1, mxM_vecs, &inc1 );
510 void CheMPS2::Davidson::MxMafterDeflation(){
514 if ( problem_type ==
'E' ){
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 ];
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 ];
531 for (
int ivec = 0; ivec < NUM_VEC_KEEP; ivec++ ){
532 mxM_rhs[ ivec ] = ddot_( &veclength, Hvecs[ ivec ], &inc1, RHS, &inc1 );
virtual ~Davidson()
Destructor.
int GetNumMultiplications() const
Get the number of matrix vector multiplications which have been performed.
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.
char FetchInstruction(double **pointers)
The iterator to converge the ground state vector.