votca 2024.2-dev
Loading...
Searching...
No Matches
gyration.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// Third party includes
21#include <boost/format.hpp>
22
23// VOTCA includes
25
26// Local VOTCA includes
27#include "votca/xtp/gyration.h"
28#include "votca/xtp/vxc_grid.h"
29
30using namespace votca::tools;
31
32namespace votca {
33namespace xtp {
34
36 std::string key = Identify();
37 state_ = options.get(key + ".state").as<QMState>();
39 key + ".difference_to_groundstate", false);
40 gridsize_ = options.ifExistsReturnElseReturnDefault<std::string>(
41 key + ".gridsize", "medium");
42}
43
45 XTP_LOG(Log::error, log_) << "===== Running on " << OPENMP::getMaxThreads()
46 << " threads ===== " << std::flush;
47
48 const QMMolecule& Atomlist = orbitals.QMAtoms();
49 BasisSet bs;
50 bs.Load(orbitals.getDFTbasisName());
51 AOBasis basis;
52 basis.Fill(bs, Atomlist);
53 AnalyzeGeometry(Atomlist);
54
55 // setup numerical integration grid
56 Vxc_Grid grid;
57 grid.GridSetup(gridsize_, Atomlist, basis);
59
60 if (!dostateonly_) {
61 Eigen::MatrixXd DMAT_tot = orbitals.DensityMatrixFull(state_);
62 Gyrationtensor gyro = numway.IntegrateGyrationTensor(DMAT_tot);
63 Eigen::SelfAdjointEigenSolver<Eigen::Matrix3d> es;
64 es.computeDirect(gyro.gyration);
66 << TimeStamp() << " Converting to Eigenframe " << std::flush;
67 XTP_LOG(Log::error, log_) << TimeStamp() << " Reporting " << std::flush;
69
70 } else {
71 // hole density first
72 std::array<Eigen::MatrixXd, 2> DMAT =
74 Gyrationtensor gyro_hole = numway.IntegrateGyrationTensor(DMAT[0]);
75 Eigen::SelfAdjointEigenSolver<Eigen::Matrix3d> es_h;
76 es_h.computeDirect(gyro_hole.gyration);
78 << TimeStamp() << " Converting to Eigenframe " << std::flush;
79 XTP_LOG(Log::error, log_) << TimeStamp() << " Reporting " << std::flush;
80 ReportAnalysis("hole", gyro_hole, es_h);
81
82 // electron density
83 Gyrationtensor gyro_electron = numway.IntegrateGyrationTensor(DMAT[1]);
84 Eigen::SelfAdjointEigenSolver<Eigen::Matrix3d> es_e;
85 es_e.computeDirect(gyro_electron.gyration);
87 << TimeStamp() << " Converting to Eigenframe " << std::flush;
88 XTP_LOG(Log::error, log_) << TimeStamp() << " Reporting " << std::flush;
89 ReportAnalysis("electron", gyro_electron, es_e);
90 }
91 return;
92}
93
95
96 tools::Elements elements;
97 double mass = 0.0;
98 Eigen::Vector3d centroid = Eigen::Vector3d::Zero();
99 Eigen::Matrix3d gyration = Eigen::Matrix3d::Zero();
100 for (const QMAtom& atom : atoms) {
101 double m = elements.getMass(atom.getElement());
102 const Eigen::Vector3d& pos = atom.getPos();
103 mass += m;
104 centroid += m * pos;
105 gyration += m * pos * pos.transpose();
106 }
107 centroid /= mass;
108 gyration /= mass;
109 gyration -= centroid * centroid.transpose();
110 Gyrationtensor gyro;
111 gyro.mass = mass;
112 gyro.centroid = centroid;
113 gyro.gyration = gyration;
114 Eigen::SelfAdjointEigenSolver<Eigen::Matrix3d> es;
115 es.computeDirect(gyro.gyration);
116 ReportAnalysis("geometry", gyro, es);
117}
118
120 std::string label, const Gyrationtensor& gyro,
121 const Eigen::SelfAdjointEigenSolver<Eigen::Matrix3d>& es) {
122
124 << "---------------- " << label << " ----------------" << std::flush;
126 << (boost::format(" Norm = %1$9.4f ") % (gyro.mass))
127 << std::flush;
128
130 << (boost::format(" Centroid x = %1$9.4f Ang") %
131 (gyro.centroid.x() * tools::conv::bohr2ang))
132 << std::flush;
134 << (boost::format(" Centroid y = %1$9.4f Ang") %
135 (gyro.centroid.y() * tools::conv::bohr2ang))
136 << std::flush;
138 << (boost::format(" Centroid y = %1$9.4f Ang") %
139 (gyro.centroid.z() * tools::conv::bohr2ang))
140 << std::flush;
141
144 << (boost::format(" Gyration Tensor xx = %1$9.4f Ang^2") %
145 (gyro.gyration(0, 0) * RA2))
146 << std::flush;
148 << (boost::format(" Gyration Tensor xy = %1$9.4f Ang^2") %
149 (gyro.gyration(0, 1) * RA2))
150 << std::flush;
152 << (boost::format(" Gyration Tensor xz = %1$9.4f Ang^2") %
153 (gyro.gyration(0, 2) * RA2))
154 << std::flush;
156 << (boost::format(" Gyration Tensor yy = %1$9.4f Ang^2") %
157 (gyro.gyration(1, 1) * RA2))
158 << std::flush;
160 << (boost::format(" Gyration Tensor yz = %1$9.4f Ang^2") %
161 (gyro.gyration(1, 2) * RA2))
162 << std::flush;
164 << (boost::format(" Gyration Tensor zz = %1$9.4f Ang^2") %
165 (gyro.gyration(2, 2) * RA2))
166 << std::flush;
167
169 << (boost::format(" Gyration Tensor D1 = %1$9.4f Ang^2") %
170 (es.eigenvalues()[0] * RA2))
171 << std::flush;
173 << (boost::format(" Gyration Tensor D2 = %1$9.4f Ang^2") %
174 (es.eigenvalues()[1] * RA2))
175 << std::flush;
177 << (boost::format(" Gyration Tensor D3 = %1$9.4f Ang^2") %
178 (es.eigenvalues()[2] * RA2))
179 << std::flush;
180
182 << (boost::format(" Radius of Gyration = %1$9.4f Ang") %
183 (std::sqrt(es.eigenvalues().sum()) * tools::conv::bohr2ang))
184 << std::flush;
185
187 << (boost::format(" Tensor EF Axis 1 1 = %1$9.4f ") %
188 es.eigenvectors().col(0).x())
189 << std::flush;
191 << (boost::format(" Tensor EF Axis 1 2 = %1$9.4f ") %
192 es.eigenvectors().col(0).y())
193 << std::flush;
195 << (boost::format(" Tensor EF Axis 1 3 = %1$9.4f ") %
196 es.eigenvectors().col(0).z())
197 << std::flush;
199 << (boost::format(" Tensor EF Axis 2 1 = %1$9.4f ") %
200 es.eigenvectors().col(1).x())
201 << std::flush;
203 << (boost::format(" Tensor EF Axis 2 2 = %1$9.4f ") %
204 es.eigenvectors().col(1).y())
205 << std::flush;
207 << (boost::format(" Tensor EF Axis 2 3 = %1$9.4f ") %
208 es.eigenvectors().col(1).z())
209 << std::flush;
211 << (boost::format(" Tensor EF Axis 3 1 = %1$9.4f ") %
212 es.eigenvectors().col(2).x())
213 << std::flush;
215 << (boost::format(" Tensor EF Axis 3 2 = %1$9.4f ") %
216 es.eigenvectors().col(2).y())
217 << std::flush;
219 << (boost::format(" Tensor EF Axis 3 3 = %1$9.4f ") %
220 es.eigenvectors().col(2).z())
221 << std::flush;
222}
223
224} // namespace xtp
225} // namespace votca
information about an element
Definition elements.h:42
double getMass(std::string name)
Returns the mass of each atom in a.u.
Definition elements.cc:62
class to manage program options with xml serialization functionality
Definition property.h:55
Property & get(const std::string &key)
get existing property
Definition property.cc:79
T as() const
return value as type
Definition property.h:283
T ifExistsReturnElseReturnDefault(const std::string &key, T defaultvalue) const
Definition property.h:321
Container to hold Basisfunctions for all atoms.
Definition aobasis.h:42
void Fill(const BasisSet &bs, const QMMolecule &atoms)
Definition aobasis.cc:85
void Load(const std::string &name)
Definition basisset.cc:149
void AnalyzeGeometry(const QMMolecule &atoms)
Definition gyration.cc:94
void ReportAnalysis(std::string label, const Gyrationtensor &gyro, const Eigen::SelfAdjointEigenSolver< Eigen::Matrix3d > &es)
Definition gyration.cc:119
void Initialize(tools::Property &options)
Definition gyration.cc:35
void AnalyzeDensity(const Orbitals &orbitals)
Definition gyration.cc:44
Gyrationtensor IntegrateGyrationTensor(const Eigen::MatrixXd &density_matrix)
container for molecular orbitals
Definition orbitals.h:46
Eigen::MatrixXd DensityMatrixFull(const QMState &state) const
Definition orbitals.cc:140
const QMMolecule & QMAtoms() const
Definition orbitals.h:172
const std::string & getDFTbasisName() const
Definition orbitals.h:202
std::array< Eigen::MatrixXd, 2 > DensityMatrixExcitedState(const QMState &state) const
Definition orbitals.cc:302
container for QM atoms
Definition qmatom.h:37
Identifier for QMstates. Strings like S1 are converted into enum +zero indexed int.
Definition qmstate.h:132
std::string ToLongString() const
Definition qmstate.cc:130
Timestamp returns the current time as a string Example: cout << TimeStamp()
Definition logger.h:224
void GridSetup(const std::string &type, const QMMolecule &atoms, const AOBasis &basis)
Definition vxc_grid.cc:229
#define XTP_LOG(level, log)
Definition logger.h:40
const double bohr2ang
Definition constants.h:49
Index getMaxThreads()
Definition eigen.h:128
base class for all analysis tools
Definition basebead.h:33