votca 2026-dev
Loading...
Searching...
No Matches
bse_initialization.h
Go to the documentation of this file.
1/*
2 * Copyright 2009-2026 The VOTCA Development Team
3 * (http://www.votca.org)
4 *
5 * Licensed under the Apache License, Version 2.0 (the "License")
6 *
7 * You may not use this file except in compliance with the License.
8 * You may obtain a copy of the License at
9 *
10 * http://www.apache.org/licenses/LICENSE-2.0
11 *
12 * Unless required by applicable law or agreed to in writing, software
13 * distributed under the License is distributed on an "AS IS" BASIS,
14 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
15 * See the License for the specific language governing permissions and
16 * limitations under the License.
17 *
18 */
19
20#pragma once
21
22#include <Eigen/Dense>
23#include <algorithm>
24#include <cmath>
25#include <stdexcept>
26#include <vector>
27
28namespace votca {
29namespace xtp {
30
47inline Eigen::MatrixXd BuildFullBSEXRankedInitialGuess(
48 const Eigen::VectorXd& adiag, const Eigen::VectorXd& bdiag, Index nroots) {
49
50 if (adiag.size() != bdiag.size()) {
51 throw std::runtime_error("BuildFullBSEXRankedInitialGuess: size mismatch.");
52 }
53
54 const Index n = adiag.size();
55 const Index nguess = std::min<Index>(n, std::max<Index>(4 * nroots, 8));
56
57 struct RankedMode {
58 double omega;
59 double a;
60 Index idx;
61 };
62
63 std::vector<RankedMode> ranked;
64 ranked.reserve(n);
65
66 for (Index i = 0; i < n; ++i) {
67 const double a = adiag(i);
68 const double b = bdiag(i);
69
70 const double disc = std::max(0.0, (a - b) * (a + b));
71 const double omega = std::sqrt(disc);
72
73 ranked.push_back({omega, a, i});
74 }
75
76 std::sort(ranked.begin(), ranked.end(),
77 [](const RankedMode& lhs, const RankedMode& rhs) {
78 if (lhs.omega != rhs.omega) {
79 return lhs.omega < rhs.omega;
80 }
81 return lhs.a < rhs.a;
82 });
83
84 Eigen::MatrixXd guess = Eigen::MatrixXd::Zero(2 * n, nguess);
85
86 // X-only basis vectors
87 for (Index col = 0; col < nguess; ++col) {
88 const Index i = ranked[col].idx;
89 guess(i, col) = 1.0;
90 }
91
92 return guess;
93}
94
95} // namespace xtp
96} // namespace votca
Charge transport classes.
Definition ERIs.h:28
Eigen::MatrixXd BuildFullBSEXRankedInitialGuess(const Eigen::VectorXd &adiag, const Eigen::VectorXd &bdiag, Index nroots)
Build an initial guess for full-BSE Davidson diagonalization.
Provides a means for comparing floating point numbers.
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26