19#ifndef VOTCA_XTP_DIPOLEDIPOLEINTERACTION_H
20#define VOTCA_XTP_DIPOLEDIPOLEINTERACTION_H
28class DipoleDipoleInteraction;
35struct traits<
votca::xtp::DipoleDipoleInteraction>
36 :
public Eigen::internal::traits<Eigen::MatrixXd> {};
44 :
public Eigen::EigenBase<DipoleDipoleInteraction> {
57 const std::vector<PolarSegment>& segs)
61 size_ += 3 * seg.size();
80 operator bool()
const {
99 template <
typename Vtype>
100 Eigen::Product<DipoleDipoleInteraction, Vtype, Eigen::AliasFreeProduct>
103 Eigen::AliasFreeProduct>(*
this, x.derived());
113 if (seg1id == seg2id) {
114 return sites_[seg1id]->getPInv()(xyz1, xyz2);
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");
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++) {
130 result.segment<3>(3 * i) += site1.
getPInv() * v.segment<3>(3 * i);
131 for (
Index j = i + 1; j < segment_size; j++) {
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);
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>> {
162 typedef typename Product<votca::xtp::DipoleDipoleInteraction, Vtype>::Scalar
165 template <
typename Dest>
168 const Vtype& v,
const Scalar& alpha) {
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>();
const DipoleDipoleInteraction & xpr_
InnerIterator(const DipoleDipoleInteraction &xpr, const Index &id)
InnerIterator & operator++()
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
votca::Index StorageIndex
std::vector< const PolarSite * > sites_
const eeInteractor & interactor_
Eigen::VectorXd multiply(const Eigen::VectorXd &v) const
Class to represent Atom/Site in electrostatic+polarization.
const Eigen::Matrix3d & getPInv() const
Mediates interaction between polar and static sites.
Eigen::Matrix3d FillTholeInteraction(const PolarSite &site1, const PolarSite &site2) const
base class for all analysis tools
static void scaleAndAddTo(Dest &dst, const votca::xtp::DipoleDipoleInteraction &op, const Vtype &v, const Scalar &alpha)
Product< votca::xtp::DipoleDipoleInteraction, Vtype >::Scalar Scalar