CheMPS2
DIIS.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 <fstream>
24 #include <string>
25 #include <sstream>
26 #include <math.h>
27 
28 #include "MyHDF5.h"
29 #include "Lapack.h"
30 #include "DIIS.h"
31 
32 using std::string;
33 using std::ifstream;
34 using std::cout;
35 using std::endl;
36 
37 CheMPS2::DIIS::DIIS(const int numVarsParamIn, const int numVarsErrorIn, const int numVecsIn){
38 
39  numVarsParam = numVarsParamIn;
40  numVarsError = numVarsErrorIn;
41  numVecs = numVecsIn;
42 
43  errorVectors = new double*[numVecs];
44  paramVectors = new double*[numVecs];
45  currentNumVecs = 0;
46 
47  lastLinco = new double[numVarsParam];
48 
49 }
50 
52 
53  for (int cnt=0; cnt<currentNumVecs; cnt++){
54  delete [] errorVectors[cnt];
55  delete [] paramVectors[cnt];
56  }
57 
58  delete [] errorVectors;
59  delete [] paramVectors;
60  delete [] lastLinco;
61 
62 }
63 
64 int CheMPS2::DIIS::getNumVarsParam() const{ return numVarsParam; }
65 
66 int CheMPS2::DIIS::getNumVarsError() const{ return numVarsError; }
67 
68 int CheMPS2::DIIS::getNumVecs() const{ return numVecs; }
69 
70 int CheMPS2::DIIS::getCurrentNumVecs() const{ return currentNumVecs; }
71 
72 double * CheMPS2::DIIS::getLastLinco(){ return lastLinco; }
73 
74 void CheMPS2::DIIS::appendNew(double * newError, double * newParam){
75 
76  if (currentNumVecs==numVecs){
77 
78  double * ptrE = errorVectors[0];
79  double * ptrP = paramVectors[0];
80 
81  for (int cnt=1; cnt<numVecs; cnt++){
82  errorVectors[cnt-1] = errorVectors[cnt];
83  paramVectors[cnt-1] = paramVectors[cnt];
84  }
85 
86  errorVectors[numVecs-1] = ptrE;
87  paramVectors[numVecs-1] = ptrP;
88 
89  } else {
90 
91  errorVectors[currentNumVecs] = new double[numVarsError];
92  paramVectors[currentNumVecs] = new double[numVarsParam];
93  currentNumVecs += 1;
94 
95  }
96 
97  int inc = 1;
98  dcopy_(&numVarsError, newError, &inc, errorVectors[currentNumVecs-1], &inc);
99  dcopy_(&numVarsParam, newParam, &inc, paramVectors[currentNumVecs-1], &inc);
100 
101 }
102 
103 void CheMPS2::DIIS::calculateParam(double * newParam){
104 
105  int lindim = currentNumVecs + 1;
106  int inc = 1;
107 
108  //The symmetric Pulay matrix
109  double * matrix = new double[ lindim * lindim ];
110  matrix[currentNumVecs*(1+lindim)] = 0.0;
111  for (int cnt=0; cnt<currentNumVecs; cnt++){
112  matrix[currentNumVecs + cnt*lindim] = 1.0 ;
113  matrix[cnt + currentNumVecs*lindim] = 1.0 ;
114  for (int cnt2=cnt; cnt2<currentNumVecs; cnt2++){
115  matrix[cnt + cnt2*lindim] = ddot_(&numVarsError, errorVectors[cnt], &inc, errorVectors[cnt2], &inc);
116  matrix[cnt2 + cnt*lindim] = matrix[cnt + cnt2*lindim];
117  }
118  }
119 
120  //Eigenvalue decomposition of the Pulay matrix
121  char jobz = 'V';
122  char uplo = 'U';
123  double * array = new double[ lindim ]; //eigs
124  int lwork = 3 * lindim;
125  double * work = new double[ lwork ];
126  int info;
127  dsyev_(&jobz, &uplo, &lindim, matrix, &lindim, array, work, &lwork, &info);
128 
129  //The desired coefficients: [vec(c), lambda] = V (1/eigs) V^T [vec(0), 1]
130  //Step 1: [vec(0), 1] --> work
131  for (int cnt=0; cnt<currentNumVecs; cnt++){ work[cnt] = 0.0; }
132  work[currentNumVecs] = 1.0;
133 
134  //Step 2: V^T [vec(0), 1] --> work+lindim
135  char trans = 'T';
136  char notra = 'N';
137  int one = 1;
138  double alpha = 1.0;
139  double beta = 0.0;
140  dgemm_(&trans, &notra, &lindim, &one, &lindim, &alpha, matrix, &lindim, work, &lindim, &beta, work+lindim, &lindim);
141 
142  //Step 3: (1/eigs) V^T [vec(0), 1] --> work+lindim
143  for (int cnt=0; cnt<lindim; cnt++){ work[lindim + cnt] /= array[cnt]; }
144 
145  //Step 4: V (1/eigs) V^T [vec(0), 1] --> work
146  dgemm_(&notra, &notra, &lindim, &one, &lindim, &alpha, matrix, &lindim, work+lindim, &lindim, &beta, work, &lindim);
147 
148  //Fill newParam and make a copy in lastLinco
149  for (int cnt=0; cnt<numVarsParam; cnt++){ newParam[cnt] = 0.0; }
150  for (int cnt=0; cnt<currentNumVecs; cnt++){ daxpy_(&numVarsParam, work+cnt, paramVectors[cnt], &inc, newParam, &inc); }
151  dcopy_(&numVarsParam, newParam, &inc, lastLinco, &inc);
152 
153  //Print out the DIIS coefficients
154  cout << " DIIS::calculateParam : coefficients (newer vectors --> older vectors) : ";
155  for (int cnt=0; cnt<currentNumVecs; cnt++){ cout << work[currentNumVecs-1-cnt] << "\t"; }
156  cout << endl;
157 
158  delete [] matrix;
159  delete [] array;
160  delete [] work;
161 
162 }
163 
164 
165 void CheMPS2::DIIS::saveDIIS(const string filename) const{
166 
167  hid_t file_id = H5Fcreate(filename.c_str(), H5F_ACC_TRUNC, H5P_DEFAULT, H5P_DEFAULT);
168  hid_t group_id = H5Gcreate(file_id, "/Data", H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
169 
170  //The number of vectors saved
171  hsize_t dimarray1 = 1;
172  hid_t dataspace_id1 = H5Screate_simple(1, &dimarray1, NULL);
173  hid_t dataset_id1 = H5Dcreate(group_id, "currentNumVecs", H5T_STD_I32LE, dataspace_id1, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
174  H5Dwrite(dataset_id1, H5T_NATIVE_INT, H5S_ALL, H5S_ALL, H5P_DEFAULT, &currentNumVecs);
175  H5Dclose(dataset_id1);
176  H5Sclose(dataspace_id1);
177 
178  //The parameter vector sizes
179  hsize_t dimarray2 = 1;
180  hid_t dataspace_id2 = H5Screate_simple(1, &dimarray2, NULL);
181  hid_t dataset_id2 = H5Dcreate(group_id, "numVarsParam", H5T_STD_I32LE, dataspace_id2, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
182  H5Dwrite(dataset_id2, H5T_NATIVE_INT, H5S_ALL, H5S_ALL, H5P_DEFAULT, &numVarsParam);
183  H5Dclose(dataset_id2);
184  H5Sclose(dataspace_id2);
185 
186  //The error vector sizes
187  hsize_t dimarray5 = 1;
188  hid_t dataspace_id5 = H5Screate_simple(1, &dimarray5, NULL);
189  hid_t dataset_id5 = H5Dcreate(group_id, "numVarsError", H5T_STD_I32LE, dataspace_id5, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
190  H5Dwrite(dataset_id5, H5T_NATIVE_INT, H5S_ALL, H5S_ALL, H5P_DEFAULT, &numVarsError);
191  H5Dclose(dataset_id5);
192  H5Sclose(dataspace_id5);
193 
194  for (int cnt=0; cnt<currentNumVecs; cnt++){
195 
196  //The error vectors
197  std::stringstream nameE;
198  nameE << "error_" << cnt;
199  hsize_t dimarray3 = numVarsError;
200  hid_t dataspace_id3 = H5Screate_simple(1, &dimarray3, NULL);
201  hid_t dataset_id3 = H5Dcreate(group_id, nameE.str().c_str(), H5T_IEEE_F64LE, dataspace_id3, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
202  H5Dwrite(dataset_id3, H5T_NATIVE_DOUBLE, H5S_ALL, H5S_ALL, H5P_DEFAULT, errorVectors[cnt]);
203  H5Dclose(dataset_id3);
204  H5Sclose(dataspace_id3);
205 
206  //The parameter vectors
207  std::stringstream nameP;
208  nameP << "param_" << cnt;
209  hsize_t dimarray4 = numVarsParam;
210  hid_t dataspace_id4 = H5Screate_simple(1, &dimarray4, NULL);
211  hid_t dataset_id4 = H5Dcreate(group_id, nameP.str().c_str(), H5T_IEEE_F64LE, dataspace_id4, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
212  H5Dwrite(dataset_id4, H5T_NATIVE_DOUBLE, H5S_ALL, H5S_ALL, H5P_DEFAULT, paramVectors[cnt]);
213  H5Dclose(dataset_id4);
214  H5Sclose(dataspace_id4);
215 
216  }
217 
218  H5Gclose(group_id);
219  H5Fclose(file_id);
220 
221 }
222 
223 void CheMPS2::DIIS::loadDIIS(const string filename){
224 
225  hid_t file_id = H5Fopen(filename.c_str(), H5F_ACC_RDONLY, H5P_DEFAULT);
226  hid_t group_id = H5Gopen(file_id, "/Data",H5P_DEFAULT);
227 
228  //The parameter vector sizes
229  int numVarsBIS;
230  hid_t dataset_id2 = H5Dopen(group_id, "numVarsParam", H5P_DEFAULT);
231  H5Dread(dataset_id2, H5T_NATIVE_INT, H5S_ALL, H5S_ALL, H5P_DEFAULT, &numVarsBIS);
232  H5Dclose(dataset_id2);
233  assert( numVarsParam==numVarsBIS );
234 
235  //The error vector sizes
236  hid_t dataset_id5 = H5Dopen(group_id, "numVarsError", H5P_DEFAULT);
237  H5Dread(dataset_id5, H5T_NATIVE_INT, H5S_ALL, H5S_ALL, H5P_DEFAULT, &numVarsBIS);
238  H5Dclose(dataset_id5);
239  assert( numVarsError==numVarsBIS );
240 
241  //The number of vectors saved
242  int currentNumVecsBIS;
243  hid_t dataset_id1 = H5Dopen(group_id, "currentNumVecs", H5P_DEFAULT);
244  H5Dread(dataset_id1, H5T_NATIVE_INT, H5S_ALL, H5S_ALL, H5P_DEFAULT, &currentNumVecsBIS);
245  H5Dclose(dataset_id1);
246  assert( currentNumVecsBIS<=numVecs );
247  assert( currentNumVecsBIS>=0 );
248 
249  //Create just enough storage for the vectors to be loaded
250  if (currentNumVecs < currentNumVecsBIS){
251  for (int cnt=currentNumVecs; cnt<currentNumVecsBIS; cnt++){
252  errorVectors[cnt] = new double[numVarsError];
253  paramVectors[cnt] = new double[numVarsParam];
254  }
255  currentNumVecs = currentNumVecsBIS;
256  }
257 
258  if (currentNumVecs > currentNumVecsBIS){
259  for (int cnt=currentNumVecs; cnt>currentNumVecsBIS; cnt--){
260  delete [] errorVectors[cnt-1];
261  delete [] paramVectors[cnt-1];
262  }
263  currentNumVecs = currentNumVecsBIS;
264  }
265 
266  for (int cnt=0; cnt<currentNumVecs; cnt++){
267 
268  //The error vectors
269  std::stringstream nameE;
270  nameE << "error_" << cnt;
271  hid_t dataset_id3 = H5Dopen(group_id, nameE.str().c_str(), H5P_DEFAULT);
272  H5Dread(dataset_id3, H5T_NATIVE_DOUBLE, H5S_ALL, H5S_ALL, H5P_DEFAULT, errorVectors[cnt]);
273  H5Dclose(dataset_id3);
274 
275  //The parameter vectors
276  std::stringstream nameP;
277  nameP << "param_" << cnt;
278  hid_t dataset_id4 = H5Dopen(group_id, nameP.str().c_str(), H5P_DEFAULT);
279  H5Dread(dataset_id4, H5T_NATIVE_DOUBLE, H5S_ALL, H5S_ALL, H5P_DEFAULT, paramVectors[cnt]);
280  H5Dclose(dataset_id4);
281 
282  }
283 
284  H5Gclose(group_id);
285  H5Fclose(file_id);
286 
287 }
288 
int getNumVarsError() const
Get the number of variables in the error vectors.
Definition: DIIS.cpp:66
DIIS(const int numVarsParamIn, const int numVarsErrorIn, const int numVecsIn)
Constructor.
Definition: DIIS.cpp:37
void appendNew(double *newError, double *newParam)
Append a new error and parameter vector.
Definition: DIIS.cpp:74
void loadDIIS(const string filename=DMRGSCF_diis_storage_name)
Load the DIIS object from disk.
Definition: DIIS.cpp:223
int getNumVarsParam() const
Get the number of variables in the parameter vectors.
Definition: DIIS.cpp:64
void saveDIIS(const string filename=DMRGSCF_diis_storage_name) const
Save the DIIS object to disk.
Definition: DIIS.cpp:165
virtual ~DIIS()
Destructor.
Definition: DIIS.cpp:51
double * getLastLinco()
Get pointer to the last linco of parameter vectors.
Definition: DIIS.cpp:72
void calculateParam(double *newParam)
Calculate the new parameter vector, based on the just appended error and parameter vectors...
Definition: DIIS.cpp:103
int getCurrentNumVecs() const
Get the current number of vectors which are used for DIIS.
Definition: DIIS.cpp:70
int getNumVecs() const
Get the max. number of parameter and error vectors which are kept.
Definition: DIIS.cpp:68