21#ifndef INCLUDE_LIBINT2_LCAO_MOLDEN_H_
22#define INCLUDE_LIBINT2_LCAO_MOLDEN_H_
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>
37#pragma GCC diagnostic push
38#pragma GCC system_header
40#pragma GCC diagnostic pop
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)
87 basis_(validate(basis)),
89 occupancies_(occupancies),
91 labels_(symmetry_labels),
93 bohr_to_angstrom_(bohr_to_angstrom),
94 coefficient_epsilon_(coefficient_epsilon) {
100 os <<
"[Molden Format]" << std::endl;
105 os <<
"[Atoms] AU" << std::endl;
108 os << std::fixed << std::setprecision(8);
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
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;
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)
158 if (dfg_is_cart_[0]) {
159 if (!dfg_is_cart_[1])
160 os <<
"[7F]" << std::endl;
162 if (!dfg_is_cart_[1]) {
163 os <<
"[5D7F]" << std::endl;
164 }
else if (f_found) {
165 os <<
"[5D10F]" << std::endl;
167 os <<
"[5D]" << std::endl;
170 if (!dfg_is_cart_[2])
171 os <<
"[9G]" << std::endl;
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))
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"))
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;
200 void write(std::ostream& os)
const {
208 void write(
const std::string& filename)
const {
209 std::ofstream os(filename);
216 double bohr_to_angstrom()
const {
217 return bohr_to_angstrom_;
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_;
233 std::vector<std::vector<long>>
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) {
249 throw std::logic_error(
250 "molden::Export cannot handle shells with l > 4");
255 throw std::logic_error(
256 "molden::Export cannot handle solid harmonics p shells");
261 if (!dfg_found[contr.l - 2]) {
262 dfg_is_cart_[contr.l - 2] = !contr.pure;
263 dfg_found[contr.l - 2] =
true;
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");
280 const std::vector<Shell>& validate(
const BasisSet& bs)
const {
281 return validate(bs.shells());
284 void initialize_bf_map() {
285 atom2shell_ = BasisSet::atom2shell(atoms_, basis_);
287 const auto nao = BasisSet::nbf(basis_);
289 assert(nao == coefs_.rows());
290 const auto shell2ao = BasisSet::compute_shell2bf(basis_);
292 for (
size_t iatom = 0; iatom < atoms_.size(); ++iatom) {
293 for (
auto ish : atom2shell_[iatom]) {
294 auto ao = shell2ao[ish];
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;
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;
306 END_FOR_SOLIDHARM_MOLDEN
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;
354 template <
typename Coeffs,
typename Occs,
typename Energies = Eigen::VectorXd>
357 const std::vector<Shell>& basis,
358 const Coeffs& coefficients,
359 const Occs& occupancies,
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)
375 os <<
"[SpaceGroup] (Number)" << std::endl;
376 os << space_group_ << std::endl;
381 os <<
"[Operators]" << std::endl;
382 os <<
"x, y, z" << std::endl;
390 os <<
"[Cell]" << std::endl;
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
415 void write(
const std::string& filename)
const {
416 std::ofstream os(filename);
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