27#include <boost/algorithm/string.hpp>
28#include <boost/format.hpp>
62 std::string fileName = options.
get(
"temporary_file").
as<std::string>();
75 std::string& el_file_name) {
85 el_file.open(el_file_name);
86 el_file <<
"$DATA" << endl;
88 for (
const std::string& element_name : UniqueElements) {
90 el_file << elementInfo.
getEleFull(element_name) << endl;
91 for (
const Shell& shell : element) {
97 el_file <<
" " << sh_idx <<
" " <<
indent(gaussian.decay());
98 el_file <<
" " <<
indent(gaussian.contraction());
115 for (
const QMAtom& atom : qmatoms) {
117 inp_file << setw(3) << atom.getElement() << setw(12)
118 << setiosflags(ios::fixed) << setprecision(6) << pos.x()
119 << setw(12) << setiosflags(ios::fixed) << setprecision(6)
120 << pos.y() << setw(12) << setiosflags(ios::fixed)
121 << setprecision(6) << pos.z() << endl;
123 inp_file <<
"* \n" << endl;
139 <<
options_.get(
"ecp").as<std::string>() << flush;
141 for (
const std::string& element_name : UniqueElements) {
144 }
catch (std::runtime_error& error) {
146 <<
"No pseudopotential for " << element_name <<
" available" << flush;
153 <<
" " << element_name << endl;
155 <<
" " << element.getNcore() << endl;
160 for (
Index i = 0; i <=
Index(element.getLmax()); i++) {
161 for (
const ECPShell& shell : element) {
162 if (
Index(shell.getL()) == i) {
169 inp_file << sh_idx <<
" " << gaussian.decay_ <<
" "
170 << gaussian.contraction_ <<
" " << gaussian.power_ << endl;
183 this->
options_.addTree(
"orca.pointcharges",
"\"background.crg\"");
192 std::ofstream crg_file;
193 std::string crg_file_name_full_ =
run_dir_ +
"/background.crg";
194 crg_file.open(crg_file_name_full_);
195 Index total_background = 0;
198 if (site->getCharge() != 0.0) {
204 crg_file << total_background << endl;
205 boost::format fmt(
"%1$+1.7f %2$+1.7f %3$+1.7f %4$+1.7f");
210 boost::str(fmt % site->getCharge() % pos.x() % pos.y() % pos.z());
211 if (site->getCharge() != 0.0) {
212 crg_file << sitestring << endl;
214 std::vector<MinimalMMCharge> split_multipoles =
SplitMultipoles(*site);
215 for (
const auto& mpoles : split_multipoles) {
218 boost::str(fmt % mpoles.q_ % pos2.x() % pos2.y() % pos2.z());
219 crg_file << multipole << endl;
233 std::vector<std::string> results;
234 std::string temp_suffix =
"/id";
236 std::ofstream inp_file;
238 inp_file.open(inp_file_name_full);
240 inp_file <<
"* xyz " <<
charge_ <<
" " <<
spin_ << endl;
247 <<
"nprocs " << threads <<
"\nend"
251 std::string el_file_name =
run_dir_ +
"/" +
"system.bas";
253 inp_file <<
"%basis\n";
254 inp_file <<
"GTOName"
257 <<
"\"system.bas\";" << endl;
258 if (
options_.exists(
"auxbasisset")) {
259 std::string aux_file_name =
run_dir_ +
"/" +
"system.aux";
260 std::string auxbasisset_name =
261 options_.get(
"auxbasisset").as<std::string>();
263 inp_file <<
"GTOAuxName"
266 <<
"\"system.aux\";" << endl;
281 if (
options_.exists(
"externalfield")) {
282 Eigen::Vector3d field =
options_.get(
"externalfield").as<Eigen::Vector3d>();
283 inp_file <<
"%scf\n ";
284 inp_file <<
" efield " << field.x() <<
", " << field.y() <<
", "
285 << field.z() <<
"\n";
287 inp_file << std::endl;
289 std::string input_options;
291 for (
const auto& prop : this->
options_.get(
"orca")) {
292 const std::string& prop_name = prop.name();
293 if (prop_name ==
"pointcharges") {
295 }
else if (prop_name !=
"method") {
301 inp_file << input_options;
307 <<
"Setting the scratch dir to " <<
scratch_dir_ + temp_suffix << flush;
317 shell_file.open(shell_file_name_full);
318 shell_file <<
"#!/usr/bin/env bash" << endl;
320 const std::string key =
"molecule";
322 auto pos = path.find(key);
324 if (pos != std::string::npos) {
332 if (
options_.get(
"initial_guess").as<std::string>() ==
"orbfile") {
333 if (!(std::filesystem::exists(
run_dir_ +
"/molA.gbw") &&
334 std::filesystem::exists(
run_dir_ +
"/molB.gbw"))) {
336 "Using guess relies on a molA.gbw and a molB.gbw file being in the "
339 shell_file <<
options_.get(
"executable").as<std::string>()
340 <<
"_mergefrag molA.gbw molB.gbw " << base_name
341 <<
".gbw > merge.log" << endl;
343 shell_file <<
options_.get(
"executable").as<std::string>() <<
" "
346 shell_file <<
options_.get(
"executable").as<std::string>() <<
"_2mkl "
347 << base_name <<
" -molden > molden.log" << endl;
356 std::vector<char*> argv = {
const_cast<char*
>(
"/bin/sh"),
357 const_cast<char*
>(
"-c"),
358 const_cast<char*
>(command.c_str()),
nullptr};
360 int rc = posix_spawn(&pid,
"/bin/sh",
nullptr,
nullptr, argv.data(),
environ);
363 throw std::runtime_error(
"posix_spawn failed: " + std::to_string(rc));
367 waitpid(pid, &status, 0);
369 if (WIFEXITED(status))
return WEXITSTATUS(status);
381 if (std::system(
nullptr)) {
410 if (
options_.get(
"initial_guess").as<std::string>() ==
"orbfile") {
411 remove((
run_dir_ +
"/" +
"molA.gbw").c_str());
412 remove((
run_dir_ +
"/" +
"molB.gbw").c_str());
413 remove((
run_dir_ +
"/" +
"dimer.gbw").c_str());
418 for (
const std::string& substring : tok_cleanup) {
419 if (substring ==
"inp") {
421 remove(file_name.c_str());
424 if (substring ==
"bas") {
425 std::string file_name =
run_dir_ +
"/system.bas";
426 remove(file_name.c_str());
429 if (substring ==
"log") {
431 remove(file_name.c_str());
434 if (substring ==
"gbw") {
436 remove(file_name.c_str());
439 if (substring ==
"ges") {
440 std::string file_name =
run_dir_ +
"/system.ges";
441 remove(file_name.c_str());
443 if (substring ==
"prop") {
444 std::string file_name =
run_dir_ +
"/system.prop";
445 remove(file_name.c_str());
460 std::ifstream input_file(log_file_name_full);
466 std::string::size_type charge_pos = line.find(
"CHELPG Charges");
468 if (charge_pos != std::string::npos) {
473 bool hasAtoms = result.
size() > 0;
474 while (nfields == 4) {
475 Index atom_id = boost::lexical_cast<Index>(row.at(0));
476 std::string atom_type = row.at(1);
477 double atom_charge = boost::lexical_cast<double>(row.at(3));
479 nfields =
Index(row.size());
483 throw std::runtime_error(
484 "Getting charges failed. Mismatch in elemts:" +
490 StaticSite(atom_id, atom_type, Eigen::Vector3d::Zero());
503 bool has_pol =
false;
505 Eigen::Matrix3d pol = Eigen::Matrix3d::Zero();
510 std::string::size_type pol_pos = line.find(
"POLARIZABILITY TENSOR");
512 if (pol_pos != std::string::npos) {
515 for (
Index i = 0; i < 10; i++) {
517 if (line.find(
"The raw cartesian tensor (atomic units)") !=
522 throw std::runtime_error(
523 "Could not find cartesian polarization tensor");
527 for (
Index i = 0; i < 3; i++) {
529 std::vector<double> values =
531 if (values.size() != 3) {
532 throw std::runtime_error(
"polarization line " + line +
533 " cannot be parsed");
536 row << values[0], values[1], values[2];
544 throw std::runtime_error(
"Could not find polarization in logfile");
550 bool found_success =
false;
563 std::map<Index, double> energies;
564 std::map<Index, double> energies_beta;
566 std::map<Index, double> occupancy;
567 std::map<Index, double> occupancy_beta;
570 std::string orca_version;
571 Index orca_major_version = 0;
574 Index number_of_electrons = 0;
575 Index number_of_electrons_beta = 0;
576 Index number_of_virtuals = 0;
577 Index number_of_virtuals_beta = 0;
579 std::vector<std::string> results;
581 std::ifstream input_file(log_file_name_full);
583 if (input_file.fail()) {
585 <<
"File " << log_file_name_full <<
" not found " << flush;
589 <<
"Reading basic ORCA output from " << log_file_name_full << flush;
602 std::string::size_type version_pos = line.find(
"Program Version");
603 if (version_pos != std::string::npos) {
605 orca_version = results[2];
606 boost::trim(orca_version);
609 orca_major_version = boost::lexical_cast<Index>(results[0]);
611 <<
"ORCA Major Version " << orca_major_version << flush;
614 std::string::size_type energy_pos = line.find(
"FINAL SINGLE");
615 if (energy_pos != std::string::npos) {
617 std::string energy = results[4];
619 orbitals.
setQMEnergy(boost::lexical_cast<double>(energy));
626 std::string::size_type HFX_pos = line.find(
"Fraction HF Exchange ScalHFX");
628 if (HFX_pos != std::string::npos) {
630 double ScaHFX = boost::lexical_cast<double>(results.back());
634 <<
"DFT with " << ScaHFX <<
" of HF exchange!" << flush;
637 std::string::size_type dim_pos = line.find(
"Basis Dimension");
638 if (dim_pos != std::string::npos) {
643 levels = boost::lexical_cast<Index>(dim);
648 std::string::size_type OE_pos = line.find(
"ORBITAL ENERGIES");
649 if (OE_pos != std::string::npos) {
651 number_of_electrons = 0;
656 line.find(
"SPIN UP ORBITALS") == std::string::npos) {
658 "Expected to read an open-shell system but found no spin orbitals");
661 if (line.find(
"E(Eh)") == std::string::npos) {
663 <<
"Warning: Orbital Energies not found in log file" << flush;
665 for (
Index i = 0; i < levels; i++) {
666 if (number_of_virtuals == 10 && orca_major_version > 5) {
670 std::string no = results[0];
672 Index levelnumber = boost::lexical_cast<Index>(no);
673 if (levelnumber != i) {
675 "something weird is going on"
678 std::string oc = results[1];
680 double occ = boost::lexical_cast<double>(oc);
683 if (occ == 2 || occ == 1) {
684 number_of_electrons++;
686 }
else if (occ == 0) {
687 number_of_virtuals++;
690 std::string
e = results[2];
692 energies[i] = boost::lexical_cast<double>(
e);
697 number_of_electrons_beta = 0;
698 number_of_virtuals_beta = 0;
702 for (
Index i = 0; i < levels; i++) {
703 if (number_of_virtuals == 10 && orca_major_version > 5) {
707 std::string no = results[0];
709 Index levelnumber = boost::lexical_cast<Index>(no);
710 if (levelnumber != i) {
712 <<
"Have a look at the orbital energies "
713 "something weird is going on"
716 std::string oc = results[1];
718 double occ = boost::lexical_cast<double>(oc);
721 number_of_electrons_beta++;
722 occupancy_beta[i] = occ;
723 }
else if (occ == 0) {
724 number_of_virtuals_beta++;
725 occupancy_beta[i] = occ;
728 "Encountered spin down orbital with occupancy != 0 or 1");
730 std::string
e = results[2];
732 energies_beta[i] = boost::lexical_cast<double>(
e);
737 std::string::size_type success =
738 line.find(
"* SUCCESS *");
739 if (success != std::string::npos) {
740 found_success =
true;
750 <<
"Alpha electrons: " << number_of_electrons << flush;
751 Index occupied_levels = number_of_electrons;
752 Index unoccupied_levels = levels - occupied_levels;
754 <<
"Occupied Alpha levels: " << occupied_levels << flush;
756 <<
"Unoccupied levels: " << unoccupied_levels << flush;
767 for (
Index i = 0; i < number_of_electrons + number_of_virtuals; i++) {
773 <<
"Beta electrons: " << number_of_electrons_beta << flush;
775 <<
"Occupied Beta levels: " << number_of_electrons_beta << flush;
777 <<
"Unoccupied Beta levels: " << levels - number_of_electrons_beta
783 for (
Index i = 0; i < number_of_electrons_beta + number_of_virtuals_beta;
791 return found_success;
795 std::string::size_type coordinates_pos =
796 line.find(
"CARTESIAN COORDINATES (ANGSTROEM)");
798 using Atom =
typename T::Atom_Type;
800 if (coordinates_pos != std::string::npos) {
802 bool has_QMAtoms = mol.size() > 0;
810 while (nfields == 4) {
811 string atom_type = row.at(0);
812 double x = boost::lexical_cast<double>(row.at(1));
813 double y = boost::lexical_cast<double>(row.at(2));
814 double z = boost::lexical_cast<double>(row.at(3));
816 nfields =
Index(row.size());
817 Eigen::Vector3d pos(x, y, z);
819 if (has_QMAtoms ==
false) {
820 mol.push_back(
Atom(atom_id, atom_type, pos));
822 Atom& pAtom = mol.at(atom_id);
834 if (input_file.fail()) {
843 std::string::size_type error = line.find(
"FATAL ERROR ENCOUNTERED");
844 if (error != std::string::npos) {
846 "look in the log file may help."
851 "mpirun detected that one or more processes exited with non-zero "
853 if (error != std::string::npos) {
855 <<
"ORCA had an mpi problem, maybe your openmpi version is not good."
868 std::vector<double> coefficients;
870 if (basis_size == 0) {
872 "Basis size not set, calculator does not parse log file first");
881 "Basisset names should be set before reading the molden file.");
885 std::string file_name =
run_dir_ +
"/" +
889 std::ifstream molden_file(file_name);
890 if (!molden_file.good()) {
891 throw std::runtime_error(
892 "Could not find the molden input file for the MO coefficients.\nIf you "
893 "have run the orca calculation manually or use data from an old\n"
894 "calculation, make sure that besides the .gbw file a .molden.input is\n"
895 "present. If not, convert the .gbw file to a .molden.input file with\n"
896 "the orca_2mkl tool from orca.\nAn example, if you have a benzene.gbw "
897 "file run:\n orca_2mkl benzene -molden\n");
915 std::stringstream ssnumber;
921 ssnumber << setiosflags(ios::fixed) << setprecision(15) << std::scientific
923 std::string snumber = ssnumber.str();
928 std::stringstream stream;
929 std::string section = key.substr(key.find(
".") + 1);
930 stream <<
"%" << section;
932 stream <<
" " <<
options_.get(key).as<std::string>() <<
"\n";
935 <<
options_.get(key).as<std::string>() <<
"\n"
944 std::vector<std::string> words = values.
ToVector();
945 return ((words.size() <= 1) ?
true :
false);
949 std::stringstream stream;
950 std::string opt = (
options_.get(
"optimize").as<bool>()) ?
"Opt" :
"";
952 std::string user_method =
953 (orca.
exists(
"method")) ? orca.
get(
"method").
as<std::string>() :
"";
954 std::string convergence =
"";
955 if (!orca.
exists(
"scf")) {
957 options_.get(
"convergence_tightness").as<std::string>());
963 << user_method <<
"\n";
969 std::map<std::string, std::string> votca_to_orca;
971 votca_to_orca[
"XC_HYB_GGA_XC_B1LYP"] =
"B1LYP";
972 votca_to_orca[
"XC_HYB_GGA_XC_B3LYP"] =
"B3LYP";
973 votca_to_orca[
"XC_HYB_GGA_XC_PBEH"] =
"PBE0";
974 votca_to_orca[
"XC_GGA_C_PBE"] =
"PBE";
975 votca_to_orca[
"XC_GGA_X_PBE"] =
"PBE";
977 std::string votca_functional =
978 options_.get(
"functional").as<std::vector<std::string>>()[0];
980 std::string orca_functional;
981 if (votca_to_orca.count(votca_functional)) {
982 orca_functional = votca_to_orca.at(votca_functional);
983 }
else if (
options_.exists(
"orca." + votca_functional)) {
985 options_.get(
"orca." + votca_functional).as<std::string>();
987 throw std::runtime_error(
988 "Cannot translate " + votca_functional +
989 " to orca functional names. Please add a <" + votca_functional +
990 "> to your orca input options with the functional orca should use.");
992 return orca_functional;
void push_back(const T &atom)
const T & at(Index index) const
std::vector< std::string > FindUniqueElements() const
void setPos(const Eigen::Vector3d &r)
const Element & getElement(std::string element_type) const
void Load(const std::string &name)
Container to hold ECPs for all atoms.
std::vector< std::string > Fill(const ECPBasisSet &bs, QMMolecule &atoms)
void Load(const std::string &name)
const ECPElement & getElement(std::string element_type) const
void parseMoldenFile(const std::string &filename, Orbitals &orbitals) const
void setBasissetInfo(const std::string &basisset_name, const std::string &aux_basisset_name="")
Container for molecular orbitals and derived one-particle data.
void setScaHFX(double ScaHFX)
Store the fraction of exact exchange associated with the functional.
double getDFTTotalEnergy() const
Return the stored total DFT energy.
const tools::EigenSystem & MOs_beta() const
Return read-only access to beta-spin molecular orbitals.
void setNumberOfAlphaElectrons(Index electrons)
Store the total number of alpha electrons.
void setNumberOfBetaElectrons(Index electrons)
Store the total number of beta electrons.
void setECPName(const std::string &ECP)
Store the effective core potential label.
void setXCGrid(std::string grid)
Store the numerical XC grid quality label.
void setNumberOfOccupiedLevels(Index occupied_levels)
Index getBasisSetSize() const
Return the number of AO basis functions in the DFT basis.
void setQMEnergy(double qmenergy)
Store the total DFT energy.
const tools::EigenSystem & MOs() const
Return read-only access to alpha/restricted molecular orbitals.
void setQMpackage(const std::string &qmpackage)
Store the name of the QM package that produced these orbitals.
const QMMolecule & QMAtoms() const
Return read-only access to the molecular geometry.
void setNumberOfOccupiedLevelsBeta(Index occupied_levels_beta)
Store the number of occupied beta-spin orbitals.
void setChargeAndSpin(Index charge, Index spin)
const std::string & getECPName() const
Return the effective core potential label.
const std::string & getDFTbasisName() const
Return the DFT basis-set name.
void SetupDftBasis(std::string basis_name)
Build and attach the DFT AO basis from the stored molecular geometry.
bool isOpenShell() const
Report whether the stored state corresponds to an open-shell system.
void setXCFunctionalName(std::string functionalname)
void GetCoordinates(T &mol, std::string &line, std::ifstream &input_file) const
void WriteCoordinates(std::ofstream &inp_file, const QMMolecule &)
void WriteECP(std::ofstream &inp_file, const QMMolecule &)
void WriteBasisset(const QMMolecule &qmatoms, std::string &bs_name, std::string &el_file_name)
std::string GetOrcaFunctionalName() const
std::string indent(const double &number)
std::string getPackageName() const override
Eigen::Matrix3d GetPolarizability() const override
void WriteChargeOption() override
std::map< std::string, std::string > convergence_map_
std::string CreateInputSection(const std::string &key) const
StaticSegment GetCharges() const override
std::string WriteMethod() const
void WriteBackgroundCharges()
bool ParseLogFile(Orbitals &orbitals) override
void ParseSpecificOptions(const tools::Property &options) final
bool ParseMOsFile(Orbitals &orbitals) override
bool KeywordIsSingleLine(const std::string &key) const
bool WriteInputFile(const Orbitals &orbitals) override
std::string log_file_name_
std::vector< std::string > GetLineAndSplit(std::ifstream &input_file, const std::string separators) const
std::vector< std::unique_ptr< StaticSite > > externalsites_
std::vector< MinimalMMCharge > SplitMultipoles(const StaticSite &site) const
std::string mo_file_name_
std::string basisset_name_
std::string shell_file_name_
std::string input_file_name_
Class to represent Atom/Site in electrostatic.
const std::string & getElement() const
#define XTP_LOG(level, log)
Charge transport classes.
std::string EnumToString(L l)
ClassicalSegment< StaticSite > StaticSegment
int run_command_spawn(const std::string &command)
Provides a means for comparing floating point numbers.