iterative-solver 0.0
helper-implementation.h
1#ifndef LINEARALGEBRA_SRC_MOLPRO_LINALG_ITERATIVESOLVER_HELPER_IMPLEMENTATION_H_
2#define LINEARALGEBRA_SRC_MOLPRO_LINALG_ITERATIVESOLVER_HELPER_IMPLEMENTATION_H_
3#include <Eigen/Dense>
4#include <cmath>
5#include <cstddef>
6#include <molpro/Profiler.h>
7#include <molpro/lapacke.h>
8#include <molpro/linalg/itsolv/helper.h>
9
10namespace molpro::linalg::itsolv {
11
12template <typename value_type>
13std::list<SVD<value_type>> svd_eigen_jacobi(size_t nrows, size_t ncols, const array::Span<value_type>& m,
14 double threshold) {
15 auto mat = Eigen::Map<const Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic>>(m.data(), nrows, ncols);
16 auto svd = Eigen::JacobiSVD<Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic>, Eigen::NoQRPreconditioner>(
17 mat, Eigen::ComputeThinV);
18 auto svd_system = std::list<SVD<value_type>>{};
19 auto sv = svd.singularValues();
20 for (int i = int(ncols) - 1; i >= 0; --i) {
21 if (std::abs(sv(i)) < threshold) { // TODO: This seems to discard values ABOVE the threshold, not below it. it's
22 auto t = SVD<value_type>{}; // also not scaling this threshold relative to the max singular value - find out why
23 t.value = sv(i);
24 t.v.reserve(ncols);
25 for (size_t j = 0; j < ncols; ++j) {
26 t.v.emplace_back(svd.matrixV()(j, i));
27 }
28 svd_system.emplace_back(std::move(t));
29 }
30 }
31 return svd_system;
32}
33
34template <typename value_type>
35std::list<SVD<value_type>> svd_eigen_bdcsvd(size_t nrows, size_t ncols, const array::Span<value_type>& m,
36 double threshold) {
37 auto mat = Eigen::Map<const Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic>>(m.data(), nrows, ncols);
38 auto svd = Eigen::BDCSVD<Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic>>(mat, Eigen::ComputeThinV);
39 auto svd_system = std::list<SVD<value_type>>{};
40 auto sv = svd.singularValues();
41 for (int i = int(ncols) - 1; i >= 0; --i) {
42 if (std::abs(sv(i)) < threshold) {
43 auto t = SVD<value_type>{};
44 t.value = sv(i);
45 t.v.reserve(ncols);
46 for (size_t j = 0; j < ncols; ++j) {
47 t.v.emplace_back(svd.matrixV()(j, i));
48 }
49 svd_system.emplace_back(std::move(t));
50 }
51 }
52 return svd_system;
53}
54
55#if defined HAVE_CBLAS
56template <typename value_type>
57std::list<SVD<value_type>> svd_lapacke_dgesdd(size_t nrows, size_t ncols, const array::Span<value_type>& mat,
58 double threshold) {
59 int info;
60 int m = nrows;
61 int n = ncols;
62 int sdim = std::min(m, n);
63 std::vector<double> sv(sdim), u(nrows * nrows), v(ncols * ncols);
64 info = LAPACKE_dgesdd(LAPACK_ROW_MAJOR, 'A', int(nrows), int(ncols), const_cast<double*>(mat.data()), int(ncols),
65 sv.data(), u.data(), int(nrows), v.data(), int(ncols));
66 auto svd_system = std::list<SVD<value_type>>{};
67 for (int i = int(ncols) - 1; i >= 0; --i) {
68 if (std::abs(sv[i]) < threshold) {
69 auto t = SVD<value_type>{};
70 t.value = sv[i];
71 t.v.reserve(ncols);
72 for (size_t j = 0; j < ncols; ++j) {
73 t.v.emplace_back(v[i * ncols + j]);
74 }
75 svd_system.emplace_back(std::move(t));
76 }
77 }
78 return svd_system;
79}
80
81template <typename value_type>
82std::list<SVD<value_type>> svd_lapacke_dgesvd(size_t nrows, size_t ncols, const array::Span<value_type>& mat,
83 double threshold) {
84 int info;
85 int m = nrows;
86 int n = ncols;
87 int sdim = std::min(m, n);
88 std::vector<double> sv(sdim), u(nrows * nrows), v(ncols * ncols);
89 double superb[sdim - 1];
90 info = LAPACKE_dgesvd(LAPACK_ROW_MAJOR, 'N', 'A', int(nrows), int(ncols), const_cast<double*>(mat.data()), int(ncols),
91 sv.data(), u.data(), int(nrows), v.data(), int(ncols), superb);
92 auto svd_system = std::list<SVD<value_type>>{};
93 for (int i = int(ncols) - 1; i >= 0; --i) {
94 if (std::abs(sv[i]) < threshold) {
95 auto t = SVD<value_type>{};
96 t.value = sv[i];
97 t.v.reserve(ncols);
98 for (size_t j = 0; j < ncols; ++j) {
99 t.v.emplace_back(v[i * ncols + j]);
100 }
101 svd_system.emplace_back(std::move(t));
102 }
103 }
104 return svd_system;
105}
106
107#endif
108
109#ifdef MOLPRO
110extern "C" int dsyev_c(char, char, int, double*, int, double*);
111#endif
112
122inline int eigensolver_lapacke_dsyev(const std::vector<double>& matrix, std::vector<double>& eigenvectors,
123 std::vector<double>& eigenvalues, const size_t dimension) {
124
125 // validate input
126 if (eigenvectors.size() != matrix.size()) {
127 throw std::runtime_error("Matrix of eigenvectors and input matrix are not the same size!");
128 }
129
130 if (eigenvectors.size() != dimension * dimension || eigenvalues.size() != dimension) {
131 throw std::runtime_error("Size of eigenvectors/eigenvlaues do not match dimension!");
132 }
133
134 // magic letters
135 // static const char compute_eigenvalues = 'N';
136 static const char compute_eigenvalues_eigenvectors = 'V';
137 // static const char store_upper_triangle = 'U';
138 static const char store_lower_triangle = 'L';
139
140 // copy input matrix (lapack overwrites)
141 std::copy(matrix.begin(), matrix.end(), eigenvectors.begin());
142
143 // set lapack vars
144 lapack_int status;
145 lapack_int leading_dimension = dimension;
146 lapack_int order = dimension;
147
148 // call to lapack
149#ifdef MOLPRO
150 status = dsyev_c(compute_eigenvalues_eigenvectors, store_lower_triangle, order, eigenvectors.data(),
151 leading_dimension, eigenvalues.data());
152#else
153 status = LAPACKE_dsyev(LAPACK_COL_MAJOR, compute_eigenvalues_eigenvectors, store_lower_triangle, order,
154 eigenvectors.data(), leading_dimension, eigenvalues.data());
155#endif
156
157 return status;
158}
159
167inline std::list<SVD<double>> eigensolver_lapacke_dsyev(size_t dimension, std::vector<double>& matrix) {
168 std::vector<double> eigvecs(dimension * dimension);
169 std::vector<double> eigvals(dimension);
170
171 // call to lapack
172 int success = eigensolver_lapacke_dsyev(matrix, eigvecs, eigvals, dimension);
173 if (success < 0) {
174 throw std::invalid_argument("Invalid argument of eigensolver_lapacke_dsyev: ");
175 }
176 if (success > 0) {
177 throw std::runtime_error("Lapacke_dsyev (eigensolver) failed to converge. "
178 " elements of an intermediate tridiagonal form did not converge to zero.");
179 }
180
181 auto eigensystem = std::list<SVD<double>>{};
182
183 // populate eigensystem
184 for (int i = dimension - 1; i >= 0;
185 i--) { // note: flipping this axis gives parity with results of eigen::jacobiSVD
186 auto temp_eigenproblem = SVD<double>{};
187 temp_eigenproblem.value = eigvals[i];
188 for (size_t j = 0; j < dimension; j++) {
189 temp_eigenproblem.v.emplace_back(eigvecs[j + (dimension * i)]);
190 }
191 eigensystem.emplace_back(temp_eigenproblem);
192 }
193
194 return eigensystem;
195}
196
205inline std::list<SVD<double>> eigensolver_lapacke_dsyev(size_t dimension,
207 // TODO: this should be the other way around, eigensolver_lapacke_dsyev should take a span by default and this should
208 // wrap it with a vector
209 std::vector<double> v;
210 v.insert(v.begin(), matrix.begin(), matrix.end());
211 return eigensolver_lapacke_dsyev(dimension, v);
212}
213
221template <typename value_type>
222size_t get_rank(std::vector<value_type> eigenvalues, value_type threshold) {
223 if (eigenvalues.size() == 0) {
224 return 0;
225 }
226 value_type max = *max_element(eigenvalues.begin(), eigenvalues.end());
227 value_type threshold_scaled = threshold * max;
228 size_t count =
229 std::count_if(eigenvalues.begin(), eigenvalues.end(), [&](auto const& val) { return val >= threshold_scaled; });
230 return count;
231}
232
240template <typename value_type>
241size_t get_rank(std::list<SVD<value_type>> svd_system, value_type threshold) {
242 // compute max
243 value_type max_value = 0;
244 std::list<SVD<double>>::iterator it;
245 for (it = svd_system.begin(); it != svd_system.end(); it++) {
246 if (it->value > max_value) {
247 max_value = it->value;
248 }
249 }
250 // scale threshold
251 value_type threshold_scaled = threshold * max_value;
252
253 size_t rank = 0;
254 // get rank
255 for (it = svd_system.begin(); it != svd_system.end(); it++) {
256 if (it->value > threshold_scaled) {
257 rank += 1;
258 }
259 }
260 return rank;
261}
262
263template <typename value_type, typename std::enable_if_t<!is_complex<value_type>{}, std::nullptr_t>>
264std::list<SVD<value_type>> svd_system(size_t nrows, size_t ncols, const array::Span<value_type>& m, double threshold,
265 bool hermitian, bool reduce_to_rank) {
266 std::list<SVD<value_type>> svds;
267 assert(m.size() == nrows * ncols);
268 if (m.empty())
269 return {};
270 if (hermitian) {
271 assert(nrows == ncols);
272 svds = eigensolver_lapacke_dsyev(nrows, m);
273 for (auto s = svds.begin(); s != svds.end();)
274 if (s->value > threshold)
275 s = svds.erase(s);
276 else
277 ++s;
278 } else {
279 //#if defined HAVE_LAPACKE
280 // if (nrows > 16)
281 // return svd_lapacke_dgesdd<value_type>(nrows, ncols, m, threshold);
282 // return svd_lapacke_dgesvd<value_type>(nrows, ncols, m, threshold);
283 //#endif
284 svds = svd_eigen_jacobi<value_type>(nrows, ncols, m, threshold);
285 // return svd_eigen_bdcsvd<value_type>(nrows, ncols, m, threshold);
286 }
287
288 // reduce to rank
289 if (reduce_to_rank) {
290 int rank = get_rank(svds, threshold);
291 for (int i = ncols; i > rank; i--) {
292 svds.pop_back();
293 }
294 }
295 return svds;
296}
297
298template <typename value_type, typename std::enable_if_t<is_complex<value_type>{}, int>>
299std::list<SVD<value_type>> svd_system(size_t nrows, size_t ncols, const array::Span<value_type>& m, double threshold,
300 bool hermitian, bool reduce_to_rank) {
301 assert(false); // Complex not implemented here
302 return {};
303}
304
305template <typename value_type>
306void printMatrix(const std::vector<value_type>& m, size_t rows, size_t cols, std::string title, std::ostream& s) {
307 s << title << "\n"
308 << Eigen::Map<const Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic>>(m.data(), rows, cols) << std::endl;
309}
310
311template <typename value_type, typename std::enable_if_t<is_complex<value_type>{}, int>>
312void eigenproblem(std::vector<value_type>& eigenvectors, std::vector<value_type>& eigenvalues,
313 const std::vector<value_type>& matrix, const std::vector<value_type>& metric, size_t dimension,
314 bool hermitian, double svdThreshold, int verbosity, bool condone_complex) {
315 assert(false); // Complex not implemented here
316}
317
318template <typename value_type, typename std::enable_if_t<!is_complex<value_type>{}, std::nullptr_t>>
319void eigenproblem(std::vector<value_type>& eigenvectors, std::vector<value_type>& eigenvalues,
320 const std::vector<value_type>& matrix, const std::vector<value_type>& metric, size_t dimension,
321 bool hermitian, double svdThreshold, int verbosity, bool condone_complex) {
322 auto prof = molpro::Profiler::single();
323 prof->start("itsolv::eigenproblem");
324 Eigen::Map<const Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>> HrowMajor(
325 matrix.data(), dimension, dimension);
326 Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic> H(dimension, dimension);
327 H = HrowMajor;
328 Eigen::Map<const Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic>> S(metric.data(), dimension, dimension);
329 Eigen::MatrixXcd subspaceEigenvectors; // FIXME templating
330 Eigen::VectorXcd subspaceEigenvalues; // FIXME templating
331 // Eigen::GeneralizedEigenSolver<Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic>> s(H, S);
332
333 // initialisation of variables
334 Eigen::VectorXd singularValues;
335 Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor> matrixV;
336 Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor> matrixU;
337 std::vector<double> eigvecs;
338 std::vector<double> eigvals;
339 int rank;
340
341 // if the matrix is hermitian, we can use lapacke_dsyev
342 if (hermitian) {
343 eigvecs.resize(dimension * dimension);
344 eigvals.resize(dimension);
345 int success = eigensolver_lapacke_dsyev(metric, eigvecs, eigvals, dimension);
346 if (success != 0) {
347 throw std::runtime_error("Eigensolver did not converge");
348 }
349 singularValues = Eigen::Map<Eigen::VectorXd>(eigvals.data(), dimension);
350 matrixV = Eigen::Map<Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic, Eigen::ColMajor>>(
351 eigvecs.data(), dimension, dimension);
352 matrixU = Eigen::Map<Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic, Eigen::ColMajor>>(
353 eigvecs.data(), dimension, dimension);
354 rank = get_rank(eigvals, svdThreshold);
355 } else {
356 Eigen::JacobiSVD<Eigen::MatrixXd> svd(S, Eigen::ComputeThinU | Eigen::ComputeThinV);
357 singularValues = svd.singularValues();
358 matrixV = svd.matrixV();
359 matrixU = svd.matrixU();
360 rank = svd.rank();
361 }
362
363 // svd.setThreshold(svdThreshold);
364 // molpro::cout << "singular values of overlap " << svd.singularValues().transpose() << std::endl;
365 // auto Hbar = svd.solve(H);
366 if (verbosity > 1 && rank < S.cols())
367 molpro::cout << "SVD rank " << rank << " in subspace of dimension " << S.cols() << std::endl;
368 if (verbosity > 2 && rank < S.cols())
369 molpro::cout << "singular values " << singularValues.transpose() << std::endl;
370 auto svmh = singularValues.head(rank);
371 for (auto k = 0; k < rank; k++)
372 svmh(k) = svmh(k) > 1e-14 ? 1 / std::sqrt(svmh(k)) : 0;
373 auto Hbar =
374 (svmh.asDiagonal()) * (matrixU.leftCols(rank).adjoint()) * H * matrixV.leftCols(rank) * (svmh.asDiagonal());
375 // std::cout << "\n\nHbar: \n" << Hbar << "\n\n";
376 // molpro::cout << "S\n"<<S<<std::endl;
377 // molpro::cout << "S singular values"<<(Eigen::DiagonalMatrix<value_type, Eigen::Dynamic,
378 // Eigen::Dynamic>(svd.singularValues().head(svd.rank())))<<std::endl; molpro::cout << "S inverse singular
379 // values"<<Eigen::DiagonalMatrix<value_type,
380 // Eigen::Dynamic>(svd.singularValues().head(svd.rank())).inverse()<<std::endl; molpro::cout << "S singular
381 // values"<<sv<<std::endl; molpro::cout << "H\n"<<H<<std::endl; molpro::cout << "Hbar\n"<<Hbar<<std::endl;
382 Eigen::EigenSolver<Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic>> s(Hbar);
383 // molpro::cout << "s.eigenvectors()\n"<<s.eigenvectors()<<std::endl;
384 subspaceEigenvalues = s.eigenvalues();
385 if (s.eigenvalues().imag().norm() < 1e-10) { // real eigenvalues
386 // molpro::cout << "eigenvalues near-enough real" << std::endl;
387 subspaceEigenvalues = subspaceEigenvalues.real();
388 subspaceEigenvectors = s.eigenvectors();
389 // complex eigenvectors need to be rotated
390 // assume that they come in consecutive pairs
391 for (int i = 0; i < subspaceEigenvectors.cols(); i++) {
392 if (subspaceEigenvectors.col(i).imag().norm() > 1e-10) {
393 int j = i + 1;
394 if (std::abs(subspaceEigenvalues(i) - subspaceEigenvalues(j)) < 1e-10 and
395 subspaceEigenvectors.col(j).imag().norm() > 1e-10) {
396 subspaceEigenvectors.col(j) = subspaceEigenvectors.col(i).imag() / subspaceEigenvectors.col(i).imag().norm();
397 subspaceEigenvectors.col(i) = subspaceEigenvectors.col(i).real() / subspaceEigenvectors.col(i).real().norm();
398 }
399 }
400 }
401 subspaceEigenvectors = matrixV.leftCols(rank) * svmh.asDiagonal() * subspaceEigenvectors;
402 } else { // complex eigenvectors
403// molpro::cout << "eigenvalues not near-enough real"<<std::endl;
404// molpro::cout << "s.eigenvalues() "<< s.eigenvalues().transpose()<<std::endl;
405#ifdef __INTEL_COMPILER
406 molpro::cout << "Hbar\n" << Hbar << std::endl;
407 molpro::cout << "Eigenvalues\n" << s.eigenvalues() << std::endl;
408 molpro::cout << "Eigenvectors\n" << s.eigenvectors() << std::endl;
409 throw std::runtime_error("Intel compiler does not support working with complex eigen3 entities properly");
410#endif
411 subspaceEigenvectors = matrixV.leftCols(rank) * svmh.asDiagonal() * s.eigenvectors();
412 // std::cout << "subspaceEigenvectors\n" << subspaceEigenvectors << std::endl;
413 }
414
415 {
416 // sort
417 auto eigval = subspaceEigenvalues;
418 auto eigvec = subspaceEigenvectors;
419 std::vector<Eigen::Index> map;
420 for (Eigen::Index k = 0; k < Hbar.cols(); k++) {
421 Eigen::Index ll;
422 for (ll = 0; std::count(map.begin(), map.end(), ll) != 0; ll++)
423 ;
424 for (Eigen::Index l = 0; l < Hbar.cols(); l++) {
425 if (std::count(map.begin(), map.end(), l) == 0) {
426 if (eigval(l).real() < eigval(ll).real())
427 ll = l;
428 }
429 }
430 map.push_back(ll);
431 subspaceEigenvalues(k) = eigval(ll);
432 // molpro::cout << "new sorted eigenvalue "<<k<<", "<<ll<<", "<<eigval(ll)<<std::endl;
433 // molpro::cout << eigvec.col(ll)<<std::endl;
434 subspaceEigenvectors.col(k) = eigvec.col(ll);
435 double maxcomp =0;
436 for (Eigen::Index l = 0; l < Hbar.cols(); l++) {
437 if (std::abs(subspaceEigenvectors.col(k)[l].real()) > std::abs(subspaceEigenvectors.col(k)[maxcomp].real()))
438 maxcomp = l;
439 }
440 if (subspaceEigenvectors.col(k)[maxcomp].real() < 0)
441 subspaceEigenvectors.col(k) = - subspaceEigenvectors.col(k);
442 }
443 }
444
445 // TODO: Need to address the case of near-zero eigenvalues (as below for non-hermitian case) and clean-up
446 // non-hermitian case
447
448 // molpro::cout << "sorted eigenvalues\n"<<subspaceEigenvalues<<std::endl;
449 // molpro::cout << "sorted eigenvectors\n"<<subspaceEigenvectors<<std::endl;
450 // molpro::cout << "hermitian="<<hermitian<<std::endl;
451 if (!hermitian) {
452 Eigen::MatrixXcd ovlTimesVec(subspaceEigenvectors.cols(), subspaceEigenvectors.rows()); // FIXME templating
453 for (auto repeat = 0; repeat < 3; ++repeat)
454 for (Eigen::Index k = 0; k < subspaceEigenvectors.cols(); k++) {
455 if (std::abs(subspaceEigenvalues(k)) <
456 1e-12) { // special case of zero eigenvalue -- make some real non-zero vector definitely in the null space
457 subspaceEigenvectors.col(k).real() += double(0.3256897) * subspaceEigenvectors.col(k).imag();
458 subspaceEigenvectors.col(k).imag().setZero();
459 }
460 if (hermitian)
461 for (Eigen::Index l = 0; l < k; l++) {
462 // auto ovl =
463 // (subspaceEigenvectors.col(l).adjoint() * m_subspaceOverlap * subspaceEigenvectors.col(k))(
464 // 0, 0); (ovlTimesVec.row(l) * subspaceEigenvectors.col(k))(0,0);
465 // ovlTimesVec.row(l).dot(subspaceEigenvectors.col(k));
466 // auto norm =
467 // (subspaceEigenvectors.col(l).adjoint() * subspaceOverlap * subspaceEigenvectors.col(l))(
468 // 0,
469 // 0);
470 // molpro::cout << "k="<<k<<", l="<<l<<", ovl="<<ovl<<" norm="<<norm<<std::endl;
471 // molpro::cout << subspaceEigenvectors.col(k).transpose()<<std::endl;
472 // molpro::cout << subspaceEigenvectors.col(l).transpose()<<std::endl;
473 subspaceEigenvectors.col(k) -= subspaceEigenvectors.col(l) * // ovl;// / norm;
474 ovlTimesVec.row(l).dot(subspaceEigenvectors.col(k));
475 // molpro::cout<<"immediately after projection " << k<<l<<" "<<
476 // (subspaceEigenvectors.col(l).adjoint() * subspaceOverlap * subspaceEigenvectors.col(k))( 0,
477 // 0)<<std::endl;
478 }
479 // for (Eigen::Index l = 0; l < k; l++) molpro::cout<<"after projection loop " << k<<l<<" "<<
480 // (subspaceEigenvectors.col(l).adjoint() * subspaceOverlap * subspaceEigenvectors.col(k))( 0,
481 // 0)<<std::endl; molpro::cout <<
482 // "eigenvector"<<std::endl<<subspaceEigenvectors.col(k).adjoint()<<std::endl;
483 auto ovl =
484 // (subspaceEigenvectors.col(k).adjoint() * subspaceOverlap *
485 // subspaceEigenvectors.col(k))(0,0);
486 subspaceEigenvectors.col(k).adjoint().dot(S * subspaceEigenvectors.col(k));
487 subspaceEigenvectors.col(k) /= std::sqrt(ovl.real());
488 ovlTimesVec.row(k) = subspaceEigenvectors.col(k).adjoint() * S;
489 // for (Eigen::Index l = 0; l < k; l++)
490 // molpro::cout<<"after normalisation " << k<<l<<" "<< (subspaceEigenvectors.col(l).adjoint() *
491 // subspaceOverlap * subspaceEigenvectors.col(k))( 0, 0)<<std::endl; molpro::cout <<
492 // "eigenvector"<<std::endl<<subspaceEigenvectors.col(k).adjoint()<<std::endl;
493 // phase
494 Eigen::Index lmax = 0;
495 for (Eigen::Index l = 0; l < subspaceEigenvectors.rows(); l++) {
496 if (std::abs(subspaceEigenvectors(l, k)) > std::abs(subspaceEigenvectors(lmax, k)))
497 lmax = l;
498 }
499 if (subspaceEigenvectors(lmax, k).real() < 0)
500 subspaceEigenvectors.col(k) = -subspaceEigenvectors.col(k);
501 // for (Eigen::Index l = 0; l < k; l++)
502 // molpro::cout << k<<l<<" "<<
503 // (subspaceEigenvectors.col(l).adjoint() * subspaceOverlap *
504 // subspaceEigenvectors.col(k))( 0, 0)<<std::endl;
505 }
506 }
507 // if (!hermitian) {
508 // molpro::cout << "eigenvalues"<<std::endl<<subspaceEigenvalues<<std::endl;
509 // molpro::cout << "eigenvectors" << std::endl << subspaceEigenvectors << std::endl;
510 // }
511 if (condone_complex) {
512 for (Eigen::Index root = 0; root < Hbar.cols(); ++root) {
513 if (subspaceEigenvalues(root).imag() != 0) {
514 // molpro::cout << "complex eigenvalues: " << subspaceEigenvalues(root) << ", " <<
515 // subspaceEigenvalues(root + 1)
516 // << std::endl;
517 subspaceEigenvalues(root) = subspaceEigenvalues(root + 1) = subspaceEigenvalues(root).real();
518 subspaceEigenvectors.col(root) = subspaceEigenvectors.col(root).real();
519 subspaceEigenvectors.col(root + 1) = subspaceEigenvectors.col(root + 1).imag();
520 ++root;
521 }
522 }
523 }
524 if ((subspaceEigenvectors - subspaceEigenvectors.real()).norm() > 1e-10 or
525 (subspaceEigenvalues - subspaceEigenvalues.real()).norm() > 1e-10) {
526 throw std::runtime_error("unexpected complex solution found");
527 }
528 eigenvectors.resize(dimension * Hbar.cols());
529 eigenvalues.resize(Hbar.cols());
530 // if constexpr (std::is_class<value_type>::value) {
531 Eigen::Map<Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic>>(eigenvectors.data(), dimension, Hbar.cols()) =
532 subspaceEigenvectors.real();
533 Eigen::Map<Eigen::Matrix<value_type, Eigen::Dynamic, 1>> ev(eigenvalues.data(), Hbar.cols());
534 ev = subspaceEigenvalues.real();
535
536 // } else {
537 // Eigen::Map<Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic>>(m_evec_xx.data(), dimension, dimension)
538 // =
539 // subspaceEigenvectors;
540 // Eigen::Map<Eigen::Matrix<value_type, Eigen::Dynamic, 1>>(eigenvalues.data(), dimension) = subspaceEigenvalues;
541 // }
542 prof->stop();
543}
544
545template <typename value_type, typename std::enable_if_t<is_complex<value_type>{}, int>>
546void solve_LinearEquations(std::vector<value_type>& solution, std::vector<value_type>& eigenvalues,
547 const std::vector<value_type>& matrix, const std::vector<value_type>& metric,
548 const std::vector<value_type>& rhs, const size_t dimension, size_t nroot,
549 double augmented_hessian, double svdThreshold, int verbosity) {
550 assert(false); // Complex not implemented here
551}
552
553template <typename value_type, typename std::enable_if_t<!is_complex<value_type>{}, std::nullptr_t>>
554void solve_LinearEquations(std::vector<value_type>& solution, std::vector<value_type>& eigenvalues,
555 const std::vector<value_type>& matrix, const std::vector<value_type>& metric,
556 const std::vector<value_type>& rhs, const size_t dimension, size_t nroot,
557 double augmented_hessian, double svdThreshold, int verbosity) {
558 const Eigen::Index nX = dimension;
559 solution.resize(nX * nroot);
560 // std::cout << "augmented_hessian "<<augmented_hessian<<std::endl;
561 if (augmented_hessian > 0) { // Augmented hessian
562 Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic> subspaceMatrix;
563 Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic> subspaceOverlap;
564 subspaceMatrix.conservativeResize(nX + 1, nX + 1);
565 subspaceOverlap.conservativeResize(nX + 1, nX + 1);
566 subspaceMatrix.block(0, 0, nX, nX) =
567 Eigen::Map<const Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic>>(matrix.data(), nX, nX);
568 subspaceOverlap.block(0, 0, nX, nX) =
569 Eigen::Map<const Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic>>(metric.data(), nX, nX);
570 eigenvalues.resize(nroot);
571 for (size_t root = 0; root < nroot; root++) {
572 for (Eigen::Index i = 0; i < nX; i++) {
573 subspaceMatrix(i, nX) = subspaceMatrix(nX, i) = -augmented_hessian * rhs[i + nX * root];
574 subspaceOverlap(i, nX) = subspaceOverlap(nX, i) = 0;
575 }
576 subspaceMatrix(nX, nX) = 0;
577 subspaceOverlap(nX, nX) = 1;
578 // std::cout << "subspace augmented hessian subspaceMatrix\n"<<subspaceMatrix<<std::endl;
579 // std::cout << "subspace augmented hessian subspaceOverlap\n"<<subspaceOverlap<<std::endl;
580
581 Eigen::GeneralizedEigenSolver<Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic>> s(subspaceMatrix,
582 subspaceOverlap);
583 auto eval = s.eigenvalues();
584 auto evec = s.eigenvectors();
585 Eigen::Index imax = 0;
586 for (Eigen::Index i = 0; i < nX + 1; i++)
587 if (eval(i).real() < eval(imax).real())
588 imax = i;
589 eigenvalues[root] = eval(imax).real();
590 auto Solution = evec.col(imax).real().head(nX) / (augmented_hessian * evec.real()(nX, imax));
591 for (auto k = 0; k < nX; k++)
592 solution[k + nX * root] = Solution(k);
593 // std::cout << "subspace augmented hessian solution\n"<<Solution<<std::endl;
594 }
595 } else { // straight solution of linear equations
596 Eigen::Map<const Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>> subspaceMatrixR(
597 matrix.data(), nX, nX);
598 Eigen::Map<const Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>> RHS_R(rhs.data(), nX,
599 nroot);
600 Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic> subspaceMatrix = subspaceMatrixR;
601 Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic> RHS = RHS_R;
602 Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic> Solution;
603// std::cout << "solve_LinearEquations RHS_R\n"<<RHS_R<<std::endl;
604// for (size_t i=0; i<RHS_R.cols()*RHS_R.rows(); ++i)
605// std::cout << " "<<RHS_R.data()[i];
606// std::cout << std::endl;
607// std::cout << "solve_LinearEquations RHS\n"<<RHS<<std::endl;
608// for (size_t i=0; i<RHS.cols()*RHS.rows(); ++i)
609// std::cout << " "<<RHS.data()[i];
610// std::cout << std::endl;
611 Solution = subspaceMatrix.householderQr().solve(RHS);
612 // std::cout << "subspace linear equations solution\n"<<Solution<<std::endl;
613 for (size_t root = 0; root < nroot; root++)
614 for (auto k = 0; k < nX; k++)
615 solution[k + nX * root] = Solution(k, root);
616 }
617}
618
619template <typename value_type, typename std::enable_if_t<!is_complex<value_type>{}, std::nullptr_t>>
620void solve_DIIS(std::vector<value_type>& solution, const std::vector<value_type>& matrix, const size_t dimension,
621 double svdThreshold, int verbosity) {
622 auto nAug = dimension + 1;
623 // auto nQ = dimension - 1;
624 solution.resize(dimension);
625 // if (nQ > 0) {
626 Eigen::VectorXd Rhs(nAug), Coeffs(nAug);
627 Eigen::MatrixXd BAug(nAug, nAug);
628 // Eigen::Matrix<value_type, Eigen::Dynamic, 1> Rhs(nQ), Coeffs(nQ);
629 // Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic> B(nQ, nQ);
630 //
631 Eigen::Map<const Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic>> subspaceMatrix(matrix.data(), dimension,
632 dimension);
633 BAug.block(0, 0, dimension, dimension) = subspaceMatrix;
634 for (size_t i = 0; i < dimension; ++i) {
635 BAug(dimension, i) = BAug(i, dimension) = -1;
636 Rhs(i) = 0;
637 }
638 BAug(dimension, dimension) = 0;
639 Rhs(dimension) = -1;
640 //
641 // molpro::cout << "BAug:" << std::endl << BAug << std::endl;
642 // molpro::cout << "Rhs:" << std::endl << Rhs << std::endl;
643
644 // invert the system, determine extrapolation coefficients.
645 Eigen::JacobiSVD<Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic>> svd(BAug, Eigen::ComputeThinU |
646 Eigen::ComputeThinV);
647
648 // std::cout << "svd thresholds " << svdThreshold << "," << svd.singularValues().maxCoeff() << std::endl;
649 // std::cout << "singular values " << svd.singularValues().transpose() << std::endl;
650 svd.setThreshold(svdThreshold * svd.singularValues().maxCoeff() * 0);
651 // molpro::cout << "svdThreshold "<<svdThreshold<<std::endl;
652 // molpro::cout << "U\n"<<svd.matrixU()<<std::endl;
653 // molpro::cout << "V\n"<<svd.matrixV()<<std::endl;
654 // molpro::cout << "singularValues\n"<<svd.singularValues()<<std::endl;
655 Coeffs = svd.solve(Rhs).head(dimension);
656 // Coeffs = BAug.fullPivHouseholderQr().solve(Rhs);
657 // molpro::cout << "Coeffs "<<Coeffs.transpose()<<std::endl;
658 if (verbosity > 1)
659 molpro::cout << "Combination of iteration vectors: " << Coeffs.transpose() << std::endl;
660 for (size_t k = 0; k < (size_t)Coeffs.rows(); k++) {
661 if (std::isnan(std::abs(Coeffs(k)))) {
662 molpro::cout << "B:" << std::endl << BAug << std::endl;
663 molpro::cout << "Rhs:" << std::endl << Rhs << std::endl;
664 molpro::cout << "Combination of iteration vectors: " << Coeffs.transpose() << std::endl;
665 throw std::overflow_error("NaN detected in DIIS submatrix solution");
666 }
667 solution[k] = Coeffs(k);
668 }
669}
670} // namespace molpro::linalg::itsolv
671
672#endif // LINEARALGEBRA_SRC_MOLPRO_LINALG_ITERATIVESOLVER_HELPER_IMPLEMENTATION_H_
Non-owning container taking a pointer to the data buffer and its size and exposing routines for itera...
Definition: Span.h:28
bool empty() const
Definition: Span.h:76
iterator begin()
Definition: Span.h:66
size_type size() const
Definition: Span.h:74
iterator end()
Definition: Span.h:70
iterator data()
Definition: Span.h:63
static std::shared_ptr< Profiler > single()
4-parameter interpolation of a 1-dimensional function given two points for which function values and ...
Definition: helper.h:10
std::list< SVD< value_type > > svd_eigen_bdcsvd(size_t nrows, size_t ncols, const array::Span< value_type > &m, double threshold)
Definition: helper-implementation.h:35
std::list< SVD< value_type > > svd_system(size_t nrows, size_t ncols, const array::Span< value_type > &m, double threshold, bool hermitian=false, bool reduce_to_rank=false)
Performs singular value decomposition and returns SVD objects for singular values less than threshold...
Definition: helper-implementation.h:264
void solve_LinearEquations(std::vector< value_type > &solution, std::vector< value_type > &eigenvalues, const std::vector< value_type > &matrix, const std::vector< value_type > &metric, const std::vector< value_type > &rhs, size_t dimension, size_t nroot, double augmented_hessian, double svdThreshold, int verbosity)
Definition: helper-implementation.h:546
size_t get_rank(std::vector< value_type > eigenvalues, value_type threshold)
Definition: helper-implementation.h:222
void solve_DIIS(std::vector< value_type > &solution, const std::vector< value_type > &matrix, size_t dimension, double svdThreshold, int verbosity=0)
Definition: helper-implementation.h:620
void eigenproblem(std::vector< value_type > &eigenvectors, std::vector< value_type > &eigenvalues, const std::vector< value_type > &matrix, const std::vector< value_type > &metric, size_t dimension, bool hermitian, double svdThreshold, int verbosity, bool condone_complex)
Definition: helper-implementation.h:312
int eigensolver_lapacke_dsyev(const std::vector< double > &matrix, std::vector< double > &eigenvectors, std::vector< double > &eigenvalues, const size_t dimension)
Definition: helper-implementation.h:122
void printMatrix(const std::vector< value_type > &, size_t rows, size_t cols, std::string title="", std::ostream &s=molpro::cout)
Definition: helper-implementation.h:306
std::list< SVD< value_type > > svd_eigen_jacobi(size_t nrows, size_t ncols, const array::Span< value_type > &m, double threshold)
Definition: helper-implementation.h:13
Stores a singular value and corresponding left and right singular vectors.
Definition: helper.h:19
value_type value
Definition: helper.h:21