1#ifndef LINEARALGEBRA_SRC_MOLPRO_LINALG_ITERATIVESOLVER_HELPER_IMPLEMENTATION_H_
2#define LINEARALGEBRA_SRC_MOLPRO_LINALG_ITERATIVESOLVER_HELPER_IMPLEMENTATION_H_
6#include <molpro/Profiler.h>
7#include <molpro/lapacke.h>
8#include <molpro/linalg/itsolv/helper.h>
12template <
typename value_type>
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) {
25 for (
size_t j = 0; j < ncols; ++j) {
26 t.v.emplace_back(svd.matrixV()(j, i));
34template <
typename value_type>
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) {
46 for (
size_t j = 0; j < ncols; ++j) {
47 t.v.emplace_back(svd.matrixV()(j, i));
56template <
typename value_type>
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>{};
72 for (
size_t j = 0; j < ncols; ++j) {
73 t.v.emplace_back(v[i * ncols + j]);
81template <
typename value_type>
82std::list<SVD<value_type>> svd_lapacke_dgesvd(
size_t nrows,
size_t ncols,
const array::Span<value_type>& mat,
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>{};
98 for (
size_t j = 0; j < ncols; ++j) {
99 t.v.emplace_back(v[i * ncols + j]);
110extern "C" int dsyev_c(
char,
char,
int,
double*,
int,
double*);
123 std::vector<double>& eigenvalues,
const size_t dimension) {
126 if (eigenvectors.size() != matrix.size()) {
127 throw std::runtime_error(
"Matrix of eigenvectors and input matrix are not the same size!");
130 if (eigenvectors.size() != dimension * dimension || eigenvalues.size() != dimension) {
131 throw std::runtime_error(
"Size of eigenvectors/eigenvlaues do not match dimension!");
136 static const char compute_eigenvalues_eigenvectors =
'V';
138 static const char store_lower_triangle =
'L';
141 std::copy(matrix.begin(), matrix.end(), eigenvectors.begin());
145 lapack_int leading_dimension = dimension;
146 lapack_int order = dimension;
150 status = dsyev_c(compute_eigenvalues_eigenvectors, store_lower_triangle, order, eigenvectors.data(),
151 leading_dimension, eigenvalues.data());
153 status = LAPACKE_dsyev(LAPACK_COL_MAJOR, compute_eigenvalues_eigenvectors, store_lower_triangle, order,
154 eigenvectors.data(), leading_dimension, eigenvalues.data());
168 std::vector<double> eigvecs(dimension * dimension);
169 std::vector<double> eigvals(dimension);
174 throw std::invalid_argument(
"Invalid argument of eigensolver_lapacke_dsyev: ");
177 throw std::runtime_error(
"Lapacke_dsyev (eigensolver) failed to converge. "
178 " elements of an intermediate tridiagonal form did not converge to zero.");
181 auto eigensystem = std::list<SVD<double>>{};
184 for (
int i = dimension - 1; i >= 0;
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)]);
191 eigensystem.emplace_back(temp_eigenproblem);
209 std::vector<double> v;
210 v.insert(v.begin(), matrix.
begin(), matrix.
end());
221template <
typename value_type>
222size_t get_rank(std::vector<value_type> eigenvalues, value_type threshold) {
223 if (eigenvalues.size() == 0) {
226 value_type max = *max_element(eigenvalues.begin(), eigenvalues.end());
227 value_type threshold_scaled = threshold * max;
229 std::count_if(eigenvalues.begin(), eigenvalues.end(), [&](
auto const& val) { return val >= threshold_scaled; });
240template <
typename value_type>
243 value_type max_value = 0;
244 std::list<SVD<double>>::iterator it;
246 if (it->value > max_value) {
247 max_value = it->value;
251 value_type threshold_scaled = threshold * max_value;
256 if (it->value > threshold_scaled) {
263template <
typename value_type,
typename std::enable_if_t<!is_complex<value_type>{}, std::
nullptr_t>>
265 bool hermitian,
bool reduce_to_rank) {
266 std::list<SVD<value_type>> svds;
267 assert(m.
size() == nrows * ncols);
271 assert(nrows == ncols);
273 for (
auto s = svds.begin(); s != svds.end();)
274 if (s->value > threshold)
284 svds = svd_eigen_jacobi<value_type>(nrows, ncols, m, threshold);
289 if (reduce_to_rank) {
290 int rank =
get_rank(svds, threshold);
291 for (
int i = ncols; i > rank; i--) {
298template <
typename value_type,
typename std::enable_if_t<is_complex<value_type>{},
int>>
300 bool hermitian,
bool reduce_to_rank) {
305template <
typename value_type>
306void printMatrix(
const std::vector<value_type>& m,
size_t rows,
size_t cols, std::string title, std::ostream& s) {
308 << Eigen::Map<const Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic>>(m.data(), rows, cols) << std::endl;
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) {
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) {
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);
328 Eigen::Map<const Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic>> S(metric.data(), dimension, dimension);
329 Eigen::MatrixXcd subspaceEigenvectors;
330 Eigen::VectorXcd subspaceEigenvalues;
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;
343 eigvecs.resize(dimension * dimension);
344 eigvals.resize(dimension);
347 throw std::runtime_error(
"Eigensolver did not converge");
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);
356 Eigen::JacobiSVD<Eigen::MatrixXd> svd(S, Eigen::ComputeThinU | Eigen::ComputeThinV);
357 singularValues = svd.singularValues();
358 matrixV = svd.matrixV();
359 matrixU = svd.matrixU();
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;
374 (svmh.asDiagonal()) * (matrixU.leftCols(rank).adjoint()) * H * matrixV.leftCols(rank) * (svmh.asDiagonal());
382 Eigen::EigenSolver<Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic>> s(Hbar);
384 subspaceEigenvalues = s.eigenvalues();
385 if (s.eigenvalues().imag().norm() < 1e-10) {
387 subspaceEigenvalues = subspaceEigenvalues.real();
388 subspaceEigenvectors = s.eigenvectors();
391 for (
int i = 0; i < subspaceEigenvectors.cols(); i++) {
392 if (subspaceEigenvectors.col(i).imag().norm() > 1e-10) {
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();
401 subspaceEigenvectors = matrixV.leftCols(rank) * svmh.asDiagonal() * subspaceEigenvectors;
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");
411 subspaceEigenvectors = matrixV.leftCols(rank) * svmh.asDiagonal() * s.eigenvectors();
417 auto eigval = subspaceEigenvalues;
418 auto eigvec = subspaceEigenvectors;
419 std::vector<Eigen::Index> map;
420 for (Eigen::Index k = 0; k < Hbar.cols(); k++) {
422 for (ll = 0; std::count(map.begin(), map.end(), ll) != 0; ll++)
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())
431 subspaceEigenvalues(k) = eigval(ll);
434 subspaceEigenvectors.col(k) = eigvec.col(ll);
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()))
440 if (subspaceEigenvectors.col(k)[maxcomp].real() < 0)
441 subspaceEigenvectors.col(k) = - subspaceEigenvectors.col(k);
452 Eigen::MatrixXcd ovlTimesVec(subspaceEigenvectors.cols(), subspaceEigenvectors.rows());
453 for (
auto repeat = 0; repeat < 3; ++repeat)
454 for (Eigen::Index k = 0; k < subspaceEigenvectors.cols(); k++) {
455 if (std::abs(subspaceEigenvalues(k)) <
457 subspaceEigenvectors.col(k).real() += double(0.3256897) * subspaceEigenvectors.col(k).imag();
458 subspaceEigenvectors.col(k).imag().setZero();
461 for (Eigen::Index l = 0; l < k; l++) {
473 subspaceEigenvectors.col(k) -= subspaceEigenvectors.col(l) *
474 ovlTimesVec.row(l).dot(subspaceEigenvectors.col(k));
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;
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)))
499 if (subspaceEigenvectors(lmax, k).real() < 0)
500 subspaceEigenvectors.col(k) = -subspaceEigenvectors.col(k);
511 if (condone_complex) {
512 for (Eigen::Index root = 0; root < Hbar.cols(); ++root) {
513 if (subspaceEigenvalues(root).imag() != 0) {
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();
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");
528 eigenvectors.resize(dimension * Hbar.cols());
529 eigenvalues.resize(Hbar.cols());
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();
545template <
typename value_type,
typename std::enable_if_t<is_complex<value_type>{},
int>>
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) {
553template <
typename value_type,
typename std::enable_if_t<!is_complex<value_type>{}, std::
nullptr_t>>
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);
561 if (augmented_hessian > 0) {
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;
576 subspaceMatrix(nX, nX) = 0;
577 subspaceOverlap(nX, nX) = 1;
581 Eigen::GeneralizedEigenSolver<Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic>> s(subspaceMatrix,
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())
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);
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,
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;
611 Solution = subspaceMatrix.householderQr().solve(RHS);
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);
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;
624 solution.resize(dimension);
626 Eigen::VectorXd Rhs(nAug), Coeffs(nAug);
627 Eigen::MatrixXd BAug(nAug, nAug);
631 Eigen::Map<const Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic>> subspaceMatrix(matrix.data(), 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;
638 BAug(dimension, dimension) = 0;
645 Eigen::JacobiSVD<Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic>> svd(BAug, Eigen::ComputeThinU |
646 Eigen::ComputeThinV);
650 svd.setThreshold(svdThreshold * svd.singularValues().maxCoeff() * 0);
655 Coeffs = svd.solve(Rhs).head(dimension);
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");
667 solution[k] = Coeffs(k);
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