votca 2024.1-dev
Loading...
Searching...
No Matches
gaussianwriter.cc
Go to the documentation of this file.
1/*
2 * Copyright 2009-2021 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
21#include "votca/xtp/basisset.h"
22#include <boost/algorithm/string.hpp>
23#include <boost/format.hpp>
24#include <sstream>
26
27namespace votca {
28namespace xtp {
29
30/*
31 * This function converts VOTCA's L enum to the gaussian equivalent.
32 * Gaussian uses a minus sign to indicate spherical shells, but -1 for sp.
33 */
35 switch (l) {
36 case L::S:
37 return 0;
38 case L::P:
39 return 1;
40 default:
41 return -1 * static_cast<Index>(l);
42 }
43}
44
46 const Orbitals& orbitals) const {
48 Eigen::MatrixXd moCoefficients = orbitals.MOs().eigenvectors();
49 reorder.reorderOrbitals(moCoefficients, orbitals.getDftBasis());
50
51 // put the reordered mos in a string
52 std::stringstream mos_string;
53
54 int temp_int = 1;
55 for (Index i = 0; i < moCoefficients.rows(); ++i) {
56 for (Index j = 0; j < moCoefficients.cols(); ++j) {
57 mos_string << boost::format("%16.8e") % moCoefficients(j, i);
58 if (temp_int % 5 == 0) {
59 mos_string << "\n";
60 }
61 temp_int++;
62 }
63 }
64 mos_string << ((temp_int - 1) % 5 == 0 ? "" : "\n");
65
66 return mos_string.str();
67}
68
70 const QMState& state,
71 bool diff2gs) const {
73
74 Eigen::MatrixXd density;
75 if (diff2gs) {
76 density = orbitals.DensityMatrixWithoutGS(state);
77 } else {
78 density = orbitals.DensityMatrixFull(state);
79 }
80
81 reorder.reorderOperator(density, orbitals.getDftBasis());
82
83 // put the reordered mos in a string
84 std::stringstream density_string;
85
86 int temp_int = 1;
87
88 for (Index i = 0; i < density.rows(); ++i) {
89 for (Index j = 0; j <= i; ++j) {
90 density_string << boost::format("%16.8e") % density(i, j);
91 if (temp_int % 5 == 0) {
92 density_string << "\n";
93 }
94 temp_int++;
95 }
96 }
97 density_string << ((temp_int - 1) % 5 == 0 ? "" : "\n");
98
99 return density_string.str();
100}
101
102void GaussianWriter::WriteFile(const std::string& basename,
103 const Orbitals& orbitals, const QMState state,
104 bool diff2gs) const {
105 if (!orbitals.hasDFTbasisName()) {
106 throw std::runtime_error(".orb file does not contain a basisset name");
107 }
108
109 AOBasis basis = orbitals.getDftBasis();
110
111 std::ofstream outFile(basename + ".fchk");
112
113 if (outFile.is_open()) {
115 << "Start writing to " << (basename + ".fchk") << std::flush;
116 int temp_int;
117 // job description
118 outFile << basename << ", fchk created by VOTCA-XTP\n";
119 outFile << "SP RHF " << orbitals.getDFTbasisName() << "\n";
120
121 // clang-format off
122 outFile << boost::format("%-43s%-2s%15d\n") % "Number of atoms" % "I" % orbitals.QMAtoms().size();
123 outFile << boost::format("%-43s%-2s%15d\n") % "Charge" % "I" % 0;
124 outFile << boost::format("%-43s%-2s%15d\n") % "Multiplicity" % "I" % 1;
125 outFile << boost::format("%-43s%-2s%15d\n") % "Number of electrons" % "I" % (2*orbitals.getNumberOfAlphaElectrons());
126 outFile << boost::format("%-43s%-2s%15d\n") % "Number of alpha electrons" % "I" % orbitals.getNumberOfAlphaElectrons();
127 outFile << boost::format("%-43s%-2s%15d\n") % "Number of beta electrons" % "I" % orbitals.getNumberOfAlphaElectrons();
128 outFile << boost::format("%-43s%-2s%15d\n") % "Number of basis functions" % "I" % basis.AOBasisSize();
129 outFile << boost::format("%-43s%-2s%15d\n") % "Number of independent functions" % "I" % basis.AOBasisSize();
130 // clang-format on
131 // ATOMIC NUMBERS
132 outFile << boost::format("%-43s%-2s N= %10d\n") % "Atomic numbers" % "I" %
133 orbitals.QMAtoms().size();
134 temp_int = 1;
135 for (const auto& atom : orbitals.QMAtoms()) {
136 outFile << boost::format("%12d") % atom.getElementNumber();
137 if (temp_int % 6 == 0) {
138 outFile << "\n";
139 }
140 temp_int++;
141 }
142 outFile << ((temp_int - 1) % 6 == 0 ? "" : "\n");
143 // NUCLEAR CHARGES
144 outFile << boost::format("%-43s%-2s N= %10d\n") % "Nuclear charges" % "R" %
145 orbitals.QMAtoms().size();
146 temp_int = 1;
147 for (const auto& atom : orbitals.QMAtoms()) {
148 outFile << boost::format("%16.8e") % (double)atom.getNuccharge();
149 if (temp_int % 5 == 0) {
150 outFile << "\n";
151 }
152 temp_int++;
153 }
154 outFile << ((temp_int - 1) % 5 == 0 ? "" : "\n");
155 // CURRENT CARTESIAN COORDINATES
156 outFile << boost::format("%-43s%-2s N= %10d\n") %
157 "Current cartesian coordinates" % "R" %
158 (3 * orbitals.QMAtoms().size());
159 temp_int = 1;
160 for (const auto& atom : orbitals.QMAtoms()) {
161 for (int i = 0; i < 3; ++i) {
162 outFile << boost::format("%16.8e") % atom.getPos()(i);
163 if (temp_int % 5 == 0) {
164 outFile << "\n";
165 }
166 temp_int++;
167 }
168 }
169 outFile << ((temp_int - 1) % 5 == 0 ? "" : "\n");
170 // NUMBER OF PRIMITIVE SHELLS
171 outFile << boost::format("%-43s%-2s%15d\n") % "Number of primitive shells" %
172 "I" % basis.getNumberOfPrimitives();
173 // NUMBER OF CONTRACTED SHELLS
174 outFile << boost::format("%-43s%-2s%15d\n") %
175 "Number of contracted shells" % "I" % basis.getNumofShells();
176 // PURE/CARTESIAN D
177 outFile << boost::format("%-43s%-2s%15d\n") % "Pure/Cartesian d shells " %
178 "I" % 0;
179 // PURE/CARTESIAN F
180 outFile << boost::format("%-43s%-2s%15d\n") % "Pure/Cartesian f shells " %
181 "I" % 0;
182 // HIGHEST ANGULAR MOMENTUM
183 outFile << boost::format("%-43s%-2s%15d\n") % "Highest angular momentum " %
184 "I" % basis.getMaxL();
185 // HIGHEST ANGULAR MOMENTUM
186 outFile << boost::format("%-43s%-2s%15d\n") %
187 "Largest degree of contraction " % "I" % basis.getMaxNprim();
188 // SHELL TYPES
189 outFile << boost::format("%-43s%-2s N= %10d\n") % "Shell types" % "I" %
190 (basis.getNumofShells());
191 temp_int = 1;
192 for (const auto& shell : basis) {
193 outFile << boost::format("%12d") % toGaussianL(shell.getL());
194 if (temp_int % 6 == 0) {
195 outFile << "\n";
196 }
197 temp_int++;
198 }
199 outFile << ((temp_int - 1) % 6 == 0 ? "" : "\n");
200 // NR OF PRIMITIVES PER SHELL
201 outFile << boost::format("%-43s%-2s N= %10d\n") %
202 "Number of primitives per shell" % "I" %
203 (basis.getNumofShells());
204 temp_int = 1;
205 for (const AOShell& shell : basis) {
206 outFile << boost::format("%12d") % shell.getSize();
207 if (temp_int % 6 == 0) {
208 outFile << "\n";
209 }
210 temp_int++;
211 }
212 outFile << ((temp_int - 1) % 6 == 0 ? "" : "\n");
213 // SHELL TO ATOM MAP
214 outFile << boost::format("%-43s%-2s N= %10d\n") % "Shell to atom map" %
215 "I" % (basis.getNumofShells());
216 temp_int = 1;
217 for (const AOShell& shell : basis) {
218 // Gaussian indices start at 1, hence the + 1
219 outFile << boost::format("%12d") % (shell.getAtomIndex() + 1);
220 if (temp_int % 6 == 0) {
221 outFile << "\n";
222 }
223 temp_int++;
224 }
225 outFile << ((temp_int - 1) % 6 == 0 ? "" : "\n");
226 // PRIMITIVE EXPONENTS
227 outFile << boost::format("%-43s%-2s N= %10d\n") % "Primitive exponents" %
228 "R" % basis.getNumberOfPrimitives();
229 temp_int = 1;
230 for (const AOShell& shell : basis) {
231 for (const AOGaussianPrimitive& prim : shell) {
232 outFile << boost::format("%16.8e") % prim.getDecay();
233 if (temp_int % 5 == 0) {
234 outFile << "\n";
235 }
236 temp_int++;
237 }
238 }
239 outFile << ((temp_int - 1) % 5 == 0 ? "" : "\n");
240 // CONTRACTION COEFFICIENTS
241 outFile << boost::format("%-43s%-2s N= %10d\n") %
242 "Contraction coefficients" % "R" %
243 basis.getNumberOfPrimitives();
244 temp_int = 1;
245 for (const AOShell& shell : basis) {
246 for (const AOGaussianPrimitive& prim : shell) {
247 outFile << boost::format("%16.8e") % prim.getContraction();
248 if (temp_int % 5 == 0) {
249 outFile << "\n";
250 }
251 temp_int++;
252 }
253 }
254 outFile << ((temp_int - 1) % 5 == 0 ? "" : "\n");
255 // SHELL COORDINATES
256 outFile << boost::format("%-43s%-2s N= %10d\n") %
257 "Coordinates of each shell" % "R" %
258 (3 * basis.getNumofShells());
259 temp_int = 1;
260 for (const AOShell& shell : basis) {
261 for (int i = 0; i < 3; ++i) {
262 outFile << boost::format("%16.8e") % shell.getPos()(i);
263 if (temp_int % 5 == 0) {
264 outFile << "\n";
265 }
266 temp_int++;
267 }
268 }
269 outFile << ((temp_int - 1) % 5 == 0 ? "" : "\n");
270 // TOTAL ENERGY
271 outFile << boost::format("%-43s%-2s%22.15e\n") % "Total Energy" % "R" %
272 orbitals.getDFTTotalEnergy();
273 // ALPHA ORBITAL ENERGIES
274 outFile << boost::format("%-43s%-2s N= %10d\n") %
275 "Alpha Orbital Energies" % "R" %
276 orbitals.MOs().eigenvalues().size();
277 temp_int = 1;
278 for (Index i = 0; i < orbitals.MOs().eigenvalues().size(); ++i) {
279 outFile << boost::format("%16.8e") % orbitals.MOs().eigenvalues()[i];
280 if (temp_int % 5 == 0) {
281 outFile << "\n";
282 }
283 temp_int++;
284 }
285 outFile << ((temp_int - 1) % 5 == 0 ? "" : "\n");
286 // ALPHA ORBITAL ENERGIES
287 outFile << boost::format("%-43s%-2s N= %10d\n") % "Alpha MO coefficients" %
288 "R" %
289 (orbitals.MOs().eigenvalues().size() *
290 orbitals.MOs().eigenvalues().size());
291 outFile << reorderedMOCoefficients(orbitals);
292 // DENSITY MATRIX
293 outFile << boost::format("%-43s%-2s N= %10d\n") % "Total SCF Density" %
294 "R" %
295 ((orbitals.MOs().eigenvalues().size() *
296 (orbitals.MOs().eigenvalues().size() - 1)) /
297 2 +
298 orbitals.MOs().eigenvalues().size());
299 outFile << densityMatrixToString(orbitals, state, diff2gs);
300
301 XTP_LOG(Log::error, log_) << "Done writing \n" << std::flush;
302 }
303}
304
305} // namespace xtp
306} // namespace votca
const Eigen::VectorXd & eigenvalues() const
Definition eigensystem.h:30
const Eigen::MatrixXd & eigenvectors() const
Definition eigensystem.h:33
Container to hold Basisfunctions for all atoms.
Definition aobasis.h:42
Index getMaxL() const
Definition aobasis.cc:29
Index AOBasisSize() const
Definition aobasis.h:46
Index getNumberOfPrimitives() const
Definition aobasis.h:58
Index getNumofShells() const
Definition aobasis.h:56
Index getMaxNprim() const
Definition aobasis.cc:37
std::array< Index, 49 > gaussianMultipliers
std::string densityMatrixToString(const Orbitals &orbitals, const QMState &state, bool diff2gs) const
std::array< Index, 49 > gaussianOrder
std::string reorderedMOCoefficients(const Orbitals &orbitals) const
void WriteFile(const std::string &basename, const Orbitals &orbitals, const QMState state=QMState(QMStateType::statetype::Gstate, 0, false), bool diff2gs=false) const
void reorderOrbitals(Eigen::MatrixXd &moCoefficients, const AOBasis &basis)
Definition orbreorder.cc:74
void reorderOperator(Eigen::MatrixXd &Matrixoperator, const AOBasis &basis)
container for molecular orbitals
Definition orbitals.h:46
double getDFTTotalEnergy() const
Definition orbitals.h:193
Index getNumberOfAlphaElectrons() const
Definition orbitals.h:93
Eigen::MatrixXd DensityMatrixWithoutGS(const QMState &state) const
Definition orbitals.cc:112
Eigen::MatrixXd DensityMatrixFull(const QMState &state) const
Definition orbitals.cc:140
const tools::EigenSystem & MOs() const
Definition orbitals.h:122
const QMMolecule & QMAtoms() const
Definition orbitals.h:172
bool hasDFTbasisName() const
Definition orbitals.h:198
const std::string & getDFTbasisName() const
Definition orbitals.h:202
const AOBasis & getDftBasis() const
Definition orbitals.h:208
Identifier for QMstates. Strings like S1 are converted into enum +zero indexed int.
Definition qmstate.h:132
#define XTP_LOG(level, log)
Definition logger.h:40
base class for all analysis tools
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26