LIBINT 2.7.2
molden.h
1/*
2 * Copyright (C) 2004-2021 Edward F. Valeev
3 *
4 * This file is part of Libint.
5 *
6 * Libint is free software: you can redistribute it and/or modify
7 * it under the terms of the GNU Lesser General Public License as published by
8 * the Free Software Foundation, either version 3 of the License, or
9 * (at your option) any later version.
10 *
11 * Libint is distributed in the hope that it will be useful,
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 * GNU Lesser General Public License for more details.
15 *
16 * You should have received a copy of the GNU Lesser General Public License
17 * along with Libint. If not, see <http://www.gnu.org/licenses/>.
18 *
19 */
20
21#ifndef INCLUDE_LIBINT2_LCAO_MOLDEN_H_
22#define INCLUDE_LIBINT2_LCAO_MOLDEN_H_
23
24#include <cmath>
25#include <iomanip>
26#include <ostream>
27#include <string>
28#include <vector>
29
30#include <libint2/atom.h>
31#include <libint2/basis.h>
32#include <libint2/cgshell_ordering.h>
33#include <libint2/chemistry/elements.h>
34#include <libint2/shell.h>
35#include <libint2/shgshell_ordering.h>
36
37#pragma GCC diagnostic push
38#pragma GCC system_header
39#include <Eigen/Core>
40#pragma GCC diagnostic pop
41
42namespace libint2 {
43namespace molden {
44
47class Export {
48 public:
76 template <typename ShellSequence, typename Coeffs, typename Occs,
77 typename Energies = Eigen::VectorXd>
78 Export(const std::vector<Atom>& atoms, const ShellSequence& basis,
79 const Coeffs& coefficients, const Occs& occupancies,
80 const Energies& energies = Energies(),
81 const std::vector<std::string>& symmetry_labels =
82 std::vector<std::string>(),
83 const std::vector<bool>& spincases = std::vector<bool>(),
84 const double bohr_to_angstrom = constants::codata_2018::bohr_to_angstrom,
85 double coefficient_epsilon = 5e-11)
86 : atoms_(atoms),
87 basis_(validate(basis)),
88 coefs_(coefficients),
89 occupancies_(occupancies),
90 energies_(energies),
91 labels_(symmetry_labels),
92 spins_(spincases),
93 bohr_to_angstrom_(bohr_to_angstrom),
94 coefficient_epsilon_(coefficient_epsilon) {
95 initialize_bf_map();
96 }
97
99 void write_prologue(std::ostream& os) const {
100 os << "[Molden Format]" << std::endl;
101 }
102
104 void write_atoms(std::ostream& os) const {
105 os << "[Atoms] AU" << std::endl;
106
107 os.fill(' ');
108 os << std::fixed << std::setprecision(8);
109 auto iatom = 0;
110 for (const auto& atom : atoms_) {
111 auto Z = atom.atomic_number;
112 os << std::setw(4) << libint2::chemistry::get_element_info().at(Z - 1).symbol
113 << std::setw(6) << (iatom + 1) << std::setw(6) << Z
114 << std::setw(14) << atom.x
115 << std::setw(14) << atom.y
116 << std::setw(14) << atom.z
117 << std::endl;
118 ++iatom;
119 }
120 }
121
124 void write_basis(std::ostream& os) const {
125 bool f_found = false;
126 os << "[GTO]" << std::endl;
127 for (size_t iatom = 0; iatom < atoms_.size(); ++iatom) {
128 os << std::setw(4) << (iatom + 1) << std::setw(4) << 0 << std::endl;
129 for (auto ish : atom2shell_[iatom]) {
130 const Shell& sh = basis_.at(ish);
131 if (sh.contr.size() == 1) {
132 const auto& contr = sh.contr[0];
133 const auto l = contr.l;
134 assert(l <= 4); // only up to g functions are supported
135 if (l == 3)
136 f_found = true;
137 const auto nprim = contr.coeff.size();
138 os << std::setw(4) << Shell::am_symbol(contr.l) << std::setw(6)
139 << nprim << std::setw(6) << "1.00" << std::endl;
140 for (int iprim = 0; iprim < nprim; ++iprim) {
141 os << std::scientific << std::uppercase << std::setprecision(10);
142 os << std::setw(20) << sh.alpha[iprim] << std::setw(20)
143 << sh.coeff_normalized(0, iprim) << std::endl;
144 } // end loop over primitives
145 } // end if ncontraction == 1
146 else {
147 assert(false); // Not implemented
148 }
149 } // end loop over shells on center
150 // format calls for a blank line here
151 os << std::endl;
152 } // end loop over centers
153
154 // write solid harmonic/cartesian tags
155 {
156 // Molden default is cartesians throughout
157 // dfg_is_cart_ is set to true even if there are no shells of a given type
158 if (dfg_is_cart_[0]) { // cartesian d
159 if (!dfg_is_cart_[1]) // solid harmonic f
160 os << "[7F]" << std::endl;
161 } else { // solid harmonic d
162 if (!dfg_is_cart_[1]) { // solid harmonic f
163 os << "[5D7F]" << std::endl;
164 } else if (f_found) { // cartesian f
165 os << "[5D10F]" << std::endl;
166 } else { // no f functions
167 os << "[5D]" << std::endl;
168 }
169 }
170 if (!dfg_is_cart_[2]) // solid harmonic g
171 os << "[9G]" << std::endl;
172 }
173 }
174
176 void write_lcao(std::ostream& os) const {
177 os << "[MO]" << std::endl;
178 for (int imo = 0; imo < coefs_.cols(); ++imo) {
179 os << std::fixed << std::setprecision(10);
180 os << std::setw(8) << "Sym= " << (labels_.empty() ? "A" : labels_.at(imo))
181 << std::endl
182 << std::setw(8) << "Ene= " << std::setw(16)
183 << (energies_.rows() == 0 ? 0.0 : energies_(imo)) << std::endl
184 << std::setw(8) << "Spin= "
185 << (spins_.empty() ? "Alpha" : (spins_.at(imo) ? "Alpha" : "Beta"))
186 << std::endl
187 << std::setw(8) << "Occup= " << occupancies_(imo) << std::endl;
188 os << std::scientific << std::uppercase << std::setprecision(10);
189 for (int iao = 0; iao < coefs_.rows(); ++iao) {
190 const auto C_ao_mo = coefs_(ao_map_[iao], imo);
191 if (std::abs(C_ao_mo) >= coefficient_epsilon_) {
192 os << std::setw(6) << (iao + 1) << " " << std::setw(16)
193 << C_ao_mo << std::endl;
194 }
195 } // end loop over AOs
196 } // end loop over MOs
197 }
198
200 void write(std::ostream& os) const {
201 write_prologue(os);
202 write_atoms(os);
203 write_basis(os);
204 write_lcao(os);
205 }
206
208 void write(const std::string& filename) const {
209 std::ofstream os(filename);
210 write_prologue(os);
211 write_atoms(os);
212 write_basis(os);
213 write_lcao(os);
214 }
215
216 double bohr_to_angstrom() const {
217 return bohr_to_angstrom_;
218 }
219
220 private:
221 const std::vector<Atom>& atoms_;
222 const std::vector<Shell>& basis_;
223 Eigen::MatrixXd coefs_;
224 Eigen::VectorXd occupancies_;
225 Eigen::VectorXd energies_;
226 std::vector<std::string> labels_;
227 std::vector<bool> spins_;
228 double bohr_to_angstrom_;
229 double coefficient_epsilon_;
230 mutable bool
231 dfg_is_cart_[3]; // whether {d, f, g} shells are cartesian (true) or
232 // solid harmonics (false)
233 std::vector<std::vector<long>>
234 atom2shell_; // maps atom -> shell indices in basis_
235 std::vector<long>
236 ao_map_; // maps from the AOs ordered according to Molden
237 // (atoms->shells, bf in shells ordered in the Molden order)
238 // to the AOs ordered according to basis_
239
242 const std::vector<Shell>& validate(const std::vector<Shell>& shells) const {
243 bool dfg_found[] = {false, false, false};
244 for(int i=0; i!=sizeof(dfg_is_cart_)/sizeof(bool); ++i)
245 dfg_is_cart_[i] = true;
246 for (const auto& shell : shells) {
247 for (const auto& contr : shell.contr) {
248 if (contr.l > 4)
249 throw std::logic_error(
250 "molden::Export cannot handle shells with l > 4");
251
252 switch (contr.l) {
253 case 1:
254 if (contr.pure)
255 throw std::logic_error(
256 "molden::Export cannot handle solid harmonics p shells");
257 break;
258 case 2:
259 case 3:
260 case 4: {
261 if (!dfg_found[contr.l - 2]) {
262 dfg_is_cart_[contr.l - 2] = !contr.pure;
263 dfg_found[contr.l - 2] = true;
264 }
265 if (!contr.pure ^ dfg_is_cart_[contr.l - 2])
266 throw std::logic_error(
267 "molden::Export only supports all-Cartesian or "
268 "all-solid-harmonics d/f/g shells");
269 }
270
271 default: {} // l = 0 is fine
272 }
273 }
274 }
275 return shells;
276 }
277
280 const std::vector<Shell>& validate(const BasisSet& bs) const {
281 return validate(bs.shells());
282 }
283
284 void initialize_bf_map() {
285 atom2shell_ = BasisSet::atom2shell(atoms_, basis_);
286
287 const auto nao = BasisSet::nbf(basis_);
288 ao_map_.resize(nao);
289 assert(nao == coefs_.rows());
290 const auto shell2ao = BasisSet::compute_shell2bf(basis_);
291 long ao_molden = 0;
292 for (size_t iatom = 0; iatom < atoms_.size(); ++iatom) {
293 for (auto ish : atom2shell_[iatom]) {
294 auto ao = shell2ao[ish]; // refers to order assumed by coefs
295 const auto& shell = basis_[ish];
296 const auto ncontr = shell.contr.size();
297 for (int c = 0; c != ncontr; ++c) {
298 const auto l = shell.contr[c].l;
299 const auto pure = shell.contr[c].pure;
300 if (pure) {
301 int m;
302 FOR_SOLIDHARM_MOLDEN(l, m)
303 const auto ao_in_shell = libint2::INT_SOLIDHARMINDEX(l, m);
304 ao_map_[ao_molden] = ao + ao_in_shell;
305 ++ao_molden;
306 END_FOR_SOLIDHARM_MOLDEN
307 ao += 2 * l + 1;
308 } else {
309 int i, j, k;
310 FOR_CART_MOLDEN(i, j, k, l)
311 const auto ao_in_shell = INT_CARTINDEX(l, i, j);
312 ao_map_[ao_molden] = ao + ao_in_shell;
313 ++ao_molden;
314 END_FOR_CART_MOLDEN
315 ao += INT_NCART(l);
316 }
317 } // contraction loop
318 }
319 }
320 }
321
322}; // Export
323
326class PBCExport: public Export{
327 public:
354 template <typename Coeffs, typename Occs, typename Energies = Eigen::VectorXd>
355 PBCExport(const std::vector<Atom>& atoms,
356 const std::array<Eigen::Vector3d, 3>& cell_axes,
357 const std::vector<Shell>& basis,
358 const Coeffs& coefficients,
359 const Occs& occupancies,
360 int space_group,
361 const Energies& energies = Energies(),
362 const std::vector<std::string>& symmetry_labels =
363 std::vector<std::string>(),
364 const std::vector<bool>& spincases = std::vector<bool>(),
365 const double bohr_to_angstrom = constants::codata_2018::bohr_to_angstrom)
366 : Export(atoms, basis, coefficients, occupancies, energies, symmetry_labels, spincases, bohr_to_angstrom),
367 cell_axes_(cell_axes),
368 space_group_(space_group)
369 {
370 //initialize_bf_map();
371 }
372
374 void write_space_group(std::ostream& os) const {
375 os << "[SpaceGroup] (Number)" << std::endl;
376 os << space_group_ << std::endl;
377 }
378
380 void write_operators(std::ostream& os) const {
381 os << "[Operators]" << std::endl;
382 os << "x, y, z" << std::endl;
383 }
384
386 void write_cell_axes(std::ostream& os) const {
387
388 // https://sourceforge.net/p/jmol/code/HEAD/tree/trunk/Jmol/src/org/jmol/adapter/readers/quantum/MoldenReader.java#l107
389 // suggests that [Cell] defaults to angstroms
390 os << "[Cell]" << std::endl;
391 {
392 // convert vectors to abcɑβɣ
393 const double a = cell_axes_[0].norm();
394 const double b = cell_axes_[1].norm();
395 const double c = cell_axes_[2].norm();
396 const double eps = std::sqrt(std::numeric_limits<double>::epsilon());
397 const bool nonzero_a = a >= eps;
398 const bool nonzero_b = b >= eps;
399 const bool nonzero_c = c >= eps;
400 static constexpr double right_angle = M_PI / 2;
401 const double alpha = nonzero_b && nonzero_c ? std::acos(cell_axes_[1].dot(cell_axes_[2]) / (b * c)) : right_angle;
402 const double beta = nonzero_a && nonzero_c ? std::acos(cell_axes_[0].dot(cell_axes_[2]) / (a * c)) : right_angle;
403 const double gamma = nonzero_a && nonzero_b ? std::acos(cell_axes_[0].dot(cell_axes_[1]) / (a * b)) : right_angle;
404 const double radian_to_degree = 180 / M_PI;
405 os << std::setw(12) << a * bohr_to_angstrom()
406 << std::setw(12) << b * bohr_to_angstrom()
407 << std::setw(12) << c * bohr_to_angstrom()
408 << std::setw(12) << alpha * radian_to_degree
409 << std::setw(12) << beta * radian_to_degree
410 << std::setw(12) << gamma * radian_to_degree
411 << std::endl;
412 }
413 }
414
415 void write(const std::string& filename) const {
416 std::ofstream os(filename);
417 write_prologue(os);
419 write_operators(os);
420 write_cell_axes(os);
421 write_atoms(os);
422 write_basis(os);
423 write_lcao(os);
424 }
425
426private:
428 int space_group_;
429
430}; // PBCExport
431
432} // namespace molden
433} // namespace libint2
434
435#endif // INCLUDE_LIBINT2_LCAO_MOLDEN_H_
Exports LCAO coefficients in Molden format.
Definition: molden.h:47
void write_basis(std::ostream &os) const
writes the "[GTO]" section, as well as optional Cartesian/solid harmonics keywords,...
Definition: molden.h:124
void write_atoms(std::ostream &os) const
writes the "[Atoms]" section to ostream os
Definition: molden.h:104
void write_prologue(std::ostream &os) const
writes "[Molden Format]" to ostream os
Definition: molden.h:99
void write(const std::string &filename) const
same as write(ostream), but creates new file named filename
Definition: molden.h:208
Export(const std::vector< Atom > &atoms, const ShellSequence &basis, const Coeffs &coefficients, const Occs &occupancies, const Energies &energies=Energies(), const std::vector< std::string > &symmetry_labels=std::vector< std::string >(), const std::vector< bool > &spincases=std::vector< bool >(), const double bohr_to_angstrom=constants::codata_2018::bohr_to_angstrom, double coefficient_epsilon=5e-11)
Definition: molden.h:78
void write_lcao(std::ostream &os) const
writes the "[MO]" section to ostream os
Definition: molden.h:176
void write(std::ostream &os) const
writes "prologue", atoms, basis, and LCAOs to ostream os
Definition: molden.h:200
Extension of the Molden exporter to support JMOL extensions for crystal orbitals (see https://sourcef...
Definition: molden.h:326
void write_space_group(std::ostream &os) const
writes the "[SpaceGroup]" section to ostream os
Definition: molden.h:374
PBCExport(const std::vector< Atom > &atoms, const std::array< Eigen::Vector3d, 3 > &cell_axes, const std::vector< Shell > &basis, const Coeffs &coefficients, const Occs &occupancies, int space_group, const Energies &energies=Energies(), const std::vector< std::string > &symmetry_labels=std::vector< std::string >(), const std::vector< bool > &spincases=std::vector< bool >(), const double bohr_to_angstrom=constants::codata_2018::bohr_to_angstrom)
Definition: molden.h:355
void write_cell_axes(std::ostream &os) const
writes the "[Cell]" section to ostream os
Definition: molden.h:386
void write_operators(std::ostream &os) const
writes the "[Operators]" section to ostream os
Definition: molden.h:380
Defaults definitions for various parameters assumed by Libint.
Definition: algebra.cc:24
generally-contracted Solid-Harmonic/Cartesion Gaussian Shell
Definition: shell.h:81
svector< Contraction > contr
contractions
Definition: shell.h:105
real_t coeff_normalized(size_t c, size_t p) const
Definition: shell.h:231
svector< real_t > alpha
exponents
Definition: shell.h:104