votca 2026-dev
Loading...
Searching...
No Matches
convergenceacc.h
Go to the documentation of this file.
1/*
2 * Copyright 2009-2020 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#ifndef VOTCA_XTP_CONVERGENCEACC_H
22#define VOTCA_XTP_CONVERGENCEACC_H
23
24// VOTCA includes
25#include <votca/tools/linalg.h>
26
27// Local VOTCA includes
28#include "adiis.h"
29#include "aomatrix.h"
30#include "diis.h"
31#include "logger.h"
32
33namespace votca {
34namespace xtp {
35
44 public:
47
67
69 struct SpinDensity {
70 Eigen::MatrixXd alpha;
71 Eigen::MatrixXd beta;
72
74 Eigen::MatrixXd total() const { return alpha + beta; }
75
77 Eigen::MatrixXd spin() const { return alpha - beta; }
78 };
79
83 opt_ = opt;
84 if (opt_.mode == KSmode::closed) {
85 nocclevels_ = opt_.numberofelectrons / 2;
86 } else if (opt_.mode == KSmode::open) {
87 nocclevels_ = opt_.numberofelectrons;
88 } else if (opt_.mode == KSmode::fractional) {
89 nocclevels_ = 0;
90 } else if (opt_.mode == KSmode::restricted_open) {
92 std::max(opt_.number_alpha_electrons, opt_.number_beta_electrons);
93 }
94 diis_.setHistLength(opt_.histlength);
95 }
96
97 void setLogger(Logger* log) { log_ = log; }
98
100 void PrintConfigOptions() const;
101
104 bool isConverged() const {
105 if (totE_.size() < 2) {
106 return false;
107 } else {
108 return std::abs(getDeltaE()) < opt_.Econverged &&
109 getDIIsError() < opt_.error_converged;
110 }
111 }
112
114 double getDeltaE() const {
115 if (totE_.size() < 2) {
116 return 0;
117 } else {
118 return totE_.back() - totE_[totE_.size() - 2];
119 }
120 }
121
122 void setOverlap(AOOverlap& S, double etol);
123
125 double getDIIsError() const { return diiserror_; }
126
129 bool getUseMixing() const { return usedmixing_; }
130
133 Eigen::MatrixXd Iterate(const Eigen::MatrixXd& dmat, Eigen::MatrixXd& H,
134 tools::EigenSystem& MOs, double totE);
136 tools::EigenSystem SolveFockmatrix(const Eigen::MatrixXd& H) const;
138 void Levelshift(Eigen::MatrixXd& H, const Eigen::MatrixXd& MOs_old) const;
139
142 Eigen::MatrixXd DensityMatrix(const tools::EigenSystem& MOs) const;
143
146 SpinDensity DensityMatrixSpinResolved(const tools::EigenSystem& MOs) const;
147
148 private:
150
153 Eigen::MatrixXd DensityMatrixGroundState(const Eigen::MatrixXd& MOs) const;
156 Eigen::MatrixXd DensityMatrixGroundState_unres(
157 const Eigen::MatrixXd& MOs) const;
159 Eigen::MatrixXd DensityMatrixGroundState_frac(
160 const tools::EigenSystem& MOs) const;
161
165 const Eigen::MatrixXd& MOs) const;
166
167 bool usedmixing_ = true;
168 double diiserror_ = std::numeric_limits<double>::max();
170 const AOOverlap* S_;
171
172 Eigen::MatrixXd Sminusahalf;
173 std::vector<Eigen::MatrixXd> mathist_;
174 std::vector<Eigen::MatrixXd> dmatHist_;
175 std::vector<double> totE_;
176
179 double maxerror_ = 0.0;
182};
183
184} // namespace xtp
185} // namespace votca
186
187#endif // VOTCA_XTP_CONVERGENCEACC_H
Eigen::MatrixXd DensityMatrix(const tools::EigenSystem &MOs) const
void Configure(const ConvergenceAcc::options &opt)
Eigen::MatrixXd DensityMatrixGroundState_unres(const Eigen::MatrixXd &MOs) const
std::vector< Eigen::MatrixXd > mathist_
void setLogger(Logger *log)
Attach the logger used for convergence diagnostics.
Eigen::MatrixXd Iterate(const Eigen::MatrixXd &dmat, Eigen::MatrixXd &H, tools::EigenSystem &MOs, double totE)
double getDIIsError() const
Return the DIIS commutator norm from the latest iteration.
tools::EigenSystem SolveFockmatrix(const Eigen::MatrixXd &H) const
Solve the generalized eigenvalue problem for the current Fock matrix.
void PrintConfigOptions() const
Print the active convergence-acceleration settings to the logger.
void Levelshift(Eigen::MatrixXd &H, const Eigen::MatrixXd &MOs_old) const
Apply a virtual-space level shift in the molecular-orbital basis.
double getDeltaE() const
Return the total-energy change between the two most recent SCF iterations.
Eigen::MatrixXd DensityMatrixGroundState_frac(const tools::EigenSystem &MOs) const
Construct a fractional-occupation density matrix from orbital occupations.
std::vector< double > totE_
void setOverlap(AOOverlap &S, double etol)
Precompute overlap-dependent quantities used when solving the Fock matrix.
Eigen::MatrixXd DensityMatrixGroundState(const Eigen::MatrixXd &MOs) const
std::vector< Eigen::MatrixXd > dmatHist_
KSmode
Occupation model used when constructing density matrices.
SpinDensity DensityMatrixSpinResolved(const tools::EigenSystem &MOs) const
SpinDensity DensityMatrixGroundState_restricted_open(const Eigen::MatrixXd &MOs) const
Logger is used for thread-safe output of messages.
Definition logger.h:164
Provides a means for comparing floating point numbers.
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26
Spin-resolved density matrices returned for open-shell SCF updates.
Eigen::MatrixXd total() const
Return the total density P = P^alpha + P^beta.
Eigen::MatrixXd spin() const
Return the spin density P^alpha - P^beta.