votca 2024.2-dev
Loading...
Searching...
No Matches
espfit.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// VOTCA includes
22
23// Local VOTCA includes
24#include "votca/xtp/aomatrix.h"
26#include "votca/xtp/espfit.h"
27#include "votca/xtp/grid.h"
28#include "votca/xtp/orbitals.h"
29#include "votca/xtp/vxc_grid.h"
30
31namespace votca {
32namespace xtp {
33using std::flush;
34
36 const QMState& state, std::string gridsize) {
37
38 const Eigen::MatrixXd dmat = orbitals.DensityMatrixFull(state);
39 // setting up grid
40 Grid grid;
41 grid.setupCHELPGGrid(orbitals.QMAtoms());
43 << " Done setting up CHELPG grid with "
44 << grid.size() << " points " << flush;
45
46 // Calculating nuclear potential at gridpoints
47 AOBasis basis = orbitals.getDftBasis();
48 AOOverlap overlap;
49 overlap.Fill(basis);
50 double N_comp = dmat.cwiseProduct(overlap.Matrix()).sum();
51
52 Vxc_Grid numintgrid;
53 numintgrid.GridSetup(gridsize, orbitals.QMAtoms(), basis);
54 XTP_LOG(Log::info, log_) << TimeStamp() << " Setup " << gridsize
55 << " Numerical Grid with "
56 << numintgrid.getGridSize() << " gridpoints."
57 << flush;
58 DensityIntegration<Vxc_Grid> numway(numintgrid);
59 double N = numway.IntegrateDensity(dmat);
61 << TimeStamp()
62 << " Calculated Densities at Numerical Grid, Number of electrons is " << N
63 << flush;
64
65 if (std::abs(N - N_comp) > 0.001) {
66 XTP_LOG(Log::error, log_) << "=======================" << flush;
68 << "WARNING: Calculated Densities at Numerical Grid, Number of "
69 "electrons "
70 << N << " is far away from the the real value " << N_comp
71 << ", you should increase the accuracy of the integration grid."
72 << flush;
73 N = N_comp;
75 << "WARNING: Electronnumber set to " << N << flush;
76 XTP_LOG(Log::error, log_) << "=======================" << flush;
77 }
78
80 << TimeStamp() << " Calculating ESP at CHELPG grid points" << flush;
81#pragma omp parallel for
82 for (Index i = 0; i < grid.size(); i++) {
83 grid.getGridValues()(i) =
84 numway.IntegratePotential(grid.getGridPositions()[i]);
85 }
86
87 XTP_LOG(Log::info, log_) << TimeStamp() << " Electron contribution calculated"
88 << flush;
89 double netcharge = -N;
90 if (!state.isTransition()) {
91 EvalNuclearPotential(orbitals.QMAtoms(), grid);
92 Index Znuc = 0;
93 for (const QMAtom& atom : orbitals.QMAtoms()) {
94 Znuc += atom.getNuccharge();
95 }
96 netcharge += double(Znuc);
97 }
98 netcharge = std::round(netcharge);
100 << TimeStamp() << " Netcharge constrained to " << netcharge << flush;
101 return FitPartialCharges(orbitals, grid, netcharge);
102 ;
103}
104
106
107 const std::vector<Eigen::Vector3d>& gridpoints = grid.getGridPositions();
108 Eigen::VectorXd& gridvalues = grid.getGridValues();
110 << " Calculating ESP of nuclei at CHELPG grid points"
111 << flush;
112
113 for (Index i = 0; i < Index(gridpoints.size()); i++) {
114 for (Index j = 0; j < atoms.size(); j++) {
115 const Eigen::Vector3d& posatom = atoms[j].getPos();
116 Index Znuc = atoms[j].getNuccharge();
117 double dist_j = (gridpoints[i] - posatom).norm();
118 gridvalues(i) += double(Znuc) / dist_j;
119 }
120 }
121 return;
122}
123
125 const Grid& grid, double netcharge) {
126 const QMMolecule& atomlist = orbitals.QMAtoms();
127 const Index NoOfConstraints =
128 1 + Index(regionconstraint_.size()) + Index(pairconstraint_.size());
129 const Index matrixSize = atomlist.size() + NoOfConstraints;
131 << " Setting up Matrices for fitting of size "
132 << matrixSize << " x " << matrixSize << flush;
133
134 const std::vector<Eigen::Vector3d>& gridpoints = grid.getGridPositions();
135 const Eigen::VectorXd& potential = grid.getGridValues();
136 XTP_LOG(Log::info, log_) << TimeStamp() << " Using " << atomlist.size()
137 << " Fittingcenters and " << gridpoints.size()
138 << " Gridpoints." << flush;
139
140 Eigen::MatrixXd Amat = Eigen::MatrixXd::Zero(matrixSize, matrixSize);
141 Eigen::VectorXd Bvec = Eigen::VectorXd::Zero(matrixSize);
142// setting up Amat_
143#pragma omp parallel for
144 for (Index i = 0; i < atomlist.size(); i++) {
145 for (Index j = i; j < atomlist.size(); j++) {
146 for (const auto& gridpoint : gridpoints) {
147 double dist_i = (atomlist[i].getPos() - gridpoint).norm();
148 double dist_j = (atomlist[j].getPos() - gridpoint).norm();
149
150 Amat(i, j) += 1.0 / dist_i / dist_j;
151 }
152 Amat(j, i) = Amat(i, j);
153 }
154 }
155
156 // setting up Bvec
157#pragma omp parallel for
158 for (Index i = 0; i < atomlist.size(); i++) {
159 for (Index j = 0; j < Index(gridpoints.size()); j++) {
160 double dist_i = (atomlist[i].getPos() - gridpoints[j]).norm();
161 Bvec(i) += potential(j) / dist_i;
162 }
163 }
164 // Total charge constraint
165 for (Index i = 0; i < atomlist.size() + 1; i++) {
166 Amat(i, atomlist.size()) = 1.0;
167 Amat(atomlist.size(), i) = 1.0;
168 }
169 Amat(atomlist.size(), atomlist.size()) = 0.0;
170 Bvec(atomlist.size()) = netcharge; // netcharge!!!!
171
172 // Pairconstraint
173 for (Index i = 0; i < Index(pairconstraint_.size()); i++) {
174 const std::pair<Index, Index>& pair = pairconstraint_[i];
175 Amat(pair.first, atomlist.size() + 1 + i) = 1.0;
176 Amat(atomlist.size() + 1 + i, pair.first) = 1.0;
177 Amat(pair.second, atomlist.size() + 1 + i) = -1.0;
178 Amat(atomlist.size() + 1 + i, pair.second) = -1.0;
179 }
180
181 // Regionconstraint
182 for (Index i = 0; i < Index(regionconstraint_.size()); i++) {
184 for (Index index : reg) {
185 Amat(index, atomlist.size() + i + 1 + pairconstraint_.size()) = 1.0;
186 Amat(atomlist.size() + i + 1 + pairconstraint_.size(), index) = 1.0;
187 }
188 Bvec(atomlist.size() + i + 1 + pairconstraint_.size()) = reg.value();
189 }
190
191 XTP_LOG(Log::info, log_) << TimeStamp() << " Solving linear Equation "
192 << flush;
193 Eigen::VectorXd charges;
194 if (do_svd_) {
195 Eigen::JacobiSVD<Eigen::MatrixXd> svd;
196 svd.setThreshold(conditionnumber_);
197 svd.compute(Amat, Eigen::ComputeThinU | Eigen::ComputeThinV);
198 charges = svd.solve(Bvec);
199 XTP_LOG(Log::info, log_) << TimeStamp() << " SVD Done. " << flush;
200 if ((Bvec.size() - svd.nonzeroSingularValues()) != 0) {
202 << TimeStamp() << Bvec.size() - svd.nonzeroSingularValues()
203 << " Sites could not be fitted and are set to zero." << flush;
204 }
205 } else {
206 Eigen::ColPivHouseholderQR<Eigen::MatrixXd> QR(Amat);
207 if (QR.info() != Eigen::ComputationInfo::Success) {
208 throw std::runtime_error(
209 "Espfit: Solving the constrained equation failed. Maybe try SVD.");
210 }
211 charges = QR.solve(Bvec);
213 << TimeStamp() << " Solved linear least square fit ." << flush;
214 }
215 // remove constraints from charges
216 charges.conservativeResize(atomlist.size());
217 StaticSegment seg =
218 StaticSegment(orbitals.QMAtoms().getType(), orbitals.QMAtoms().getId());
219
221 << " Sum of fitted charges: " << charges.sum() << flush;
222 for (Index i = 0; i < atomlist.size(); i++) {
223 seg.push_back(StaticSite(atomlist[i], charges(i)));
224 }
225 // get RMSE
226 double rmse = 0.0;
227 double totalPotSq = 0.0;
228 for (Index k = 0; k < Index(gridpoints.size()); k++) {
229 double temp = 0.0;
230 for (const StaticSite& atom : seg) {
231 double dist = (gridpoints[k] - atom.getPos()).norm();
232 temp += atom.getCharge() / dist;
233 }
234 rmse += (potential(k) - temp) * (potential(k) - temp);
235 totalPotSq += potential(k) * potential(k);
236 }
238 << " RMSE of fit: " << std::sqrt(rmse / double(gridpoints.size()))
239 << flush;
241 << " RRMSE of fit: " << std::sqrt(rmse / totalPotSq) << flush;
242
243 return seg;
244}
245
246} // namespace xtp
247} // namespace votca
Container to hold Basisfunctions for all atoms.
Definition aobasis.h:42
void Fill(const AOBasis &aobasis) final
const Eigen::MatrixXd & Matrix() const
Definition aomatrix.h:52
const std::string & getType() const
void push_back(const T &atom)
const Eigen::Vector3d & getPos() const
double IntegratePotential(const Eigen::Vector3d &rvector) const
double IntegrateDensity(const Eigen::MatrixXd &density_matrix)
std::vector< QMFragment< double > > regionconstraint_
Definition espfit.h:74
Logger & log_
Definition espfit.h:64
StaticSegment FitPartialCharges(const Orbitals &orbitals, const Grid &grid, double netcharge)
Definition espfit.cc:124
std::vector< std::pair< Index, Index > > pairconstraint_
Definition espfit.h:68
StaticSegment Fit2Density(const Orbitals &orbitals, const QMState &state, std::string gridsize)
Definition espfit.cc:35
double conditionnumber_
Definition espfit.h:66
void EvalNuclearPotential(const QMMolecule &atoms, Grid &grid)
Definition espfit.cc:105
void setupCHELPGGrid(const QMMolecule &Atomlist)
Definition grid.h:54
Eigen::VectorXd & getGridValues()
Definition grid.h:48
Index size()
Definition grid.h:50
const std::vector< Eigen::Vector3d > & getGridPositions() const
Definition grid.h:44
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 AOBasis & getDftBasis() const
Definition orbitals.h:208
container for QM atoms
Definition qmatom.h:37
const T & value() const
Definition qmfragment.h:59
Identifier for QMstates. Strings like S1 are converted into enum +zero indexed int.
Definition qmstate.h:132
bool isTransition() const
Definition qmstate.h:153
Class to represent Atom/Site in electrostatic.
Definition staticsite.h:37
Timestamp returns the current time as a string Example: cout << TimeStamp()
Definition logger.h:224
Index getGridSize() const
Definition vxc_grid.h:41
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
ClassicalSegment< StaticSite > StaticSegment
base class for all analysis tools
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26