votca 2024.2-dev
Loading...
Searching...
No Matches
dftengine.h
Go to the documentation of this file.
1/*
2 * Copyright 2009-2023 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_DFTENGINE_H
22#define VOTCA_XTP_DFTENGINE_H
23
24// VOTCA includes
26
27// Local VOTCA includes
28#include "ERIs.h"
29#include "convergenceacc.h"
30#include "ecpaobasis.h"
31#include "logger.h"
32#include "qmmolecule.h"
33#include "staticsite.h"
34#include "vxc_grid.h"
35#include "vxc_potential.h"
36namespace votca {
37namespace xtp {
38class Orbitals;
39
48class DFTEngine {
49 public:
50 void Initialize(tools::Property& options);
51
52 void setLogger(Logger* pLog) { pLog_ = pLog; }
53
55 std::vector<std::unique_ptr<StaticSite> >* externalsites) {
56 externalsites_ = externalsites;
57 }
58
59 bool Evaluate(Orbitals& orb);
60
63
64 std::string getDFTBasisName() const { return dftbasis_name_; };
65
66 private:
67 void Prepare(Orbitals& orb, Index numofelectrons = -1);
68
70
71 Eigen::MatrixXd OrthogonalizeGuess(const Eigen::MatrixXd& GuessMOs) const;
72 void PrintMOs(const Eigen::VectorXd& MOEnergies, Log::Level level);
73 void CalcElDipole(const Orbitals& orb) const;
74
75 std::array<Eigen::MatrixXd, 2> CalcERIs_EXX(const Eigen::MatrixXd& MOCoeff,
76 const Eigen::MatrixXd& Dmat,
77 double error) const;
78
79 Eigen::MatrixXd CalcERIs(const Eigen::MatrixXd& Dmat, double error) const;
80
81 void ConfigOrbfile(Orbitals& orb);
83 Eigen::MatrixXd McWeenyPurification(Eigen::MatrixXd& Dmat_in,
84 AOOverlap& overlap);
85
86 Mat_p_Energy SetupH0(const QMMolecule& mol) const;
88 const QMMolecule& mol,
89 const std::vector<std::unique_ptr<StaticSite> >& multipoles) const;
91 const Orbitals& extdensity) const;
92
93 Eigen::MatrixXd IntegrateExternalField(const QMMolecule& mol) const;
94
97 const Mat_p_Energy& H0, const QMMolecule& mol,
98 const Vxc_Potential<Vxc_Grid>& vxcpotential) const;
99
100 Eigen::MatrixXd AtomicGuess(const QMMolecule& mol) const;
101
102 Eigen::MatrixXd RunAtomicDFT_unrestricted(const QMAtom& uniqueAtom) const;
103
104 double NuclearRepulsion(const QMMolecule& mol) const;
105 double ExternalRepulsion(
106 const QMMolecule& mol,
107 const std::vector<std::unique_ptr<StaticSite> >& multipoles) const;
108 Eigen::MatrixXd SphericalAverageShells(const Eigen::MatrixXd& dmat,
109 const AOBasis& dftbasis) const;
110
111 void TruncateBasis(Orbitals& orb, std::vector<Index>& activeatoms,
112 Mat_p_Energy& H0,
113 Eigen::MatrixXd InitialActiveDensityMatrix,
114 Eigen::MatrixXd v_embedding,
115 Eigen::MatrixXd InitialInactiveMOs);
116
117 void TruncMOsFullBasis(Orbitals& orb, std::vector<Index> activeatoms,
118 std::vector<Index> numfuncpatom);
119 Eigen::MatrixXd InsertZeroCols(Eigen::MatrixXd MOsMatrix, Index startidx,
120 Index numofzerocols);
121 Eigen::MatrixXd InsertZeroRows(Eigen::MatrixXd MOsMatrix, Index startidx,
122 Index numofzerorows);
124
125 // basis sets
126 std::string auxbasis_name_;
127 std::string dftbasis_name_;
128 std::string ecp_name_;
132
134 // Pre-screening
136
137 // numerical integration Vxc
138 std::string grid_name_;
139
140 // AO Matrices
142
143 std::string initial_guess_;
144
145 // Convergence
149 // DIIS variables
151 // Electron repulsion integrals
153
154 // external charges
155 std::vector<std::unique_ptr<StaticSite> >* externalsites_ = nullptr;
156
157 // exchange and correlation
158 double ScaHFX_;
160
162 // integrate external density
163 std::string orbfilename_;
164 std::string gridquality_;
165 std::string state_;
166
168 QMMolecule("molecule made of atoms participating in Active region", 1);
169
170 Eigen::Vector3d extfield_ = Eigen::Vector3d::Zero();
172
176
177 // truncation
178 Eigen::MatrixXd H0_trunc_;
180 Eigen::MatrixXd v_embedding_trunc_;
184 double E_nuc_;
186 std::vector<Index> active_and_border_atoms_;
187 std::vector<Index> numfuncpatom_;
188};
189
190} // namespace xtp
191} // namespace votca
192
193#endif // VOTCA_XTP_DFTENGINE_H
class to manage program options with xml serialization functionality
Definition property.h:55
Container to hold Basisfunctions for all atoms.
Definition aobasis.h:42
Electronic ground-state via Density-Functional Theory.
Definition dftengine.h:48
Eigen::MatrixXd CalcERIs(const Eigen::MatrixXd &Dmat, double error) const
Definition dftengine.cc:174
void setExternalcharges(std::vector< std::unique_ptr< StaticSite > > *externalsites)
Definition dftengine.h:54
bool EvaluateTruncatedActiveRegion(Orbitals &trunc_orb)
std::string auxbasis_name_
Definition dftengine.h:126
std::string getDFTBasisName() const
Definition dftengine.h:64
std::string gridquality_
Definition dftengine.h:164
tools::EigenSystem ModelPotentialGuess(const Mat_p_Energy &H0, const QMMolecule &mol, const Vxc_Potential< Vxc_Grid > &vxcpotential) const
Definition dftengine.cc:188
Mat_p_Energy IntegrateExternalMultipoles(const QMMolecule &mol, const std::vector< std::unique_ptr< StaticSite > > &multipoles) const
Definition dftengine.cc:988
tools::EigenSystem IndependentElectronGuess(const Mat_p_Energy &H0) const
Definition dftengine.cc:183
Eigen::MatrixXd InitialActiveDmat_trunc_
Definition dftengine.h:179
void TruncMOsFullBasis(Orbitals &orb, std::vector< Index > activeatoms, std::vector< Index > numfuncpatom)
void setLogger(Logger *pLog)
Definition dftengine.h:52
void TruncateBasis(Orbitals &orb, std::vector< Index > &activeatoms, Mat_p_Energy &H0, Eigen::MatrixXd InitialActiveDensityMatrix, Eigen::MatrixXd v_embedding, Eigen::MatrixXd InitialInactiveMOs)
std::string dftbasis_name_
Definition dftengine.h:127
Mat_p_Energy IntegrateExternalDensity(const QMMolecule &mol, const Orbitals &extdensity) const
std::string ecp_name_
Definition dftengine.h:128
std::string grid_name_
Definition dftengine.h:138
bool EvaluateActiveRegion(Orbitals &orb)
Eigen::MatrixXd RunAtomicDFT_unrestricted(const QMAtom &uniqueAtom) const
Definition dftengine.cc:509
void Prepare(Orbitals &orb, Index numofelectrons=-1)
Definition dftengine.cc:801
std::string initial_guess_
Definition dftengine.h:143
std::string orbfilename_
Definition dftengine.h:163
Eigen::MatrixXd AtomicGuess(const QMMolecule &mol) const
Definition dftengine.cc:697
Eigen::MatrixXd InsertZeroRows(Eigen::MatrixXd MOsMatrix, Index startidx, Index numofzerorows)
Eigen::MatrixXd IntegrateExternalField(const QMMolecule &mol) const
Definition dftengine.cc:975
Eigen::Vector3d extfield_
Definition dftengine.h:170
std::array< Eigen::MatrixXd, 2 > CalcERIs_EXX(const Eigen::MatrixXd &MOCoeff, const Eigen::MatrixXd &Dmat, double error) const
Definition dftengine.cc:159
Eigen::MatrixXd H0_trunc_
Definition dftengine.h:178
void Initialize(tools::Property &options)
Definition dftengine.cc:46
std::string xc_functional_name_
Definition dftengine.h:159
void CalcElDipole(const Orbitals &orb) const
Definition dftengine.cc:150
double ExternalRepulsion(const QMMolecule &mol, const std::vector< std::unique_ptr< StaticSite > > &multipoles) const
Definition dftengine.cc:949
ConvergenceAcc::options conv_opt_
Definition dftengine.h:148
std::string active_atoms_as_string_
Definition dftengine.h:173
void ConfigOrbfile(Orbitals &orb)
Definition dftengine.cc:738
QMMolecule activemol_
Definition dftengine.h:167
void PrintMOs(const Eigen::VectorXd &MOEnergies, Log::Level level)
Definition dftengine.cc:134
std::vector< Index > numfuncpatom_
Definition dftengine.h:187
AOOverlap dftAOoverlap_
Definition dftengine.h:141
std::vector< Index > active_and_border_atoms_
Definition dftengine.h:186
Mat_p_Energy SetupH0(const QMMolecule &mol) const
Definition dftengine.cc:369
Eigen::MatrixXd OrthogonalizeGuess(const Eigen::MatrixXd &GuessMOs) const
Eigen::MatrixXd v_embedding_trunc_
Definition dftengine.h:180
Eigen::MatrixXd SphericalAverageShells(const Eigen::MatrixXd &dmat, const AOBasis &dftbasis) const
Definition dftengine.cc:918
bool Evaluate(Orbitals &orb)
Definition dftengine.cc:209
std::vector< std::unique_ptr< StaticSite > > * externalsites_
Definition dftengine.h:155
Eigen::MatrixXd McWeenyPurification(Eigen::MatrixXd &Dmat_in, AOOverlap &overlap)
Vxc_Potential< Vxc_Grid > SetupVxc(const QMMolecule &mol)
Definition dftengine.cc:881
double NuclearRepulsion(const QMMolecule &mol) const
Definition dftengine.cc:902
ConvergenceAcc conv_accelerator_
Definition dftengine.h:150
Eigen::MatrixXd InsertZeroCols(Eigen::MatrixXd MOsMatrix, Index startidx, Index numofzerocols)
Container to hold ECPs for all atoms.
Definition ecpaobasis.h:43
Takes a density matrix and and an auxiliary basis set and calculates the electron repulsion integrals...
Definition ERIs.h:35
Logger is used for thread-safe output of messages.
Definition logger.h:164
container for molecular orbitals
Definition orbitals.h:46
container for QM atoms
Definition qmatom.h:37
base class for all analysis tools
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26
Level
be loud and noisy
Definition globals.h:28