votca 2024.1-dev
Loading...
Searching...
No Matches
aopotential.h
Go to the documentation of this file.
1/*
2 * Copyright 2009-2021 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_AOPOTENTIAL_H
22#define VOTCA_XTP_AOPOTENTIAL_H
23
24// Local VOTCA includes
25#include "aobasis.h"
26#include "ecpaobasis.h"
27#include "staticsite.h"
28
29namespace votca {
30namespace xtp {
31
32class QMMolecule;
33
34// base class for 1D atomic orbital matrix types (overlap, Coulomb, ESP)
35template <class T>
37 public:
38 Index Dimension() { return aopotential_.rows(); }
39 const Eigen::Matrix<T, Eigen::Dynamic, Eigen::Dynamic>& Matrix() const {
40 return aopotential_;
41 }
42
43 protected:
44 Eigen::Matrix<T, Eigen::Dynamic, Eigen::Dynamic> Fill(
45 const AOBasis& aobasis) const;
46 virtual void FillBlock(
47 Eigen::Block<Eigen::Matrix<T, Eigen::Dynamic, Eigen::Dynamic>>& matrix,
48 const AOShell& shell_row, const AOShell& shell_col) const = 0;
49 Eigen::Matrix<T, Eigen::Dynamic, Eigen::Dynamic> aopotential_;
50};
51
52// derived class for Effective Core Potentials
53class AOECP : public AOPotential<double> {
54 public:
55 void FillPotential(const AOBasis& aobasis, const ECPAOBasis& ecp);
56
57 protected:
58 void FillBlock(Eigen::Block<Eigen::MatrixXd>&, const AOShell&,
59 const AOShell&) const final {};
60
61 private:
62};
63
64class AOMultipole : public AOPotential<double> {
65 public:
66 void FillPotential(const AOBasis& aobasis, const QMMolecule& atoms);
67 void FillPotential(const AOBasis& aobasis, const Eigen::Vector3d& r);
68 void FillPotential(
69 const AOBasis& aobasis,
70 const std::vector<std::unique_ptr<StaticSite>>& externalsites);
71
72 protected:
73 void FillBlock(Eigen::Block<Eigen::MatrixXd>& matrix,
74 const AOShell& shell_row,
75 const AOShell& shell_col) const override;
76
77 private:
78 void setSite(const StaticSite* site) { site_ = site; };
79
81};
82
83class AOPlanewave : public AOPotential<std::complex<double>> {
84 public:
85 void FillPotential(const AOBasis& aobasis,
86 const std::vector<Eigen::Vector3d>& kpoints);
87
88 protected:
89 void FillBlock(Eigen::Block<Eigen::MatrixXcd>& matrix,
90 const AOShell& shell_row,
91 const AOShell& shell_col) const override;
92
93 private:
94 void setkVector(const Eigen::Vector3d& k) { k_ = k; };
95 Eigen::Vector3d k_;
96};
97
98} // namespace xtp
99} // namespace votca
100
101#endif // VOTCA_XTP_AOPOTENTIAL_H
Container to hold Basisfunctions for all atoms.
Definition aobasis.h:42
void FillPotential(const AOBasis &aobasis, const ECPAOBasis &ecp)
Definition aoecp.cc:56
void FillBlock(Eigen::Block< Eigen::MatrixXd > &, const AOShell &, const AOShell &) const final
Definition aopotential.h:58
void FillBlock(Eigen::Block< Eigen::MatrixXd > &matrix, const AOShell &shell_row, const AOShell &shell_col) const override
void setSite(const StaticSite *site)
Definition aopotential.h:78
void FillPotential(const AOBasis &aobasis, const QMMolecule &atoms)
const StaticSite * site_
Definition aopotential.h:80
void FillPotential(const AOBasis &aobasis, const std::vector< Eigen::Vector3d > &kpoints)
void FillBlock(Eigen::Block< Eigen::MatrixXcd > &matrix, const AOShell &shell_row, const AOShell &shell_col) const override
void setkVector(const Eigen::Vector3d &k)
Definition aopotential.h:94
Eigen::Vector3d k_
Definition aopotential.h:95
const Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > & Matrix() const
Definition aopotential.h:39
virtual void FillBlock(Eigen::Block< Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > > &matrix, const AOShell &shell_row, const AOShell &shell_col) const =0
Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > Fill(const AOBasis &aobasis) const
Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > aopotential_
Definition aopotential.h:49
Container to hold ECPs for all atoms.
Definition ecpaobasis.h:43
Class to represent Atom/Site in electrostatic.
Definition staticsite.h:37
base class for all analysis tools
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26