Exports LCAO coefficients in Molden format.
More...
#include <molden.h>
|
| template<typename ShellSequence , typename Coeffs , typename Occs , typename Energies = Eigen::VectorXd> |
| | 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) |
| |
|
void | write_prologue (std::ostream &os) const |
| | writes "[Molden Format]" to ostream os
|
| |
|
void | write_atoms (std::ostream &os) const |
| | writes the "[Atoms]" section to ostream os
|
| |
|
void | write_basis (std::ostream &os) const |
| | writes the "[GTO]" section, as well as optional Cartesian/solid harmonics keywords, to ostream os
|
| |
|
void | write_lcao (std::ostream &os) const |
| | writes the "[MO]" section to ostream os
|
| |
|
void | write (std::ostream &os) const |
| | writes "prologue", atoms, basis, and LCAOs to ostream os
|
| |
|
void | write (const std::string &filename) const |
| | same as write(ostream), but creates new file named filename
|
| |
|
double | bohr_to_angstrom () const |
| |
Exports LCAO coefficients in Molden format.
◆ Export()
template<typename ShellSequence , typename Coeffs , typename Occs , typename Energies = Eigen::VectorXd>
| libint2::molden::Export::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 ) |
|
inline |
- Template Parameters
-
| ShellSequence | BasisSet or std::vector<Shell> |
| Coeffs | the type of LCAO coefficient matrix |
| Energies | the type of LCAO energy vector |
| Occs | the type of LCAO occupancy vector |
- Parameters
-
| atoms | the set of atoms (coordinates in atomic units) |
| basis | a sequence of shells; must meet Molden requirements (see below) |
| coefficients | the matrix of LCAO coefficients (columns are LCAOs, rows are AOs; AOs are ordered according to the order of shells in basis and by the ordering conventions of this Libint configuration) |
| occupancies | the vector of occupancies (size = # LCAOs) |
| energies | the vector of energies (size = # of LCAOs); the default is to assign zero to each LCAO |
| symmetry_labels | the vector of symmetry labels (size = # LCAOs); the default is to assign empty label to each LCAO |
| spincases | the vector of spin cases (size = # LCAOs; true = spin-up or m_s=1/2, false = spin-down or m_s=-1/2); the default is to assign spin-up to each LCAO |
| bohr_to_angstrom | the conversion factor from bohr to angstrom; the default is CODATA 2018 value |
| coefficient_epsilon | omit LCAO coefficients with absolute magnitude smaller than this value; set to 0 to write all coefficients (some Molden parsers, e.g. Avogadro2, require this) |
- Exceptions
-
| std::logic_error | if the basis does not conforms Molden requirements |
- Note
- Molden can only handle basis sets that:
- p (l=1) shells are Cartesian, not solid harmonics
- d, f, and g (l=2..4) shells are all Cartesian or all solid harmonics
- there are no shells with l>5
The documentation for this class was generated from the following file: