votca 2024.1-dev
Loading...
Searching...
No Matches
staticsite.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 <fstream>
22#include <string>
23
24// Third party includes
25#include <boost/format.hpp>
26
27// VOTCA includes
29
30// Local VOTCA includes
33
34namespace votca {
35namespace xtp {
36
38 // We are transforming here just quadrupoles
39 // const Eigen::VectorXd& MP = multipole_;
40 const Vector9d& MP = Q_;
41 Eigen::Matrix3d theta = Eigen::Matrix3d::Zero();
42 if (rank_ > 1) {
43 double sqr3 = std::sqrt(3);
44 theta(0, 0) = 0.5 * (-MP(4) + sqr3 * MP(7)); // theta_xx
45 theta(1, 1) = 0.5 * (-MP(4) - sqr3 * MP(7)); // theta_yy
46 theta(2, 2) = MP(4); // theta_zz
47 theta(0, 1) = theta(1, 0) = 0.5 * sqr3 * MP(8); // theta_xy = theta_yx
48 theta(0, 2) = theta(2, 0) = 0.5 * sqr3 * MP(5); // theta_xz = theta_zx
49 theta(1, 2) = theta(2, 1) = 0.5 * sqr3 * MP(6); // theta_yz = theta_zy
50 }
51 return theta;
52}
53
55 const Eigen::Matrix3d& quad_cart) {
56 Eigen::VectorXd quadrupole_polar = Eigen::VectorXd::Zero(5);
57 const double sqr3 = std::sqrt(3);
58 quadrupole_polar(0) = quad_cart(2, 2);
59 quadrupole_polar(1) = (2. / sqr3) * quad_cart(0, 2);
60 quadrupole_polar(2) = (2. / sqr3) * quad_cart(1, 2);
61 quadrupole_polar(3) = (1. / sqr3) * (quad_cart(0, 0) - quad_cart(1, 1));
62 quadrupole_polar(4) = (2. / sqr3) * quad_cart(0, 1);
63 return quadrupole_polar;
64}
65
66void StaticSite::Rotate(const Eigen::Matrix3d& R,
67 const Eigen::Vector3d& refPos) {
68 Eigen::Vector3d dir = pos_ - refPos;
69 dir = R * dir;
70 pos_ = refPos + dir; // Rotated Position
71 if (rank_ > 0) {
72 const Eigen::Vector3d temp = R * Q_.segment<3>(1);
73 Q_.segment<3>(1) = temp;
74 }
75 if (rank_ > 1) {
76 Eigen::Matrix3d cartesianquad = CalculateCartesianMultipole();
77 Eigen::Matrix3d rotated = R * cartesianquad * R.transpose();
78 Q_.segment<5>(4) = CalculateSphericalMultipole(rotated);
79 }
80 return;
81}
82
83void StaticSite::Translate(const Eigen::VectorXd& shift) {
84 pos_ += shift;
85 return;
86}
87
88std::string StaticSite::writepolarization() const {
90 double default_pol = 1; // default is alway 1A^3
91 try {
92 default_pol =
94 } catch (const std::runtime_error&) {
95 ;
96 }
97 return (boost::format(" P %1$+1.7f\n") % default_pol).str();
98}
99
100std::string StaticSite::WriteMpsLine(std::string unit) const {
101 double conv_pos = 1.;
102 if (unit == "angstrom") {
103 conv_pos = tools::conv::bohr2ang;
104 } else if (unit == "bohr") {
105 conv_pos = 1.;
106 } else {
107 throw std::runtime_error(
108 " StaticSite::WriteMpsLine: Unit conversion not known");
109 }
110 std::string output = "";
111 output += (boost::format(" %1$2s %2$+1.7f %3$+1.7f %4$+1.7f Rank %5$d\n") %
112 element_ % (pos_(0) * conv_pos) % (pos_(1) * conv_pos) %
113 (pos_(2) * conv_pos) % rank_)
114 .str();
115 output += (boost::format(" %1$+1.7f\n") % getCharge()).str();
116 if (rank_ > 0) {
117 // Dipole z x y
118 output += (boost::format(" %1$+1.7f %2$+1.7f %3$+1.7f\n") % Q_(3) %
119 Q_(1) % Q_(2))
120 .str();
121 if (rank_ > 1) {
122 // Quadrupole 20 21c 21s 22c 22s
123 output +=
124 (boost::format(" %1$+1.7f %2$+1.7f %3$+1.7f %4$+1.7f %5$+1.7f\n") %
125 Q_(4) % Q_(5) % Q_(6) % Q_(7) % Q_(8))
126 .str();
127 }
128 }
129 // Polarizability
130 output += writepolarization();
131 return output;
132}
133
135 table.addCol<Index>("index", HOFFSET(data, id));
136 table.addCol<std::string>("type", HOFFSET(data, element));
137
138 table.addCol<double>("posX", HOFFSET(data, posX));
139 table.addCol<double>("posY", HOFFSET(data, posY));
140 table.addCol<double>("posZ", HOFFSET(data, posZ));
141
142 table.addCol<Index>("rank", HOFFSET(data, rank));
143
144 table.addCol<double>("Q00", HOFFSET(data, Q00));
145 table.addCol<double>("Q11c", HOFFSET(data, Q11c));
146 table.addCol<double>("Q11s", HOFFSET(data, Q11s));
147 table.addCol<double>("Q10", HOFFSET(data, Q10));
148 table.addCol<double>("Q20", HOFFSET(data, Q20));
149 table.addCol<double>("Q21c", HOFFSET(data, Q21c));
150 table.addCol<double>("Q21s", HOFFSET(data, Q21s));
151 table.addCol<double>("Q22c", HOFFSET(data, Q22c));
152 table.addCol<double>("Q22s", HOFFSET(data, Q22s));
153}
154
156 d.id = id_;
157 d.element = const_cast<char*>(element_.c_str());
158 d.posX = pos_[0];
159 d.posY = pos_[1];
160 d.posZ = pos_[2];
161
162 d.rank = rank_;
163
164 d.Q00 = Q_[0];
165 d.Q11c = Q_[1];
166 d.Q11s = Q_[2];
167 d.Q10 = Q_[3];
168 d.Q20 = Q_[4];
169 d.Q21c = Q_[5];
170 d.Q21s = Q_[6];
171 d.Q22c = Q_[7];
172 d.Q22s = Q_[8];
173}
174
176 id_ = d.id;
177 element_ = std::string(d.element);
178 free(d.element);
179 pos_[0] = d.posX;
180 pos_[1] = d.posY;
181 pos_[2] = d.posZ;
182
183 rank_ = d.rank;
184
185 Q_[0] = d.Q00;
186 Q_[1] = d.Q11c;
187 Q_[2] = d.Q11s;
188 Q_[3] = d.Q10;
189 Q_[4] = d.Q20;
190 Q_[5] = d.Q21c;
191 Q_[6] = d.Q21s;
192 Q_[7] = d.Q22c;
193 Q_[8] = d.Q22s;
194}
195
196} // namespace xtp
197} // namespace votca
information about an element
Definition elements.h:42
double getPolarizability(std::string name)
Returns the atomic polarisability of atom.
Definition elements.cc:97
void addCol(const std::string &name, const size_t &offset)
static Eigen::VectorXd CalculateSphericalMultipole(const Eigen::Matrix3d &quad_cart)
Definition staticsite.cc:54
Eigen::Matrix3d CalculateCartesianMultipole() const
Definition staticsite.cc:37
static void SetupCptTable(CptTable &table)
std::string WriteMpsLine(std::string unit="bohr") const
double getCharge() const
Definition staticsite.h:99
void WriteData(data &d) const
virtual std::string writepolarization() const
Definition staticsite.cc:88
void ReadData(const data &d)
void Translate(const Eigen::VectorXd &shift)
Definition staticsite.cc:83
virtual void Rotate(const Eigen::Matrix3d &R, const Eigen::Vector3d &refPos)
Definition staticsite.cc:66
Eigen::Vector3d pos_
Definition staticsite.h:133
const double bohr2ang
Definition constants.h:49
const double nm2ang
Definition constants.h:50
base class for all analysis tools
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26
Eigen::Matrix< double, 9, 1 > Vector9d
Definition eigen.h:34