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;
359 std::vector<char*> argv = {
360 const_cast<char*
>(
"/bin/sh"),
361 const_cast<char*
>(
"-c"),
362 const_cast<char*
>(command.c_str()),
366 int rc = posix_spawn(
376 throw std::runtime_error(
"posix_spawn failed: " + std::to_string(rc));
380 waitpid(pid, &status, 0);
382 if (WIFEXITED(status))
383 return WEXITSTATUS(status);
397 if (std::system(
nullptr)) {
426 if (
options_.get(
"initial_guess").as<std::string>() ==
"orbfile") {
427 remove((
run_dir_ +
"/" +
"molA.gbw").c_str());
428 remove((
run_dir_ +
"/" +
"molB.gbw").c_str());
429 remove((
run_dir_ +
"/" +
"dimer.gbw").c_str());
434 for (
const std::string& substring : tok_cleanup) {
435 if (substring ==
"inp") {
437 remove(file_name.c_str());
440 if (substring ==
"bas") {
441 std::string file_name =
run_dir_ +
"/system.bas";
442 remove(file_name.c_str());
445 if (substring ==
"log") {
447 remove(file_name.c_str());
450 if (substring ==
"gbw") {
452 remove(file_name.c_str());
455 if (substring ==
"ges") {
456 std::string file_name =
run_dir_ +
"/system.ges";
457 remove(file_name.c_str());
459 if (substring ==
"prop") {
460 std::string file_name =
run_dir_ +
"/system.prop";
461 remove(file_name.c_str());
476 std::ifstream input_file(log_file_name_full);
482 std::string::size_type charge_pos = line.find(
"CHELPG Charges");
484 if (charge_pos != std::string::npos) {
489 bool hasAtoms = result.
size() > 0;
490 while (nfields == 4) {
491 Index atom_id = boost::lexical_cast<Index>(row.at(0));
492 std::string atom_type = row.at(1);
493 double atom_charge = boost::lexical_cast<double>(row.at(3));
495 nfields =
Index(row.size());
499 throw std::runtime_error(
500 "Getting charges failed. Mismatch in elemts:" +
506 StaticSite(atom_id, atom_type, Eigen::Vector3d::Zero());
519 bool has_pol =
false;
521 Eigen::Matrix3d pol = Eigen::Matrix3d::Zero();
526 std::string::size_type pol_pos = line.find(
"POLARIZABILITY TENSOR");
528 if (pol_pos != std::string::npos) {
531 for (
Index i = 0; i < 10; i++){
533 if (line.find(
"The raw cartesian tensor (atomic units)") !=
538 throw std::runtime_error(
539 "Could not find cartesian polarization tensor");
543 for (
Index i = 0; i < 3; i++) {
545 std::vector<double> values =
547 if (values.size() != 3) {
548 throw std::runtime_error(
"polarization line " + line +
549 " cannot be parsed");
552 row << values[0], values[1], values[2];
560 throw std::runtime_error(
"Could not find polarization in logfile");
566 bool found_success =
false;
579 std::map<Index, double> energies;
580 std::map<Index, double> energies_beta;
582 std::map<Index, double> occupancy;
583 std::map<Index, double> occupancy_beta;
586 std::string orca_version;
587 Index orca_major_version = 0;
590 Index number_of_electrons = 0;
591 Index number_of_electrons_beta = 0;
592 Index number_of_virtuals = 0;
593 Index number_of_virtuals_beta = 0;
595 std::vector<std::string> results;
597 std::ifstream input_file(log_file_name_full);
599 if (input_file.fail()) {
601 <<
"File " << log_file_name_full <<
" not found " << flush;
605 <<
"Reading basic ORCA output from " << log_file_name_full << flush;
618 std::string::size_type version_pos = line.find(
"Program Version");
619 if (version_pos != std::string::npos) {
621 orca_version = results[2];
622 boost::trim(orca_version);
625 orca_major_version = boost::lexical_cast<Index>(results[0]);
627 <<
"ORCA Major Version " << orca_major_version << flush;
630 std::string::size_type energy_pos = line.find(
"FINAL SINGLE");
631 if (energy_pos != std::string::npos) {
633 std::string energy = results[4];
635 orbitals.
setQMEnergy(boost::lexical_cast<double>(energy));
642 std::string::size_type HFX_pos = line.find(
"Fraction HF Exchange ScalHFX");
644 if (HFX_pos != std::string::npos) {
646 double ScaHFX = boost::lexical_cast<double>(results.back());
650 <<
"DFT with " << ScaHFX <<
" of HF exchange!" << flush;
653 std::string::size_type dim_pos = line.find(
"Basis Dimension");
654 if (dim_pos != std::string::npos) {
659 levels = boost::lexical_cast<Index>(dim);
664 std::string::size_type OE_pos = line.find(
"ORBITAL ENERGIES");
665 if (OE_pos != std::string::npos) {
667 number_of_electrons = 0;
672 line.find(
"SPIN UP ORBITALS") == std::string::npos) {
674 "Expected to read an open-shell system but found no spin orbitals");
677 if (line.find(
"E(Eh)") == std::string::npos) {
679 <<
"Warning: Orbital Energies not found in log file" << flush;
681 for (
Index i = 0; i < levels; i++) {
682 if (number_of_virtuals == 10 && orca_major_version > 5) {
686 std::string no = results[0];
688 Index levelnumber = boost::lexical_cast<Index>(no);
689 if (levelnumber != i) {
691 "something weird is going on"
694 std::string oc = results[1];
696 double occ = boost::lexical_cast<double>(oc);
699 if (occ == 2 || occ == 1) {
700 number_of_electrons++;
702 }
else if (occ == 0) {
703 number_of_virtuals++;
706 std::string
e = results[2];
708 energies[i] = boost::lexical_cast<double>(
e);
713 number_of_electrons_beta = 0;
714 number_of_virtuals_beta = 0;
718 for (
Index i = 0; i < levels; i++) {
719 if (number_of_virtuals == 10 && orca_major_version > 5) {
723 std::string no = results[0];
725 Index levelnumber = boost::lexical_cast<Index>(no);
726 if (levelnumber != i) {
728 <<
"Have a look at the orbital energies "
729 "something weird is going on"
732 std::string oc = results[1];
734 double occ = boost::lexical_cast<double>(oc);
737 number_of_electrons_beta++;
738 occupancy_beta[i] = occ;
739 }
else if (occ == 0) {
740 number_of_virtuals_beta++;
741 occupancy_beta[i] = occ;
744 "Encountered spin down orbital with occupancy != 0 or 1");
746 std::string
e = results[2];
748 energies_beta[i] = boost::lexical_cast<double>(
e);
753 std::string::size_type success =
754 line.find(
"* SUCCESS *");
755 if (success != std::string::npos) {
756 found_success =
true;
766 <<
"Alpha electrons: " << number_of_electrons << flush;
767 Index occupied_levels = number_of_electrons;
768 Index unoccupied_levels = levels - occupied_levels;
770 <<
"Occupied Alpha levels: " << occupied_levels << flush;
772 <<
"Unoccupied levels: " << unoccupied_levels << flush;
783 for (
Index i = 0; i < number_of_electrons + number_of_virtuals; i++) {
789 <<
"Beta electrons: " << number_of_electrons_beta << flush;
791 <<
"Occupied Beta levels: " << number_of_electrons_beta << flush;
793 <<
"Unoccupied Beta levels: " << levels - number_of_electrons_beta
799 for (
Index i = 0; i < number_of_electrons_beta + number_of_virtuals_beta;
807 return found_success;
811 std::string::size_type coordinates_pos =
812 line.find(
"CARTESIAN COORDINATES (ANGSTROEM)");
814 using Atom =
typename T::Atom_Type;
816 if (coordinates_pos != std::string::npos) {
818 bool has_QMAtoms = mol.size() > 0;
826 while (nfields == 4) {
827 string atom_type = row.at(0);
828 double x = boost::lexical_cast<double>(row.at(1));
829 double y = boost::lexical_cast<double>(row.at(2));
830 double z = boost::lexical_cast<double>(row.at(3));
832 nfields =
Index(row.size());
833 Eigen::Vector3d pos(x, y, z);
835 if (has_QMAtoms ==
false) {
836 mol.push_back(
Atom(atom_id, atom_type, pos));
838 Atom& pAtom = mol.at(atom_id);
850 if (input_file.fail()) {
859 std::string::size_type error = line.find(
"FATAL ERROR ENCOUNTERED");
860 if (error != std::string::npos) {
862 "look in the log file may help."
867 "mpirun detected that one or more processes exited with non-zero "
869 if (error != std::string::npos) {
871 <<
"ORCA had an mpi problem, maybe your openmpi version is not good."
884 std::vector<double> coefficients;
886 if (basis_size == 0) {
888 "Basis size not set, calculator does not parse log file first");
897 "Basisset names should be set before reading the molden file.");
901 std::string file_name =
run_dir_ +
"/" +
905 std::ifstream molden_file(file_name);
906 if (!molden_file.good()) {
907 throw std::runtime_error(
908 "Could not find the molden input file for the MO coefficients.\nIf you "
909 "have run the orca calculation manually or use data from an old\n"
910 "calculation, make sure that besides the .gbw file a .molden.input is\n"
911 "present. If not, convert the .gbw file to a .molden.input file with\n"
912 "the orca_2mkl tool from orca.\nAn example, if you have a benzene.gbw "
913 "file run:\n orca_2mkl benzene -molden\n");
931 std::stringstream ssnumber;
937 ssnumber << setiosflags(ios::fixed) << setprecision(15) << std::scientific
939 std::string snumber = ssnumber.str();
944 std::stringstream stream;
945 std::string section = key.substr(key.find(
".") + 1);
946 stream <<
"%" << section;
948 stream <<
" " <<
options_.get(key).as<std::string>() <<
"\n";
951 <<
options_.get(key).as<std::string>() <<
"\n"
960 std::vector<std::string> words = values.
ToVector();
961 return ((words.size() <= 1) ?
true :
false);
965 std::stringstream stream;
966 std::string opt = (
options_.get(
"optimize").as<bool>()) ?
"Opt" :
"";
968 std::string user_method =
969 (orca.
exists(
"method")) ? orca.
get(
"method").
as<std::string>() :
"";
970 std::string convergence =
"";
971 if (!orca.
exists(
"scf")) {
973 options_.get(
"convergence_tightness").as<std::string>());
979 << user_method <<
"\n";
985 std::map<std::string, std::string> votca_to_orca;
987 votca_to_orca[
"XC_HYB_GGA_XC_B1LYP"] =
"B1LYP";
988 votca_to_orca[
"XC_HYB_GGA_XC_B3LYP"] =
"B3LYP";
989 votca_to_orca[
"XC_HYB_GGA_XC_PBEH"] =
"PBE0";
990 votca_to_orca[
"XC_GGA_C_PBE"] =
"PBE";
991 votca_to_orca[
"XC_GGA_X_PBE"] =
"PBE";
993 std::string votca_functional =
994 options_.get(
"functional").as<std::vector<std::string>>()[0];
996 std::string orca_functional;
997 if (votca_to_orca.count(votca_functional)) {
998 orca_functional = votca_to_orca.at(votca_functional);
999 }
else if (
options_.exists(
"orca." + votca_functional)) {
1001 options_.get(
"orca." + votca_functional).as<std::string>();
1003 throw std::runtime_error(
1004 "Cannot translate " + votca_functional +
1005 " to orca functional names. Please add a <" + votca_functional +
1006 "> to your orca input options with the functional orca should use.");
1008 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
void setScaHFX(double ScaHFX)
double getDFTTotalEnergy() const
const tools::EigenSystem & MOs_beta() const
void setNumberOfAlphaElectrons(Index electrons)
void setNumberOfBetaElectrons(Index electrons)
void setECPName(const std::string &ECP)
void setXCGrid(std::string grid)
void setNumberOfOccupiedLevels(Index occupied_levels)
Index getBasisSetSize() const
void setQMEnergy(double qmenergy)
const tools::EigenSystem & MOs() const
void setQMpackage(const std::string &qmpackage)
const QMMolecule & QMAtoms() const
void setNumberOfOccupiedLevelsBeta(Index occupied_levels_beta)
void setChargeAndSpin(Index charge, Index spin)
const std::string & getECPName() const
const std::string & getDFTbasisName() const
void SetupDftBasis(std::string basis_name)
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.