votca 2024.2-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
37 public:
39
55
57 opt_ = opt;
58 if (opt_.mode == KSmode::closed) {
60 } else if (opt_.mode == KSmode::open) {
62 } else if (opt_.mode == KSmode::fractional) {
63 nocclevels_ = 0;
64 }
66 }
67 void setLogger(Logger* log) { log_ = log; }
68
69 void PrintConfigOptions() const;
70
71 bool isConverged() const {
72 if (totE_.size() < 2) {
73 return false;
74 } else {
75 return std::abs(getDeltaE()) < opt_.Econverged &&
77 }
78 }
79
80 double getDeltaE() const {
81 if (totE_.size() < 2) {
82 return 0;
83 } else {
84 return totE_.back() - totE_[totE_.size() - 2];
85 }
86 }
87 void setOverlap(AOOverlap& S, double etol);
88
89 double getDIIsError() const { return diiserror_; }
90
91 bool getUseMixing() const { return usedmixing_; }
92
93 Eigen::MatrixXd Iterate(const Eigen::MatrixXd& dmat, Eigen::MatrixXd& H,
94 tools::EigenSystem& MOs, double totE);
95 tools::EigenSystem SolveFockmatrix(const Eigen::MatrixXd& H) const;
96 void Levelshift(Eigen::MatrixXd& H, const Eigen::MatrixXd& MOs_old) const;
97
98 Eigen::MatrixXd DensityMatrix(const tools::EigenSystem& MOs) const;
99
100 private:
102
103 Eigen::MatrixXd DensityMatrixGroundState(const Eigen::MatrixXd& MOs) const;
104 Eigen::MatrixXd DensityMatrixGroundState_unres(
105 const Eigen::MatrixXd& MOs) const;
106 Eigen::MatrixXd DensityMatrixGroundState_frac(
107 const tools::EigenSystem& MOs) const;
108
109 bool usedmixing_ = true;
110 double diiserror_ = std::numeric_limits<double>::max();
112 const AOOverlap* S_;
113
114 Eigen::MatrixXd Sminusahalf;
115 std::vector<Eigen::MatrixXd> mathist_;
116 std::vector<Eigen::MatrixXd> dmatHist_;
117 std::vector<double> totE_;
118
121 double maxerror_ = 0.0;
124};
125
126} // namespace xtp
127} // namespace votca
128
129#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)
Eigen::MatrixXd Iterate(const Eigen::MatrixXd &dmat, Eigen::MatrixXd &H, tools::EigenSystem &MOs, double totE)
tools::EigenSystem SolveFockmatrix(const Eigen::MatrixXd &H) const
void Levelshift(Eigen::MatrixXd &H, const Eigen::MatrixXd &MOs_old) const
Eigen::MatrixXd DensityMatrixGroundState_frac(const tools::EigenSystem &MOs) const
std::vector< double > totE_
void setOverlap(AOOverlap &S, double etol)
Eigen::MatrixXd DensityMatrixGroundState(const Eigen::MatrixXd &MOs) const
std::vector< Eigen::MatrixXd > dmatHist_
void setHistLength(Index length)
Definition diis.h:38
Logger is used for thread-safe output of messages.
Definition logger.h:164
base class for all analysis tools
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26