votca 2024.2-dev
Loading...
Searching...
No Matches
polarregion.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// Standard includes
21#include <iomanip>
22#include <numeric>
23
24// Local VOTCA includes
28#include "votca/xtp/qmregion.h"
30
31namespace votca {
32namespace xtp {
33
35 max_iter_ = prop.get("max_iter").as<Index>();
36 deltaD_ = prop.get("tolerance_dipole").as<double>();
37 deltaE_ = prop.get("tolerance_energy").as<double>();
38 exp_damp_ = prop.get("exp_damp").as<double>();
39}
40
42
43 if (!E_hist_.filled()) {
44 return false;
45 }
46 double Echange = E_hist_.getDiff().Etotal();
47 std::string info = "not converged";
48 bool converged = false;
49 if (std::abs(Echange) < deltaE_) {
50 info = "converged";
51 converged = true;
52 }
54 << TimeStamp() << " Region:" << this->identify() << " " << this->getId()
55 << " is " << info << " deltaE=" << Echange << std::flush;
56 return converged;
57}
58
60
61 eeInteractor eeinteractor;
62 double e = 0.0;
63#pragma omp parallel for reduction(+ : e)
64 for (Index i = 0; i < size(); ++i) {
65 for (Index j = 0; j < size(); ++j) {
66 if (i == j) {
67 continue;
68 }
69 PolarSegment& seg1 = segments_[i];
70 const PolarSegment& seg2 = segments_[j];
71 e += eeinteractor.ApplyStaticField<PolarSegment, Estatic::noE_V>(seg2,
72 seg1);
73 }
74 }
75
76 return 0.5 * e;
77}
78
80#pragma omp declare reduction(CustomPlus : eeInteractor::E_terms : omp_out += \
81 omp_in)
82
83 eeInteractor eeinteractor(exp_damp_);
84
86
87#pragma omp parallel for reduction(CustomPlus : terms)
88 for (Index i = 0; i < size(); ++i) {
89 for (Index j = 0; j < i; ++j) {
90 terms += eeinteractor.CalcPolarEnergy(segments_[i], segments_[j]);
91 }
92 }
93
94#pragma omp parallel for reduction(CustomPlus : terms)
95 for (Index i = 0; i < size(); ++i) {
96 terms.E_indu_indu() +=
98 }
99
100#pragma omp parallel for reduction(CustomPlus : terms)
101 for (Index i = 0; i < size(); ++i) {
102 for (const PolarSite& site : segments_[i]) {
103 terms.E_internal() += site.InternalEnergy();
104 }
105 }
106 return terms;
107}
108
110 double e = 0.0;
111#pragma omp parallel for reduction(+ : e)
112 for (Index i = 0; i < size(); ++i) {
113 for (const PolarSite& site : segments_[i]) {
114 e += site.deltaQ_V_ext();
115 }
116 }
117 return e;
118}
119
121 for (PolarSegment& seg : segments_) {
122 for (PolarSite& site : seg) {
123 site.Reset();
124 }
125 }
126}
127
129 const Energy_terms& e = this->E_hist_.back();
130 prop.add("E_static", std::to_string(e.Estatic() * tools::conv::hrt2ev));
131 prop.add("E_polar", std::to_string(e.Epolar() * tools::conv::hrt2ev));
132 prop.add("E_total", std::to_string(e.Etotal() * tools::conv::hrt2ev));
133}
134
136 Index dof_polarization = 0;
137 for (const PolarSegment& seg : segments_) {
138 dof_polarization += seg.size() * 3;
139 }
140 return dof_polarization;
141}
142
144 Index dof_polarization = CalcPolDoF();
145 Eigen::VectorXd initial_induced_dipoles =
146 Eigen::VectorXd::Zero(dof_polarization);
147 eeInteractor interactor(exp_damp_);
148 Index index = 0;
149 for (const PolarSegment& seg : segments_) {
150 initial_induced_dipoles.segment(index, 3 * seg.size()) =
151 interactor.Cholesky_IntraSegment(seg);
152 index += 3 * seg.size();
153 }
154 return initial_induced_dipoles;
155}
156
158 Index dof_polarization = CalcPolDoF();
159 Eigen::VectorXd last_induced_dipoles =
160 Eigen::VectorXd::Zero(dof_polarization);
161 Index index = 0;
162 for (const PolarSegment& seg : segments_) {
163 for (const PolarSite& site : seg) {
164 last_induced_dipoles.segment<3>(index) = site.Induced_Dipole();
165 index += 3;
166 }
167 }
168 return last_induced_dipoles;
169}
170
171void PolarRegion::WriteInducedDipolesToSegments(const Eigen::VectorXd& x) {
172 Index index = 0;
173 for (PolarSegment& seg : segments_) {
174 for (PolarSite& site : seg) {
175 site.setInduced_Dipole(x.segment<3>(index));
176 index += 3;
177 }
178 }
179}
180
182 const Eigen::VectorXd& initial_guess) {
183 Eigen::VectorXd b = Eigen::VectorXd::Zero(initial_guess.size());
184 Index index = 0;
185 for (PolarSegment& seg : segments_) {
186 for (const PolarSite& site : seg) {
187 auto V = site.V() + site.V_noE();
188 b.segment<3>(index) = -V;
189 index += 3;
190 }
191 }
192 eeInteractor interactor(exp_damp_);
193 DipoleDipoleInteraction A(interactor, segments_);
194 Eigen::ConjugateGradient<DipoleDipoleInteraction, Eigen::Lower | Eigen::Upper,
195 Eigen::DiagonalPreconditioner<double>>
196 cg;
197 cg.setMaxIterations(max_iter_);
198 cg.setTolerance(deltaD_);
199 cg.compute(A);
200 Eigen::VectorXd x = cg.solveWithGuess(b, initial_guess);
201
203 << TimeStamp() << " CG: #iterations: " << cg.iterations()
204 << ", estimated error: " << cg.error() << std::flush;
205
206 if (cg.info() == Eigen::ComputationInfo::NoConvergence) {
207 info_ = false;
208 errormsg_ = "PCG iterations did not converge";
209 }
210
211 if (cg.info() == Eigen::ComputationInfo::NumericalIssue) {
212 info_ = false;
213 errormsg_ = "PCG had a numerical issue";
214 }
215 return x;
216}
217
218void PolarRegion::Evaluate(std::vector<std::unique_ptr<Region>>& regions) {
219
220 std::vector<double> energies = ApplyInfluenceOfOtherRegions(regions);
221 Energy_terms e_contrib;
222 e_contrib.E_static_ext() =
223 std::accumulate(energies.begin(), energies.end(), 0.0);
225 << " Calculated static-static and polar-static "
226 "interaction with other regions"
227 << std::flush;
228 e_contrib.E_static_static() = StaticInteraction();
230 << " Calculated static interaction in region "
231 << std::flush;
232 Index dof_polarization = CalcPolDoF();
234 << TimeStamp() << " Starting Solving for classical polarization with "
235 << dof_polarization << " degrees of freedom." << std::flush;
236
237 Eigen::VectorXd initial_induced_dipoles;
238 if (!E_hist_.filled() || segments_.size() == 1) {
239 initial_induced_dipoles = CalcInducedDipoleInsideSegments();
240 } else {
241 initial_induced_dipoles = ReadInducedDipolesFromLastIteration();
242 }
243
244 Eigen::VectorXd x; // if only one segment
245 // it is solved exactly through the initial guess
246 if (segments_.size() != 1) {
247 x = CalcInducedDipolesViaPCG(initial_induced_dipoles);
248 if (!info_) {
249 return;
250 }
251 } else {
252 x = initial_induced_dipoles;
253 }
254
256
259 << " Calculated polar interaction in region"
260 << std::flush;
261 e_contrib.E_polar_ext() = PolarEnergy_extern();
263 << " Calculated polar interaction with other regions"
264 << std::flush;
265
266 XTP_LOG(Log::info, log_) << std::setprecision(10)
267 << " Internal static energy [hrt]= "
268 << e_contrib.E_static_static() << std::flush;
269 XTP_LOG(Log::info, log_) << std::setprecision(10)
270 << " External static energy [hrt]= "
271 << e_contrib.E_static_ext() << std::flush;
273 << std::setprecision(10)
274 << " Total static energy [hrt]= " << e_contrib.Estatic() << std::flush;
275
276 XTP_LOG(Log::info, log_) << std::setprecision(10)
277 << " internal dQ-dQ energy [hrt]= "
278 << e_contrib.E_indu_indu() << std::flush;
279
280 XTP_LOG(Log::info, log_) << std::setprecision(10)
281 << " internal Q-dQ energy [hrt]= "
282 << e_contrib.E_indu_stat() << std::flush;
283
284 XTP_LOG(Log::info, log_) << std::setprecision(10)
285 << " Internal energy [hrt]= "
286 << e_contrib.E_internal() << std::flush;
287
288 XTP_LOG(Log::info, log_) << std::setprecision(10)
289 << " External polar energy [hrt]= "
290 << e_contrib.E_polar_ext() << std::flush;
291
293 << std::setprecision(10)
294 << " Total polar energy [hrt]= " << e_contrib.Epolar() << std::flush;
295
297 << std::setprecision(10) << " Total energy [hrt]= " << e_contrib.Etotal()
298 << std::flush;
299 E_hist_.push_back(e_contrib);
300 return;
301}
302
304 // QMregions always have lower ids than other regions
306 return 0.0;
307}
309 bool noE_V = true;
310 if (this->getId() < region.getId()) {
311 noE_V = false;
312 }
313
314 double e = 0;
315#pragma omp parallel for reduction(+ : e)
316 for (Index i = 0; i < Index(segments_.size()); i++) {
317 PolarSegment& pseg1 = segments_[i];
318 double e_thread = 0.0;
319 eeInteractor ee;
320 for (const PolarSegment& pseg2 : region) {
321 if (noE_V) {
323 ee.ApplyInducedField<Estatic::noE_V>(pseg2, pseg1);
324 } else {
325 e_thread += ee.ApplyStaticField<PolarSegment, Estatic::V>(pseg2, pseg1);
326 e_thread += ee.ApplyInducedField<Estatic::V>(pseg2, pseg1);
327 }
328 e += e_thread;
329 }
330 }
331 return e;
332}
333
335 // Static regions always have higher ids than other regions
336
337 double e = 0.0;
338#pragma omp parallel for reduction(+ : e)
339 for (Index i = 0; i < Index(segments_.size()); i++) {
340 double e_thread = 0.0;
341 PolarSegment& pseg = segments_[i];
342 eeInteractor ee;
343 for (const StaticSegment& sseg : region) {
344 e_thread += ee.ApplyStaticField<StaticSegment, Estatic::V>(sseg, pseg);
345 }
346 e += e_thread;
347 }
348
349 return e;
350}
351
355
359
360} // namespace xtp
361} // namespace votca
class to manage program options with xml serialization functionality
Definition property.h:55
Property & add(const std::string &key, const std::string &value)
add a new property to structure
Definition property.cc:108
Property & get(const std::string &key)
get existing property
Definition property.cc:79
T as() const
return value as type
Definition property.h:283
void addInternalPolarContrib(const eeInteractor::E_terms &induction_terms)
Index size() const override
Definition mmregion.h:45
std::vector< PolarSegment > segments_
Definition mmregion.h:86
void WriteToCpt(CheckpointWriter &w) const override
Definition mmregion.cc:45
void ReadFromCpt(CheckpointReader &r) override
Definition mmregion.cc:58
Eigen::VectorXd CalcInducedDipoleInsideSegments() const
void AppendResult(tools::Property &prop) const override
void WriteInducedDipolesToSegments(const Eigen::VectorXd &x)
void ReadFromCpt(CheckpointReader &r) override
double InteractwithQMRegion(const QMRegion &region) override
void Initialize(const tools::Property &prop) override
bool Converged() const override
eeInteractor::E_terms PolarEnergy() const
Eigen::VectorXd CalcInducedDipolesViaPCG(const Eigen::VectorXd &initial_guess)
Eigen::VectorXd ReadInducedDipolesFromLastIteration() const
std::string identify() const override
Definition polarregion.h:48
hist< Energy_terms > E_hist_
Definition polarregion.h:86
void WriteToCpt(CheckpointWriter &w) const override
double InteractwithPolarRegion(const PolarRegion &region) override
void Evaluate(std::vector< std::unique_ptr< Region > > &regions) override
void Reset() override
double InteractwithStaticRegion(const StaticRegion &region) override
double PolarEnergy_extern() const
Class to represent Atom/Site in electrostatic+polarization.
Definition polarsite.h:36
void ApplyQMFieldToPolarSegments(std::vector< PolarSegment > &segments) const
Definition qmregion.cc:377
std::vector< double > ApplyInfluenceOfOtherRegions(std::vector< std::unique_ptr< Region > > &regions)
Definition region.cc:32
Index getId() const
Definition region.h:79
std::string errormsg_
Definition region.h:91
Logger & log_
Definition region.h:100
Timestamp returns the current time as a string Example: cout << TimeStamp()
Definition logger.h:224
Mediates interaction between polar and static sites.
double ApplyStaticField(const T &segment1, PolarSegment &segment2) const
double ApplyInducedField(const PolarSegment &segment1, PolarSegment &segment2) const
E_terms CalcPolarEnergy(const S1 &segment1, const S2 &segment2) const
Eigen::VectorXd Cholesky_IntraSegment(const PolarSegment &seg) const
double CalcPolarEnergy_IntraSegment(const PolarSegment &seg) const
#define XTP_LOG(level, log)
Definition logger.h:40
const double hrt2ev
Definition constants.h:53
ClassicalSegment< PolarSite > PolarSegment
ClassicalSegment< StaticSite > StaticSegment
base class for all analysis tools
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26