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