votca 2024.1-dev
Loading...
Searching...
No Matches
dipoledipoleinteraction.h
Go to the documentation of this file.
1/*
2 * Copyright 2009-2020 The VOTCA Development Team (http://www.votca.org)
3 *
4 * Licensed under the Apache License, Version 2.0 (the "License");
5 * you may not use this file except in compliance with the License.
6 * You may obtain a copy of the License at
7 *
8 * http://www.apache.org/licenses/LICENSE-2.0
9 *
10 * Unless required by applicable law or agreed to in writing, software
11 * distributed under the License is distributed on an "AS IS" BASIS,
12 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
13 * See the License for the specific language governing permissions and
14 * limitations under the License.
15 *
16 */
17
18#pragma once
19#ifndef VOTCA_XTP_DIPOLEDIPOLEINTERACTION_H
20#define VOTCA_XTP_DIPOLEDIPOLEINTERACTION_H
21
22// Local VOTCA includes
23#include "eeinteractor.h"
24#include "eigen.h"
25
26namespace votca {
27namespace xtp {
28class DipoleDipoleInteraction;
29}
30} // namespace votca
31namespace Eigen {
32namespace internal {
33// MatrixReplacement looks-like a SparseMatrix, so let's inherits its traits:
34template <>
35struct traits<votca::xtp::DipoleDipoleInteraction>
36 : public Eigen::internal::traits<Eigen::MatrixXd> {};
37} // namespace internal
38} // namespace Eigen
39
40namespace votca {
41namespace xtp {
42
44 : public Eigen::EigenBase<DipoleDipoleInteraction> {
45 public:
46 // Required typedefs, constants, and method:
47 using Scalar = double;
48 using RealScalar = double;
50 enum {
51 ColsAtCompileTime = Eigen::Dynamic,
52 MaxColsAtCompileTime = Eigen::Dynamic,
53 IsRowMajor = false
54 };
55
57 const std::vector<PolarSegment>& segs)
58 : interactor_(interactor) {
59 size_ = 0;
60 for (const PolarSegment& seg : segs) {
61 size_ += 3 * seg.size();
62 }
63 sites_.reserve(size_ / 3);
64 for (const PolarSegment& seg : segs) {
65 for (const PolarSite& site : seg) {
66 sites_.push_back(&site);
67 }
68 }
69 }
70
72 public:
74 : xpr_(xpr), id_(id) {};
75
77 row_++;
78 return *this;
79 }
80 operator bool() const {
81 return row_ < xpr_.size_;
82 } // DO not use the size method, it returns linear dimension*linear
83 // dimension i.e size_^2
84 double value() const { return xpr_(row_, id_); }
85 Index row() const { return row_; }
86 Index col() const { return id_; }
87 Index index() const { return row(); }
88
89 private:
91 const Index id_;
93 };
94
95 Index rows() const { return this->size_; }
96 Index cols() const { return this->size_; }
97 Index outerSize() const { return this->size_; }
98
99 template <typename Vtype>
100 Eigen::Product<DipoleDipoleInteraction, Vtype, Eigen::AliasFreeProduct>
101 operator*(const Eigen::MatrixBase<Vtype>& x) const {
102 return Eigen::Product<DipoleDipoleInteraction, Vtype,
103 Eigen::AliasFreeProduct>(*this, x.derived());
104 }
105
106 // this is not a fast method
107 double operator()(const Index i, const Index j) const {
108 Index seg1id = Index(i / 3);
109 Index xyz1 = Index(i % 3);
110 Index seg2id = Index(j / 3);
111 Index xyz2 = Index(j % 3);
112
113 if (seg1id == seg2id) {
114 return sites_[seg1id]->getPInv()(xyz1, xyz2);
115 } else {
116 const PolarSite& site1 = *sites_[seg1id];
117 const PolarSite& site2 = *sites_[seg2id];
118 return interactor_.FillTholeInteraction(site1, site2)(xyz1, xyz2);
119 }
120 };
121
122 Eigen::VectorXd multiply(const Eigen::VectorXd& v) const {
123 assert(v.size() == size_ &&
124 "input vector has the wrong size for multiply with operator");
125 const Index segment_size = Index(sites_.size());
126 Eigen::VectorXd result = Eigen::VectorXd::Zero(size_);
127#pragma omp parallel for schedule(dynamic) reduction(+ : result)
128 for (Index i = 0; i < segment_size; i++) {
129 const PolarSite& site1 = *sites_[i];
130 result.segment<3>(3 * i) += site1.getPInv() * v.segment<3>(3 * i);
131 for (Index j = i + 1; j < segment_size; j++) {
132 const PolarSite& site2 = *sites_[j];
133 Eigen::Matrix3d block = interactor_.FillTholeInteraction(site1, site2);
134 result.segment<3>(3 * i) += block * v.segment<3>(3 * j);
135 result.segment<3>(3 * j) += block.transpose() * v.segment<3>(3 * i);
136 }
137 }
138
139 return result;
140 }
141
142 private:
144 std::vector<const PolarSite*> sites_;
146};
147} // namespace xtp
148} // namespace votca
149
150namespace Eigen {
151
152namespace internal {
153
154// replacement of the mat*vect operation
155template <typename Vtype>
156struct generic_product_impl<votca::xtp::DipoleDipoleInteraction, Vtype,
157 DenseShape, DenseShape, GemvProduct>
158 : generic_product_impl_base<
159 votca::xtp::DipoleDipoleInteraction, Vtype,
160 generic_product_impl<votca::xtp::DipoleDipoleInteraction, Vtype>> {
161
162 typedef typename Product<votca::xtp::DipoleDipoleInteraction, Vtype>::Scalar
164
165 template <typename Dest>
166 static void scaleAndAddTo(Dest& dst,
168 const Vtype& v, const Scalar& alpha) {
169 // returns dst = alpha * op * v
170 // alpha must be 1 here
171 assert(alpha == Scalar(1) && "scaling is not implemented");
172 EIGEN_ONLY_USED_FOR_DEBUG(alpha);
173 Eigen::VectorXd temp = op.multiply(v);
174 dst = temp.cast<Scalar>(); // tumbleweed fix do not delete
175 }
176};
177
178} // namespace internal
179} // namespace Eigen
180
181#endif // VOTCA_XTP_DIPOLEDIPOLEINTERACTION_H
InnerIterator(const DipoleDipoleInteraction &xpr, const Index &id)
DipoleDipoleInteraction(const eeInteractor &interactor, const std::vector< PolarSegment > &segs)
double operator()(const Index i, const Index j) const
Eigen::Product< DipoleDipoleInteraction, Vtype, Eigen::AliasFreeProduct > operator*(const Eigen::MatrixBase< Vtype > &x) const
std::vector< const PolarSite * > sites_
Eigen::VectorXd multiply(const Eigen::VectorXd &v) const
Class to represent Atom/Site in electrostatic+polarization.
Definition polarsite.h:36
const Eigen::Matrix3d & getPInv() const
Definition polarsite.h:54
Mediates interaction between polar and static sites.
Eigen::Matrix3d FillTholeInteraction(const PolarSite &site1, const PolarSite &site2) const
base class for all analysis tools
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26
static void scaleAndAddTo(Dest &dst, const votca::xtp::DipoleDipoleInteraction &op, const Vtype &v, const Scalar &alpha)