votca 2026-dev
Loading...
Searching...
No Matches
uks_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
22#ifndef VOTCA_XTP_UKS_CONVERGENCEACC_H
23#define VOTCA_XTP_UKS_CONVERGENCEACC_H
24
25// VOTCA includes
26#include <votca/tools/linalg.h>
27
28#include "votca/xtp/adiis.h"
29#include "votca/xtp/aomatrix.h"
31#include "votca/xtp/diis.h"
32#include "votca/xtp/logger.h"
33
34namespace votca {
35namespace xtp {
36
38 public:
41
42 struct SpinDensity {
43 Eigen::MatrixXd alpha;
44 Eigen::MatrixXd beta;
45
46 Eigen::MatrixXd total() const { return alpha + beta; }
47 };
48
49 struct SpinFock {
50 Eigen::MatrixXd alpha;
51 Eigen::MatrixXd beta;
52 };
53
54 void Configure(const options& opt_alpha, const options& opt_beta);
55 void setLogger(Logger* log);
56 void setOverlap(AOOverlap& S, double etol);
57
59 const tools::EigenSystem& MOs_beta) const;
60
62 tools::EigenSystem& MOs_alpha,
63 tools::EigenSystem& MOs_beta, double totE);
64
65 tools::EigenSystem SolveFockmatrix(const Eigen::MatrixXd& H) const;
66
67 bool isConverged() const;
68 double getDIIsError() const { return diiserror_; }
69 double getDeltaE() const;
70 bool getUseMixing() const { return usedmixing_; }
71
72 private:
73 Eigen::MatrixXd DensityMatrixGroundState_unres(const Eigen::MatrixXd& MOs,
74 Index nocclevels) const;
75
76 void Levelshift(Eigen::MatrixXd& H, const Eigen::MatrixXd& MOs_old,
77 const options& opt, Index nocclevels) const;
78
79 Eigen::MatrixXd BuildErrorMatrix(const Eigen::MatrixXd& dmat,
80 const Eigen::MatrixXd& H) const;
81
82 double CombinedError(const Eigen::MatrixXd& err_alpha,
83 const Eigen::MatrixXd& err_beta) const;
84
87
88 AOOverlap* S_ = nullptr;
89 Logger* log_ = nullptr;
90 Eigen::MatrixXd Sminusahalf;
91
94
95 std::vector<Eigen::MatrixXd> mathist_alpha_;
96 std::vector<Eigen::MatrixXd> mathist_beta_;
97 std::vector<Eigen::MatrixXd> dmatHist_alpha_;
98 std::vector<Eigen::MatrixXd> dmatHist_beta_;
99 std::vector<double> totE_;
100
103
104 double diiserror_ = 1.0;
105 double maxerror_ = -1.0;
107 bool usedmixing_ = true;
108};
109
110} // namespace xtp
111} // namespace votca
112
113#endif // VOTCA_XTP_UKS_CONVERGENCEACC_H
KSmode
Occupation model used when constructing density matrices.
Logger is used for thread-safe output of messages.
Definition logger.h:164
double CombinedError(const Eigen::MatrixXd &err_alpha, const Eigen::MatrixXd &err_beta) const
tools::EigenSystem SolveFockmatrix(const Eigen::MatrixXd &H) const
ConvergenceAcc::KSmode KSmode
SpinDensity DensityMatrix(const tools::EigenSystem &MOs_alpha, const tools::EigenSystem &MOs_beta) const
void Levelshift(Eigen::MatrixXd &H, const Eigen::MatrixXd &MOs_old, const options &opt, Index nocclevels) const
std::vector< Eigen::MatrixXd > mathist_beta_
void setOverlap(AOOverlap &S, double etol)
void Configure(const options &opt_alpha, const options &opt_beta)
Eigen::MatrixXd BuildErrorMatrix(const Eigen::MatrixXd &dmat, const Eigen::MatrixXd &H) const
ConvergenceAcc::options options
SpinDensity Iterate(const SpinDensity &dmat, SpinFock &H, tools::EigenSystem &MOs_alpha, tools::EigenSystem &MOs_beta, double totE)
std::vector< Eigen::MatrixXd > dmatHist_alpha_
std::vector< Eigen::MatrixXd > dmatHist_beta_
Eigen::MatrixXd DensityMatrixGroundState_unres(const Eigen::MatrixXd &MOs, Index nocclevels) const
std::vector< Eigen::MatrixXd > mathist_alpha_
Provides a means for comparing floating point numbers.
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26