votca 2024-dev
Loading...
Searching...
No Matches
polarsite.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_POLARSITE_H
22#define VOTCA_XTP_POLARSITE_H
23
24// Local VOTCA includes
25#include "eigen.h"
26#include "staticsite.h"
27
28namespace votca {
29namespace xtp {
30
36class PolarSite final : public StaticSite {
37
38 public:
39 // delete these two functions because we do not want to be able to read
40 // StaticSite::data but PolarSite::data
41 void WriteData(StaticSite::data& d) const = delete;
42 void ReadData(StaticSite::data& d) = delete;
43
44 PolarSite(Index id, std::string element, Eigen::Vector3d pos);
45 PolarSite(Index id, std::string element)
46 : PolarSite(id, element, Eigen::Vector3d::Zero()){};
47
48 ~PolarSite() final = default;
49
50 void setpolarization(const Eigen::Matrix3d& pol) final;
51
52 Eigen::Matrix3d getpolarization() const { return pinv_.inverse(); }
53
54 const Eigen::Matrix3d& getPInv() const { return pinv_; }
55
56 // MULTIPOLES DEFINITION
57 Eigen::Vector3d getDipole() const final;
58
59 double getSqrtInvEigenDamp() const { return eigendamp_invsqrt_; }
60
61 void Rotate(const Eigen::Matrix3d& R, const Eigen::Vector3d& ref_pos) final {
62 StaticSite::Rotate(R, ref_pos);
63 pinv_ = R.transpose() * pinv_ * R;
64 }
65
66 const Eigen::Vector3d& V() const { return V_; }
67
68 Eigen::Vector3d& V() { return V_; }
69
70 const Eigen::Vector3d& V_noE() const { return V_noE_; }
71
72 Eigen::Vector3d& V_noE() { return V_noE_; }
73
74 void Reset() {
75 V_.setZero();
76 V_noE_.setZero();
77 }
78
79 double deltaQ_V_ext() const { return induced_dipole_.dot(V_); }
80
81 double InternalEnergy() const {
82 return 0.5 * induced_dipole_.transpose() * pinv_ * induced_dipole_;
83 }
84
85 const Eigen::Vector3d& Induced_Dipole() const { return induced_dipole_; }
86 void setInduced_Dipole(const Eigen::Vector3d& induced_dipole) {
87 induced_dipole_ = induced_dipole;
88 }
89
90 struct data {
92 char* element;
93 double posX;
94 double posY;
95 double posZ;
96
98
99 double Q00;
100 double Q11c;
101 double Q11s;
102 double Q10;
103 double Q20;
104 double Q21c;
105 double Q21s;
106 double Q22c;
107 double Q22s;
108
109 double Vx;
110 double Vy;
111 double Vz;
112
113 double Vx_noE;
114 double Vy_noE;
115 double Vz_noE;
116
117 double pxx;
118 double pxy;
119 double pxz;
120 double pyy;
121 double pyz;
122 double pzz;
123
124 double d_x_ind;
125 double d_y_ind;
126 double d_z_ind;
127 };
128 // do not move up has to be below data definition
129 PolarSite(const data& d);
130
131 double DipoleChange() const;
132
133 static void SetupCptTable(CptTable& table);
134 void WriteData(data& d) const;
135 void ReadData(const data& d);
136
137 std::string identify() const final { return "polarsite"; }
138
139 friend std::ostream& operator<<(std::ostream& out, const PolarSite& site) {
140 out << site.getId() << " " << site.getElement() << " " << site.getRank();
141 out << " " << site.getPos().transpose() << " "
142 << site.Induced_Dipole().transpose() << "\n";
143 return out;
144 }
145
146 private:
147 std::string writepolarization() const final;
148
149 // PolarSite has two external fields,
150 // the first is used for interaction with regions, which are further out, i.e.
151 // the interaction energy with it is included in the polar region energy
152 Eigen::Vector3d V_ = Eigen::Vector3d::Zero();
153 // the second is used for interaction with regions, which are further inside,
154 // i.e. the interaction energy with it is included in the other region's
155 // energy
156 Eigen::Vector3d V_noE_ = Eigen::Vector3d::Zero();
157
158 Eigen::Vector3d induced_dipole_ = Eigen::Vector3d::Zero();
159 Eigen::Matrix3d pinv_ = Eigen::Matrix3d::Zero();
160 double eigendamp_invsqrt_ = 0.0;
161};
162
163} // namespace xtp
164} // namespace votca
165
166#endif // VOTCA_XTP_POLARSITE_H
Class to represent Atom/Site in electrostatic+polarization.
Definition polarsite.h:36
Eigen::Vector3d & V_noE()
Definition polarsite.h:72
double DipoleChange() const
double deltaQ_V_ext() const
Definition polarsite.h:79
double getSqrtInvEigenDamp() const
Definition polarsite.h:59
double InternalEnergy() const
Definition polarsite.h:81
void setInduced_Dipole(const Eigen::Vector3d &induced_dipole)
Definition polarsite.h:86
const Eigen::Matrix3d & getPInv() const
Definition polarsite.h:54
const Eigen::Vector3d & Induced_Dipole() const
Definition polarsite.h:85
static void SetupCptTable(CptTable &table)
Definition polarsite.cc:75
friend std::ostream & operator<<(std::ostream &out, const PolarSite &site)
Definition polarsite.h:139
PolarSite(Index id, std::string element)
Definition polarsite.h:45
const Eigen::Vector3d & V_noE() const
Definition polarsite.h:70
void WriteData(StaticSite::data &d) const =delete
void Rotate(const Eigen::Matrix3d &R, const Eigen::Vector3d &ref_pos) final
Definition polarsite.h:61
Eigen::Vector3d V_noE_
Definition polarsite.h:156
Eigen::Vector3d induced_dipole_
Definition polarsite.h:158
void setpolarization(const Eigen::Matrix3d &pol) final
Definition polarsite.cc:58
Eigen::Matrix3d getpolarization() const
Definition polarsite.h:52
std::string identify() const final
Definition polarsite.h:137
Eigen::Matrix3d pinv_
Definition polarsite.h:159
const Eigen::Vector3d & V() const
Definition polarsite.h:66
std::string writepolarization() const final
Definition polarsite.cc:66
Eigen::Vector3d getDipole() const final
Definition polarsite.cc:54
Eigen::Vector3d & V()
Definition polarsite.h:68
void ReadData(StaticSite::data &d)=delete
~PolarSite() final=default
Eigen::Vector3d V_
Definition polarsite.h:152
Class to represent Atom/Site in electrostatic.
Definition staticsite.h:37
const std::string & getElement() const
Definition staticsite.h:79
Index getId() const
Definition staticsite.h:77
const Eigen::Vector3d & getPos() const
Definition staticsite.h:80
Index getRank() const
Definition staticsite.h:78
virtual void Rotate(const Eigen::Matrix3d &R, const Eigen::Vector3d &refPos)
Definition staticsite.cc:66
base class for all analysis tools
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26