39 numVarsParam = numVarsParamIn;
40 numVarsError = numVarsErrorIn;
43 errorVectors =
new double*[numVecs];
44 paramVectors =
new double*[numVecs];
47 lastLinco =
new double[numVarsParam];
53 for (
int cnt=0; cnt<currentNumVecs; cnt++){
54 delete [] errorVectors[cnt];
55 delete [] paramVectors[cnt];
58 delete [] errorVectors;
59 delete [] paramVectors;
76 if (currentNumVecs==numVecs){
78 double * ptrE = errorVectors[0];
79 double * ptrP = paramVectors[0];
81 for (
int cnt=1; cnt<numVecs; cnt++){
82 errorVectors[cnt-1] = errorVectors[cnt];
83 paramVectors[cnt-1] = paramVectors[cnt];
86 errorVectors[numVecs-1] = ptrE;
87 paramVectors[numVecs-1] = ptrP;
91 errorVectors[currentNumVecs] =
new double[numVarsError];
92 paramVectors[currentNumVecs] =
new double[numVarsParam];
98 dcopy_(&numVarsError, newError, &inc, errorVectors[currentNumVecs-1], &inc);
99 dcopy_(&numVarsParam, newParam, &inc, paramVectors[currentNumVecs-1], &inc);
105 int lindim = currentNumVecs + 1;
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];
123 double * array =
new double[ lindim ];
124 int lwork = 3 * lindim;
125 double * work =
new double[ lwork ];
127 dsyev_(&jobz, &uplo, &lindim, matrix, &lindim, array, work, &lwork, &info);
131 for (
int cnt=0; cnt<currentNumVecs; cnt++){ work[cnt] = 0.0; }
132 work[currentNumVecs] = 1.0;
140 dgemm_(&trans, ¬ra, &lindim, &one, &lindim, &alpha, matrix, &lindim, work, &lindim, &beta, work+lindim, &lindim);
143 for (
int cnt=0; cnt<lindim; cnt++){ work[lindim + cnt] /= array[cnt]; }
146 dgemm_(¬ra, ¬ra, &lindim, &one, &lindim, &alpha, matrix, &lindim, work+lindim, &lindim, &beta, work, &lindim);
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);
154 cout <<
" DIIS::calculateParam : coefficients (newer vectors --> older vectors) : ";
155 for (
int cnt=0; cnt<currentNumVecs; cnt++){ cout << work[currentNumVecs-1-cnt] <<
"\t"; }
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);
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, ¤tNumVecs);
175 H5Dclose(dataset_id1);
176 H5Sclose(dataspace_id1);
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);
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);
194 for (
int cnt=0; cnt<currentNumVecs; cnt++){
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);
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);
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);
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 );
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 );
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, ¤tNumVecsBIS);
245 H5Dclose(dataset_id1);
246 assert( currentNumVecsBIS<=numVecs );
247 assert( currentNumVecsBIS>=0 );
250 if (currentNumVecs < currentNumVecsBIS){
251 for (
int cnt=currentNumVecs; cnt<currentNumVecsBIS; cnt++){
252 errorVectors[cnt] =
new double[numVarsError];
253 paramVectors[cnt] =
new double[numVarsParam];
255 currentNumVecs = currentNumVecsBIS;
258 if (currentNumVecs > currentNumVecsBIS){
259 for (
int cnt=currentNumVecs; cnt>currentNumVecsBIS; cnt--){
260 delete [] errorVectors[cnt-1];
261 delete [] paramVectors[cnt-1];
263 currentNumVecs = currentNumVecsBIS;
266 for (
int cnt=0; cnt<currentNumVecs; cnt++){
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);
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);
int getNumVarsError() const
Get the number of variables in the error vectors.
DIIS(const int numVarsParamIn, const int numVarsErrorIn, const int numVecsIn)
Constructor.
void appendNew(double *newError, double *newParam)
Append a new error and parameter vector.
void loadDIIS(const string filename=DMRGSCF_diis_storage_name)
Load the DIIS object from disk.
int getNumVarsParam() const
Get the number of variables in the parameter vectors.
void saveDIIS(const string filename=DMRGSCF_diis_storage_name) const
Save the DIIS object to disk.
virtual ~DIIS()
Destructor.
double * getLastLinco()
Get pointer to the last linco of parameter vectors.
void calculateParam(double *newParam)
Calculate the new parameter vector, based on the just appended error and parameter vectors...
int getCurrentNumVecs() const
Get the current number of vectors which are used for DIIS.
int getNumVecs() const
Get the max. number of parameter and error vectors which are kept.