votca 2024-dev
Loading...
Searching...
No Matches
eeinteractor.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_EEINTERACTOR_H
22#define VOTCA_XTP_EEINTERACTOR_H
23
24// Local VOTCA includes
25#include "classicalsegment.h"
26#include "eigen.h"
27
28namespace votca {
29namespace xtp {
30
31enum Estatic : bool {
32 noE_V = true,
33 V = false,
34};
35
40 public:
41 explicit eeInteractor() = default;
42 explicit eeInteractor(double expdamping) : expdamping_(expdamping){};
43
44 Eigen::Matrix3d FillTholeInteraction(const PolarSite& site1,
45 const PolarSite& site2) const;
46
47 Eigen::VectorXd Cholesky_IntraSegment(const PolarSegment& seg) const;
48
49 template <class T, enum Estatic>
50 double ApplyStaticField(const T& segment1, PolarSegment& segment2) const;
51
52 template <enum Estatic>
53 double ApplyInducedField(const PolarSegment& segment1,
54 PolarSegment& segment2) const;
55
56 template <class S1, class S2>
57 double CalcStaticEnergy(const S1& segment1, const S2& segment2) const;
58
59 class E_terms {
60 public:
61 E_terms& operator+=(const E_terms& right) {
62 this->data_ += right.data_;
63 return *this;
64 }
65
66 E_terms operator+(E_terms right) const {
67 right.data_ += this->data_;
68 return right;
69 }
70
71 double sum() { return data_.sum(); }
72
73 double& E_indu_indu() { return data_.x(); }
74 double& E_indu_stat() { return data_.y(); }
75 double& E_internal() { return data_.z(); }
76
77 const Eigen::Vector3d& data() const { return data_; }
78
79 private:
80 Eigen::Vector3d data_ = Eigen::Vector3d::Zero();
81 };
82
83 template <class S1, class S2>
84 E_terms CalcPolarEnergy(const S1& segment1, const S2& segment2) const;
85
86 template <class S>
87 double CalcStaticEnergy_IntraSegment(const S& seg) const;
88
89 double CalcPolarEnergy_IntraSegment(const PolarSegment& seg) const;
90
91 double CalcStaticEnergy_site(const StaticSite& site1,
92 const StaticSite& site2) const;
93
94 private:
95 template <int N>
96 Eigen::Matrix<double, N, 1> VSiteA(const StaticSite& site1,
97 const StaticSite& site2) const;
98 template <enum Estatic>
99 double ApplyInducedField_site(const PolarSite& site1, PolarSite& site2) const;
100 template <enum Estatic>
101 double ApplyStaticField_site(const StaticSite& site1, PolarSite& site2) const;
102
103 double CalcPolar_stat_Energy_site(const PolarSite& site1,
104 const StaticSite& site2) const;
105
107 const StaticSite& site2) const;
108
110 const PolarSite& site2) const;
111
112 double expdamping_ = 0.39; // dimensionless
113};
114
115} // namespace xtp
116} // namespace votca
117
118#endif // VOTCA_XTP_EEINTERACTOR_H
Class to represent Atom/Site in electrostatic+polarization.
Definition polarsite.h:36
Class to represent Atom/Site in electrostatic.
Definition staticsite.h:37
const Eigen::Vector3d & data() const
E_terms & operator+=(const E_terms &right)
E_terms operator+(E_terms right) const
Mediates interaction between polar and static sites.
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
eeInteractor(double expdamping)
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