votca 2024.2-dev
Loading...
Searching...
No Matches
eeinteractor.cc
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// Local VOTCA includes
22
23namespace votca {
24namespace xtp {
25
26template <int N>
27Eigen::Matrix<double, N, 1> eeInteractor::VSiteA(
28 const StaticSite& siteA, const StaticSite& siteB) const {
29
30 const Eigen::Vector3d& posA = siteA.getPos();
31 const Eigen::Vector3d& posB = siteB.getPos();
32 Index rankB = siteB.getRank();
33 Eigen::Vector3d a =
34 posB - posA; // Vector of the distance between polar sites
35 const double R = a.norm();
36 const double fac1 = 1.0 / R;
37 a *= fac1; // unit vector pointing from A to B
38 Eigen::Matrix<double, N, 1> V;
39 V(0) = fac1 * siteB.getCharge();
40
41 if (N > 1 || rankB > 0) {
42 const double fac2 = std::pow(fac1, 2);
43 // Dipole-Charge Interaction
44 if constexpr (N > 1) {
45 V.template segment<3>(1) =
46 fac2 * a * siteB.getCharge(); // T_1alpha,00 (alpha=x,y,z)
47 }
48 // Charge-Dipole Interaction
49 if (rankB > 0) {
50 V(0) -=
51 fac2 * a.dot(siteB.Q().segment<3>(1)); // T_00,1alpha (alpha=x,y,z)
52 }
53
54 const double fac3 = std::pow(fac1, 3);
55 if (N > 1 && rankB > 0) {
56 // Dipole-Dipole Interaction
57 Eigen::Vector3d result =
58 -3 * a * (a.transpose() * siteB.Q().segment<3>(1));
59 result += siteB.Q().segment<3>(1);
60 result *= fac3;
61 V.template segment<3>(1) += result;
62 // T_1alpha,1beta // (alpha,beta=x,y,z)
63 }
64
65 if (N > 4 || rankB > 1) {
66 const double sqr3 = std::sqrt(3);
67 const AxA r(a);
68 {
69 // Quadrupole-Charge interaction
70 Eigen::Matrix<double, 1, 5> Qq;
71 Qq(0) = 1.5 * r.zz() - 0.5; // T20,00
72 Qq(1) = r.xz(); // T21c,00
73 Qq(2) = r.yz(); // T21s,000
74 Qq(3) = 0.5 * (r.xx() - r.yy()); // T22c,00
75 Qq(4) = r.xy(); // T22s,00
76 Qq.tail<4>() *= sqr3;
77 if constexpr (N > 4) {
78 V.template segment<5>(4) = fac3 * Qq.transpose() * siteB.getCharge();
79 }
80 if (rankB > 1) {
81 V(0) += fac3 * Qq * siteB.Q().segment<5>(4);
82 }
83 }
84 if ((N > 4 || rankB > 0) || (N > 1 || rankB > 1)) {
85
86 Eigen::Matrix<double, 3, 5> dQ;
87 const double fac4 = std::pow(fac1, 4);
88 // Quadrupole-Dipole Interaction
89 const Eigen::Vector3d afac = a * fac4;
90
91 dQ.col(0) = (1.5 - 7.5 * r.zz()) * afac;
92 dQ.col(0).z() += 3 * afac.z(); // T20-1beta (beta=x,y,z)
93
94 dQ.col(1) = -5 * r.xz() * afac;
95 dQ.col(1).z() += afac.x();
96 dQ.col(1).x() += afac.z(); // T21c-1beta (beta=x,y,z)
97
98 dQ.col(2) = -5 * r.yz() * afac;
99 dQ.col(2).z() += afac.y();
100 dQ.col(2).y() += afac.z();
101
102 dQ.col(3) = -2.5 * (r.xx() - r.yy()) * afac;
103 dQ.col(3).x() += afac.x();
104 dQ.col(3).y() -= afac.y(); // T22c-1beta (beta=x,y,z)
105
106 dQ.col(4) = -5 * r.xy() * afac;
107 dQ.col(4).y() += afac.x();
108 dQ.col(4).x() += afac.y(); // T22s-1beta (beta=x,y,z)
109
110 dQ.rightCols<4>() *= sqr3;
111
112 if constexpr (N > 4) {
113 V.template segment<5>(4) += dQ.transpose() * siteB.Q().segment<3>(1);
114 }
115 if (rankB > 1) {
116 V.template segment<3>(1) -= dQ * siteB.Q().segment<5>(4);
117 }
118 }
119
120 if (N > 4 && rankB > 1) {
121
122 // Quadrupole-Quadrupole Interaction
123 Eigen::Matrix<double, 5, 5> QQ;
124 QQ(0, 0) = 0.75 * ((35 * r.zz() - 30) * r.zz() + 3); // T20,20
125 double temp = 0.5 * sqr3 * (35 * r.zz() - 15);
126 QQ(1, 0) = temp * r.xz(); // T20,21c
127 QQ(2, 0) = temp * r.yz(); // T20,21s
128
129 temp = 5 * (7 * r.zz() - 1);
130 QQ(3, 0) = sqr3 * 0.25 * temp * (r.xx() - r.yy()); // T20,22c
131 QQ(4, 0) = sqr3 * 0.5 * temp * r.xy(); // T20,22s
132 QQ(1, 1) =
133 35 * r.zz() * r.xx() - 5 * (r.xx() + r.zz()) + 1; // T21c,21c
134 QQ(2, 1) = r.xy() * temp; // T21c,21s
135 QQ(3, 1) = 0.5 * r.xz() * (35 * (r.xx() - r.yy()) - 10); // T21c,22c
136 QQ(4, 1) = r.yz() * 5 * (7 * r.xx() - 1); // T21c,22s
137 QQ(2, 2) =
138 5 * (7 * r.yy() * r.zz() - (r.yy() + r.zz())) + 1; // T21s,21s
139 QQ(3, 2) = 0.5 * r.yz() * (35 * (r.xx() - r.yy()) + 10); // T21s,22c
140 QQ(4, 2) = r.xz() * 5 * (7 * r.yy() - 1); // T21s,22s
141 QQ(3, 3) = 8.75 * std::pow(r.xx() - r.yy(), 2) - 5 * (r.xx() + r.yy()) +
142 1; // T22c,22c
143 QQ(4, 3) = 17.5 * r.xy() * (r.xx() - r.yy()); // T22c,22s
144 QQ(4, 4) =
145 5 * (7 * r.xx() * r.yy() - (r.xx() + r.yy())) + 1; // T22s,22s
146
147 const double fac5 = std::pow(fac1, 5);
148 V.template segment<5>(4) +=
149 fac5 * QQ.selfadjointView<Eigen::Lower>() * siteB.Q().segment<5>(4);
150 }
151 }
152 }
153 return V;
154}
155
157 const PolarSite& site1, const PolarSite& site2) const {
158
159 const Eigen::Vector3d& posB = site2.getPos();
160 const Eigen::Vector3d& posA = site1.getPos();
161 Eigen::Vector3d a =
162 posB - posA; // Vector of the distance between polar sites
163 const double R = a.norm(); // Norm of distance vector
164 const double fac1 = 1 / R;
165 a *= fac1; // unit vector pointing from A to B
166
167 double lambda3 = std::pow(fac1, 3);
168 double lambda5 = lambda3;
169 const double au3 = expdamping_ * std::pow(R, 3) *
170 site1.getSqrtInvEigenDamp() *
171 site2.getSqrtInvEigenDamp(); // au3 is dimensionless
172 if (au3 < 40) {
173 const double exp_ua = std::exp(-au3);
174 lambda3 *= (1 - exp_ua);
175 lambda5 *= (1 - (1 + au3) * exp_ua);
176 }
177 Eigen::Matrix3d result = -3 * lambda5 * a * a.transpose();
178 result.diagonal().array() += lambda3;
179 return result; // T_1alpha,1beta (alpha,beta=x,y,z)
180}
181
182template <enum Estatic CE>
184 PolarSite& site2) const {
185 Eigen::Vector3d V = Eigen::Vector3d::Zero();
186 double e = 0.0;
187 if (site2.getRank() < 2) {
188 const Eigen::Vector4d V_full = VSiteA<4>(site2, site1);
189 V = V_full.segment<3>(1);
190 e = V_full.dot(site2.Q().head<4>());
191 } else {
192 const Eigen::Matrix<double, 9, 1> V_full = VSiteA<9>(site2, site1);
193 V = V_full.segment<3>(1);
194 e = V_full.dot(site2.Q());
195 }
196
197 if (CE == Estatic::noE_V) {
198 site2.V_noE() += V;
199 } else {
200 site2.V() += V;
201 }
202 return e;
203}
204
205template <enum Estatic CE>
207 PolarSite& site2) const {
208 const Eigen::Matrix3d tTab = FillTholeInteraction(site1, site2);
209 const Eigen::Vector3d V = tTab.transpose() * site1.Induced_Dipole();
210 double e = 0.0;
211 if (CE == Estatic::noE_V) {
212 site2.V_noE() += V;
213 } else {
214 site2.V() += V;
215 e = CalcPolar_stat_Energy_site(site1, site2);
216 }
217 return e;
218}
219
221 const StaticSite& site2) const {
222 double e = 0.0;
223 if (site2.getRank() < 2) {
224 const Eigen::Vector4d V_full = VSiteA<4>(site2, site1);
225 e = V_full.dot(site2.Q().head<4>());
226 } else {
227 const Eigen::Matrix<double, 9, 1> V_full = VSiteA<9>(site2, site1);
228 e = V_full.dot(site2.Q());
229 }
230 return e;
231}
232
234 const StaticSite& site2) const {
235 const Eigen::Vector4d V_full = VSiteA<4>(site1, site2);
236 return V_full.tail<3>().dot(site1.Induced_Dipole());
237}
238
240 const PolarSite& site1, const StaticSite& site2) const {
242 val.E_indu_stat() = CalcPolar_stat_Energy_site(site1, site2);
243 return val;
244}
245
247 const PolarSite& site1, const PolarSite& site2) const {
248 // contributions are stat-induced, induced-stat and induced induced
250 val.E_indu_indu() = site1.Induced_Dipole().transpose() *
251 FillTholeInteraction(site1, site2) *
252 site2.Induced_Dipole();
253 val.E_indu_stat() = CalcPolar_stat_Energy_site(site1, site2);
254 val.E_indu_stat() += CalcPolar_stat_Energy_site(site2, site1);
255 return val;
256}
257
258template <class T, enum Estatic CE>
259double eeInteractor::ApplyStaticField(const T& segment1,
260 PolarSegment& segment2) const {
261 double e = 0.0;
262 for (PolarSite& s2 : segment2) {
263 for (const auto& s1 : segment1) {
264 e += ApplyStaticField_site<CE>(s1, s2);
265 }
266 }
267 return e;
268}
269
270template double eeInteractor::ApplyStaticField<StaticSegment, Estatic::V>(
271 const StaticSegment& seg1, PolarSegment& seg2) const;
272template double eeInteractor::ApplyStaticField<StaticSegment, Estatic::noE_V>(
273 const StaticSegment& seg1, PolarSegment& seg2) const;
274template double eeInteractor::ApplyStaticField<PolarSegment, Estatic::V>(
275 const PolarSegment& seg1, PolarSegment& seg2) const;
276template double eeInteractor::ApplyStaticField<PolarSegment, Estatic::noE_V>(
277 const PolarSegment& seg1, PolarSegment& seg2) const;
278
279template <enum Estatic CE>
281 PolarSegment& segment2) const {
282 double e = 0;
283 for (PolarSite& s2 : segment2) {
284 for (const PolarSite& s1 : segment1) {
285 e += ApplyInducedField_site<CE>(s1, s2);
286 }
287 }
288 return e;
289}
290
291template double eeInteractor::ApplyInducedField<Estatic::V>(
292 const PolarSegment& seg1, PolarSegment& seg2) const;
293template double eeInteractor::ApplyInducedField<Estatic::noE_V>(
294 const PolarSegment& seg1, PolarSegment& seg2) const;
295
296template <class S1, class S2>
297double eeInteractor::CalcStaticEnergy(const S1& segment1,
298 const S2& segment2) const {
299 double e = 0;
300 for (const auto& s1 : segment2) {
301 for (const auto& s2 : segment1) {
302 e += CalcStaticEnergy_site(s2, s1);
303 }
304 }
305 return e;
306}
307
308template double eeInteractor::CalcStaticEnergy(const StaticSegment& seg1,
309 const PolarSegment& seg2) const;
310template double eeInteractor::CalcStaticEnergy(const StaticSegment& seg1,
311 const StaticSegment& seg2) const;
312template double eeInteractor::CalcStaticEnergy(const PolarSegment& seg1,
313 const PolarSegment& seg2) const;
314template double eeInteractor::CalcStaticEnergy(const PolarSegment& seg1,
315 const StaticSegment& seg2) const;
316template <class S>
318 double e = 0;
319 for (Index i = 0; i < seg.size(); i++) {
320 for (Index j = 0; j < i; j++) {
321 e += CalcStaticEnergy_site(seg[i], seg[j]);
322 }
323 }
324 return e;
325}
327 const PolarSegment& seg1) const;
329 const StaticSegment& seg2) const;
330
331template <class S1, class S2>
333 const S2& segment2) const {
335 for (const auto& s1 : segment2) {
336 for (const auto& s2 : segment1) {
337 e += CalcPolarEnergy_site(s2, s1);
338 }
339 }
340 return e;
341}
342
344 const PolarSegment& seg1, const PolarSegment& seg2) const;
346 const PolarSegment& seg1, const StaticSegment& seg2) const;
347
349 const PolarSegment& seg) const {
350 double e = 0;
351 for (Index i = 0; i < seg.size(); i++) {
352 for (Index j = 0; j < i; j++) {
353 e += seg[i].Induced_Dipole().transpose() *
354 FillTholeInteraction(seg[i], seg[j]) * seg[j].Induced_Dipole();
355 }
356 }
357 return e;
358}
359
361 const PolarSegment& seg) const {
362 Index size = 3 * seg.size();
363
364 Eigen::MatrixXd A = Eigen::MatrixXd::Zero(size, size);
365 for (Index i = 1; i < seg.size(); i++) {
366 for (Index j = 0; j < i; j++) {
367 A.block<3, 3>(3 * i, 3 * j) = FillTholeInteraction(seg[i], seg[j]);
368 }
369 }
370
371 for (Index i = 0; i < seg.size(); i++) {
372 A.block<3, 3>(3 * i, 3 * i) = seg[i].getPInv();
373 }
374 Eigen::VectorXd b = Eigen::VectorXd(size);
375 for (Index i = 0; i < seg.size(); i++) {
376 const Eigen::Vector3d V = seg[i].V() + seg[i].V_noE();
377 b.segment<3>(3 * i) = -V;
378 }
379
380 Eigen::LLT<Eigen::Ref<Eigen::MatrixXd>> lltOfA(A);
381 return lltOfA.solve(b);
382}
383
384} // namespace xtp
385} // namespace votca
const double & zz() const
Definition eigen.h:103
const double & yz() const
Definition eigen.h:102
const double & yy() const
Definition eigen.h:101
const double & xy() const
Definition eigen.h:99
const double & xz() const
Definition eigen.h:100
const double & xx() const
Definition eigen.h:98
Class to represent Atom/Site in electrostatic+polarization.
Definition polarsite.h:36
double getSqrtInvEigenDamp() const
Definition polarsite.h:59
const Eigen::Vector3d & Induced_Dipole() const
Definition polarsite.h:85
const Eigen::Vector3d & V_noE() const
Definition polarsite.h:70
const Eigen::Vector3d & V() const
Definition polarsite.h:66
Class to represent Atom/Site in electrostatic.
Definition staticsite.h:37
const Eigen::Vector3d & getPos() const
Definition staticsite.h:80
double getCharge() const
Definition staticsite.h:99
Index getRank() const
Definition staticsite.h:78
const Vector9d & Q() const
Definition staticsite.h:100
double ApplyStaticField_site(const StaticSite &site1, PolarSite &site2) const
double ApplyStaticField(const T &segment1, PolarSegment &segment2) const
double ApplyInducedField(const PolarSegment &segment1, PolarSegment &segment2) const
double CalcStaticEnergy_site(const StaticSite &site1, const StaticSite &site2) const
Eigen::Matrix< double, N, 1 > VSiteA(const StaticSite &site1, const StaticSite &site2) const
E_terms CalcPolarEnergy(const S1 &segment1, const S2 &segment2) const
Eigen::VectorXd Cholesky_IntraSegment(const PolarSegment &seg) const
double CalcPolar_stat_Energy_site(const PolarSite &site1, const StaticSite &site2) const
E_terms CalcPolarEnergy_site(const PolarSite &site1, const StaticSite &site2) const
Eigen::Matrix3d FillTholeInteraction(const PolarSite &site1, const PolarSite &site2) const
double CalcStaticEnergy_IntraSegment(const S &seg) const
double ApplyInducedField_site(const PolarSite &site1, PolarSite &site2) const
double CalcStaticEnergy(const S1 &segment1, const S2 &segment2) const
double CalcPolarEnergy_IntraSegment(const PolarSegment &seg) const
base class for all analysis tools
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26