votca 2026-dev
Loading...
Searching...
No Matches
orca.cc
Go to the documentation of this file.
1
2/*
3 * Copyright 2009-2024 The VOTCA Development Team
4 * (http://www.votca.org)
5 *
6 * Licensed under the Apache License, Version 2.0 (the "License")
7 *
8 * You may not use this file except in compliance with the License.
9 * You may obtain a copy of the License at
10 *
11 * http://www.apache.org/licenses/LICENSE-2.0
12 *
13 * Unless required by applicable law or agreed to in writing, software
14 * distributed under the License is distributed on an "AS IS" BASIS,
15 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
16 * See the License for the specific language governing permissions and
17 * limitations under the License.
18 *
19 */
20
21// Standard includes
22#include <cstdio>
23#include <filesystem>
24#include <iomanip>
25
26// Third party includes
27#include <boost/algorithm/string.hpp>
28#include <boost/format.hpp>
29
30// VOTCA includes
31#include <stdexcept>
32#include <string>
34#include <votca/tools/getline.h>
35
36// Local VOTCA includes
37#include "votca/tools/globals.h"
39#include "votca/xtp/basisset.h"
41#include "votca/xtp/molden.h"
42#include "votca/xtp/orbitals.h"
43
44// trying posix_spwan
45#include <spawn.h>
46#include <stdexcept>
47#include <sys/wait.h>
48#include <unistd.h>
49
50extern char** environ;
51
52// Local private VOTCA includes
53#include "orca.h"
54
55namespace votca {
56namespace xtp {
57using namespace std;
58
60
61 // Orca file names
62 std::string fileName = options.get("temporary_file").as<std::string>();
63
64 input_file_name_ = fileName + ".inp";
65 log_file_name_ = fileName + ".log";
66 shell_file_name_ = fileName + ".sh";
67 mo_file_name_ = fileName + ".gbw";
68}
69
70/* Custom basis sets are written on a per-element basis to
71 * the system.bas/aux file(s), which are then included in the
72 * Orca input file using GTOName = "system.bas/aux"
73 */
74void Orca::WriteBasisset(const QMMolecule& qmatoms, std::string& bs_name,
75 std::string& el_file_name) {
76
77 std::vector<std::string> UniqueElements = qmatoms.FindUniqueElements();
78
79 tools::Elements elementInfo;
80 BasisSet bs;
81 bs.Load(bs_name);
82 XTP_LOG(Log::error, *pLog_) << "Loaded Basis Set " << bs_name << flush;
83 ofstream el_file;
84
85 el_file.open(el_file_name);
86 el_file << "$DATA" << endl;
87
88 for (const std::string& element_name : UniqueElements) {
89 const Element& element = bs.getElement(element_name);
90 el_file << elementInfo.getEleFull(element_name) << endl;
91 for (const Shell& shell : element) {
92 el_file << xtp::EnumToString(shell.getL()) << " " << shell.getSize()
93 << endl;
94 Index sh_idx = 0;
95 for (const GaussianPrimitive& gaussian : shell) {
96 sh_idx++;
97 el_file << " " << sh_idx << " " << indent(gaussian.decay());
98 el_file << " " << indent(gaussian.contraction());
99 el_file << endl;
100 }
101 }
102 }
103 el_file << "STOP\n";
104 el_file.close();
105
106 return;
107}
108
109/* Coordinates are written in standard Element,x,y,z format to the
110 * input file.
111 */
112void Orca::WriteCoordinates(std::ofstream& inp_file,
113 const QMMolecule& qmatoms) {
114
115 for (const QMAtom& atom : qmatoms) {
116 Eigen::Vector3d pos = atom.getPos() * tools::conv::bohr2ang;
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;
122 }
123 inp_file << "* \n" << endl;
124 return;
125}
126
127/* If custom ECPs are used, they need to be specified in the input file
128 * in a section following the basis set includes.
129 */
130void Orca::WriteECP(std::ofstream& inp_file, const QMMolecule& qmatoms) {
131
132 inp_file << endl;
133 std::vector<std::string> UniqueElements = qmatoms.FindUniqueElements();
134
135 ECPBasisSet ecp;
136 ecp.Load(options_.get("ecp").as<std::string>());
137
138 XTP_LOG(Log::error, *pLog_) << "Loaded Pseudopotentials "
139 << options_.get("ecp").as<std::string>() << flush;
140
141 for (const std::string& element_name : UniqueElements) {
142 try {
143 ecp.getElement(element_name);
144 } catch (std::runtime_error& error) {
146 << "No pseudopotential for " << element_name << " available" << flush;
147 continue;
148 }
149 const ECPElement& element = ecp.getElement(element_name);
150
151 inp_file << "\n"
152 << "NewECP"
153 << " " << element_name << endl;
154 inp_file << "N_core"
155 << " " << element.getNcore() << endl;
156 inp_file << "lmax"
157 << " " << EnumToString(element.getLmax()) << endl;
158 // For Orca the order doesn't matter but let's write it in ascending order
159 // write remaining shells in ascending order s,p,d...
160 for (Index i = 0; i <= Index(element.getLmax()); i++) {
161 for (const ECPShell& shell : element) {
162 if (Index(shell.getL()) == i) {
163 // shell type, number primitives, scale factor
164 inp_file << xtp::EnumToString(shell.getL()) << " " << shell.getSize()
165 << endl;
166 Index sh_idx = 0;
167 for (const ECPGaussianPrimitive& gaussian : shell) {
168 sh_idx++;
169 inp_file << sh_idx << " " << gaussian.decay_ << " "
170 << gaussian.contraction_ << " " << gaussian.power_ << endl;
171 }
172 }
173 }
174 }
175 inp_file << "end\n "
176 << "\n"
177 << endl;
178 }
179 return;
180}
181
183 this->options_.addTree("orca.pointcharges", "\"background.crg\"");
184}
185
186/* For QM/MM the molecules in the MM environment are represented by
187 * their atomic partial charge distributions. ORCA expects them in
188 * q,x,y,z format in a separate file "background.crg"
189 */
191
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;
196
197 for (const std::unique_ptr<StaticSite>& site : externalsites_) {
198 if (site->getCharge() != 0.0) {
199 total_background++;
200 }
201 total_background += SplitMultipoles(*site).size();
202 } // counting only
203
204 crg_file << total_background << endl;
205 boost::format fmt("%1$+1.7f %2$+1.7f %3$+1.7f %4$+1.7f");
206 // now write
207 for (const std::unique_ptr<StaticSite>& site : externalsites_) {
208 Eigen::Vector3d pos = site->getPos() * tools::conv::bohr2ang;
209 string sitestring =
210 boost::str(fmt % site->getCharge() % pos.x() % pos.y() % pos.z());
211 if (site->getCharge() != 0.0) {
212 crg_file << sitestring << endl;
213 }
214 std::vector<MinimalMMCharge> split_multipoles = SplitMultipoles(*site);
215 for (const auto& mpoles : split_multipoles) {
216 Eigen::Vector3d pos2 = mpoles.pos_ * tools::conv::bohr2ang;
217 string multipole =
218 boost::str(fmt % mpoles.q_ % pos2.x() % pos2.y() % pos2.z());
219 crg_file << multipole << endl;
220 }
221 }
222
223 return;
224}
225
231bool Orca::WriteInputFile(const Orbitals& orbitals) {
232
233 std::vector<std::string> results;
234 std::string temp_suffix = "/id";
235 std::string scratch_dir_backup = scratch_dir_;
236 std::ofstream inp_file;
237 std::string inp_file_name_full = run_dir_ + "/" + input_file_name_;
238 inp_file.open(inp_file_name_full);
239 // header
240 inp_file << "* xyz " << charge_ << " " << spin_ << endl;
241 Index threads = OPENMP::getMaxThreads();
242 const QMMolecule& qmatoms = orbitals.QMAtoms();
243 // put coordinates
244 WriteCoordinates(inp_file, qmatoms);
245 // add parallelization info
246 inp_file << "%pal\n"
247 << "nprocs " << threads << "\nend"
248 << "\n"
249 << endl;
250 // basis set info
251 std::string el_file_name = run_dir_ + "/" + "system.bas";
252 WriteBasisset(qmatoms, basisset_name_, el_file_name);
253 inp_file << "%basis\n";
254 inp_file << "GTOName"
255 << " "
256 << "="
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>();
262 WriteBasisset(qmatoms, auxbasisset_name, aux_file_name);
263 inp_file << "GTOAuxName"
264 << " "
265 << "="
266 << "\"system.aux\";" << endl;
267 } // write_auxbasis set
268
269 // ECPs
270 if (options_.exists("ecp")) {
271 WriteECP(inp_file, qmatoms);
272 }
273 inp_file << "end\n "
274 << "\n"
275 << endl; // This end is for the basis set block
276 if (!externalsites_.empty()) {
278 }
279
280 // External Electric field
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";
286 inp_file << "end\n";
287 inp_file << std::endl;
288 }
289 std::string input_options;
290 // Write Orca section specified by the user
291 for (const auto& prop : this->options_.get("orca")) {
292 const std::string& prop_name = prop.name();
293 if (prop_name == "pointcharges") {
294 input_options += this->CreateInputSection("orca.pointcharges");
295 } else if (prop_name != "method") {
296 input_options += this->CreateInputSection("orca." + prop_name);
297 }
298 }
299 // Write main DFT method
300 input_options += this->WriteMethod();
301 inp_file << input_options;
302 inp_file << '\n';
303 inp_file.close();
304 // and now generate a shell script to run both jobs, if neccessary
305
307 << "Setting the scratch dir to " << scratch_dir_ + temp_suffix << flush;
308 scratch_dir_ = scratch_dir_backup + temp_suffix;
310 scratch_dir_ = scratch_dir_backup;
311 return true;
312}
313
315 ofstream shell_file;
316 std::string shell_file_name_full = run_dir_ + "/" + shell_file_name_;
317 shell_file.open(shell_file_name_full);
318 shell_file << "#!/usr/bin/env bash" << endl;
319 if (run_dir_ != ".") {
320 const std::string key = "molecule";
321 std::string path = run_dir_;
322 auto pos = path.find(key);
323
324 if (pos != std::string::npos) {
325 scratch_dir_ += "/" + path.erase(0, pos);
326 }
327 }
328 shell_file << "mkdir -p " << scratch_dir_ << endl;
329 shell_file << "export TMPDIR=" << scratch_dir_ << endl;
330 std::string base_name = mo_file_name_.substr(0, mo_file_name_.size() - 4);
331
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"))) {
335 throw runtime_error(
336 "Using guess relies on a molA.gbw and a molB.gbw file being in the "
337 "directory.");
338 }
339 shell_file << options_.get("executable").as<std::string>()
340 << "_mergefrag molA.gbw molB.gbw " << base_name
341 << ".gbw > merge.log" << endl;
342 }
343 shell_file << options_.get("executable").as<std::string>() << " "
344 << input_file_name_ << " > " << log_file_name_ << endl;
345
346 shell_file << options_.get("executable").as<std::string>() << "_2mkl "
347 << base_name << " -molden > molden.log" << endl;
348 shell_file.close();
349 return true;
350}
351
352int run_command_spawn(const std::string& command) {
353 pid_t pid;
354
355 // system() implicitly does: /bin/sh -c "command"
356 std::vector<char*> argv = {const_cast<char*>("/bin/sh"),
357 const_cast<char*>("-c"),
358 const_cast<char*>(command.c_str()), nullptr};
359
360 int rc = posix_spawn(&pid, "/bin/sh", nullptr, nullptr, argv.data(), environ);
361
362 if (rc != 0) {
363 throw std::runtime_error("posix_spawn failed: " + std::to_string(rc));
364 }
365
366 int status;
367 waitpid(pid, &status, 0);
368
369 if (WIFEXITED(status)) return WEXITSTATUS(status);
370
371 return -1; // abnormal termination
372}
373
378
379 XTP_LOG(Log::error, *pLog_) << "Running Orca job\n" << flush;
380
381 if (std::system(nullptr)) {
382
383 std::string command = "cd " + run_dir_ + "; sh " + shell_file_name_;
384 Index check = run_command_spawn(command);
385 if (check == -1) {
387 << input_file_name_ << " failed to start" << flush;
388 return false;
389 }
390 if (CheckLogFile()) {
391 XTP_LOG(Log::error, *pLog_) << "Finished Orca job" << flush;
392 return true;
393 } else {
394 XTP_LOG(Log::error, *pLog_) << "Orca job failed" << flush;
395 }
396 } else {
398 << input_file_name_ << " failed to start" << flush;
399 return false;
400 }
401
402 return true;
403}
404
409
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());
414 }
415 // cleaning up the generated files
416 if (!cleanup_.empty()) {
417 tools::Tokenizer tok_cleanup(cleanup_, ",");
418 for (const std::string& substring : tok_cleanup) {
419 if (substring == "inp") {
420 std::string file_name = run_dir_ + "/" + input_file_name_;
421 remove(file_name.c_str());
422 }
423
424 if (substring == "bas") {
425 std::string file_name = run_dir_ + "/system.bas";
426 remove(file_name.c_str());
427 }
428
429 if (substring == "log") {
430 std::string file_name = run_dir_ + "/" + log_file_name_;
431 remove(file_name.c_str());
432 }
433
434 if (substring == "gbw") {
435 std::string file_name = run_dir_ + "/" + mo_file_name_;
436 remove(file_name.c_str());
437 }
438
439 if (substring == "ges") {
440 std::string file_name = run_dir_ + "/system.ges";
441 remove(file_name.c_str());
442 }
443 if (substring == "prop") {
444 std::string file_name = run_dir_ + "/system.prop";
445 remove(file_name.c_str());
446 }
447 }
448 }
449 return;
450}
451
453
454 StaticSegment result("charges", 0);
455
456 XTP_LOG(Log::error, *pLog_) << "Parsing " << log_file_name_ << flush;
457 std::string log_file_name_full = run_dir_ + "/" + log_file_name_;
458 std::string line;
459
460 std::ifstream input_file(log_file_name_full);
461 while (input_file) {
462 tools::getline(input_file, line);
463 boost::trim(line);
464 GetCoordinates(result, line, input_file);
465
466 std::string::size_type charge_pos = line.find("CHELPG Charges");
467
468 if (charge_pos != std::string::npos) {
469 XTP_LOG(Log::error, *pLog_) << "Getting charges" << flush;
470 tools::getline(input_file, line);
471 std::vector<std::string> row = GetLineAndSplit(input_file, "\t ");
472 Index nfields = Index(row.size());
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));
478 row = GetLineAndSplit(input_file, "\t ");
479 nfields = Index(row.size());
480 if (hasAtoms) {
481 StaticSite& temp = result.at(atom_id);
482 if (temp.getElement() != atom_type) {
483 throw std::runtime_error(
484 "Getting charges failed. Mismatch in elemts:" +
485 temp.getElement() + " vs " + atom_type);
486 }
487 temp.setCharge(atom_charge);
488 } else {
489 StaticSite temp =
490 StaticSite(atom_id, atom_type, Eigen::Vector3d::Zero());
491 temp.setCharge(atom_charge);
492 result.push_back(temp);
493 }
494 }
495 }
496 }
497 return result;
498}
499
500Eigen::Matrix3d Orca::GetPolarizability() const {
501 std::string line;
502 ifstream input_file((run_dir_ + "/" + log_file_name_));
503 bool has_pol = false;
504
505 Eigen::Matrix3d pol = Eigen::Matrix3d::Zero();
506 while (input_file) {
507 tools::getline(input_file, line);
508 boost::trim(line);
509
510 std::string::size_type pol_pos = line.find("POLARIZABILITY TENSOR");
511
512 if (pol_pos != std::string::npos) {
513
514 XTP_LOG(Log::error, *pLog_) << "Getting polarizability" << flush;
515 for (Index i = 0; i < 10; i++) {
516 tools::getline(input_file, line);
517 if (line.find("The raw cartesian tensor (atomic units)") !=
518 std::string::npos) {
519 break;
520 }
521 if (i == 9) {
522 throw std::runtime_error(
523 "Could not find cartesian polarization tensor");
524 }
525 }
526
527 for (Index i = 0; i < 3; i++) {
528 tools::getline(input_file, line);
529 std::vector<double> values =
530 tools::Tokenizer(line, " ").ToVector<double>();
531 if (values.size() != 3) {
532 throw std::runtime_error("polarization line " + line +
533 " cannot be parsed");
534 }
535 Eigen::Vector3d row;
536 row << values[0], values[1], values[2];
537 pol.row(i) = row;
538 }
539
540 has_pol = true;
541 }
542 }
543 if (!has_pol) {
544 throw std::runtime_error("Could not find polarization in logfile");
545 }
546 return pol;
547}
548
550 bool found_success = false;
551 orbitals.setQMpackage(getPackageName());
552
553 orbitals.setXCGrid("xfine"); // TODO find a better approximation for the
554 // orca grid.
555 orbitals.setXCFunctionalName(options_.get("functional").as<std::string>());
556
557 XTP_LOG(Log::error, *pLog_) << "Parsing " << log_file_name_ << flush;
558 std::string log_file_name_full = run_dir_ + "/" + log_file_name_;
559 // check if LOG file is complete
560 if (!CheckLogFile()) {
561 return false;
562 }
563 std::map<Index, double> energies;
564 std::map<Index, double> energies_beta;
565
566 std::map<Index, double> occupancy;
567 std::map<Index, double> occupancy_beta;
568
569 std::string line;
570 std::string orca_version;
571 Index orca_major_version = 0;
572
573 Index levels = 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;
578
579 std::vector<std::string> results;
580
581 std::ifstream input_file(log_file_name_full);
582
583 if (input_file.fail()) {
585 << "File " << log_file_name_full << " not found " << flush;
586 return false;
587 } else {
589 << "Reading basic ORCA output from " << log_file_name_full << flush;
590 }
591 // Coordinates of the final configuration depending on whether it is an
592 // optimization or not
593
594 QMMolecule& mol = orbitals.QMAtoms();
595 orbitals.setChargeAndSpin(charge_, spin_);
596 while (input_file) {
597 tools::getline(input_file, line);
598 boost::trim(line);
599
600 GetCoordinates(mol, line, input_file);
601
602 std::string::size_type version_pos = line.find("Program Version");
603 if (version_pos != std::string::npos) {
604 results = tools::Tokenizer(line, " ").ToVector();
605 orca_version = results[2];
606 boost::trim(orca_version);
607 XTP_LOG(Log::error, *pLog_) << "ORCA Version " << orca_version << flush;
608 results = tools::Tokenizer(orca_version, ".").ToVector();
609 orca_major_version = boost::lexical_cast<Index>(results[0]);
611 << "ORCA Major Version " << orca_major_version << flush;
612 }
613
614 std::string::size_type energy_pos = line.find("FINAL SINGLE");
615 if (energy_pos != std::string::npos) {
616 results = tools::Tokenizer(line, " ").ToVector();
617 std::string energy = results[4];
618 boost::trim(energy);
619 orbitals.setQMEnergy(boost::lexical_cast<double>(energy));
620 XTP_LOG(Log::error, *pLog_) << (boost::format("QM energy[Hrt]: %4.6f ") %
621 orbitals.getDFTTotalEnergy())
622 .str()
623 << flush;
624 }
625
626 std::string::size_type HFX_pos = line.find("Fraction HF Exchange ScalHFX");
627
628 if (HFX_pos != std::string::npos) {
629 results = tools::Tokenizer(line, " ").ToVector();
630 double ScaHFX = boost::lexical_cast<double>(results.back());
631
632 orbitals.setScaHFX(ScaHFX);
634 << "DFT with " << ScaHFX << " of HF exchange!" << flush;
635 }
636
637 std::string::size_type dim_pos = line.find("Basis Dimension");
638 if (dim_pos != std::string::npos) {
639 results = tools::Tokenizer(line, " ").ToVector();
640 std::string dim =
641 results[4]; // The 4th element of results vector is the Basis Dim
642 boost::trim(dim);
643 levels = boost::lexical_cast<Index>(dim);
644 XTP_LOG(Log::info, *pLog_) << "Basis Dimension: " << levels << flush;
645 XTP_LOG(Log::info, *pLog_) << "Energy levels: " << levels << flush;
646 }
647
648 std::string::size_type OE_pos = line.find("ORBITAL ENERGIES");
649 if (OE_pos != std::string::npos) {
650
651 number_of_electrons = 0;
652 tools::getline(input_file, line);
653 tools::getline(input_file, line); // for open shell systems, this line
654 // will have "SPIN UP ORBITALS"
655 if (orbitals.isOpenShell() &&
656 line.find("SPIN UP ORBITALS") == std::string::npos) {
657 throw runtime_error(
658 "Expected to read an open-shell system but found no spin orbitals");
659 }
660 tools::getline(input_file, line);
661 if (line.find("E(Eh)") == std::string::npos) {
663 << "Warning: Orbital Energies not found in log file" << flush;
664 }
665 for (Index i = 0; i < levels; i++) {
666 if (number_of_virtuals == 10 && orca_major_version > 5) {
667 break;
668 }
669 results = GetLineAndSplit(input_file, " ");
670 std::string no = results[0];
671 boost::trim(no);
672 Index levelnumber = boost::lexical_cast<Index>(no);
673 if (levelnumber != i) {
674 XTP_LOG(Log::error, *pLog_) << "Have a look at the orbital energies "
675 "something weird is going on"
676 << flush;
677 }
678 std::string oc = results[1];
679 boost::trim(oc);
680 double occ = boost::lexical_cast<double>(oc);
681 // We only count alpha electrons, each orbital must be empty or doubly
682 // occupied
683 if (occ == 2 || occ == 1) {
684 number_of_electrons++;
685 occupancy[i] = occ;
686 } else if (occ == 0) {
687 number_of_virtuals++;
688 occupancy[i] = occ;
689 }
690 std::string e = results[2];
691 boost::trim(e);
692 energies[i] = boost::lexical_cast<double>(e);
693 }
694
695 // now read spin down energies, if needed
696 if (orbitals.isOpenShell()) {
697 number_of_electrons_beta = 0;
698 number_of_virtuals_beta = 0;
699 tools::getline(input_file, line);
700 tools::getline(input_file, line);
701 tools::getline(input_file, line);
702 for (Index i = 0; i < levels; i++) {
703 if (number_of_virtuals == 10 && orca_major_version > 5) {
704 break;
705 }
706 results = GetLineAndSplit(input_file, " ");
707 std::string no = results[0];
708 boost::trim(no);
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"
714 << flush;
715 }
716 std::string oc = results[1];
717 boost::trim(oc);
718 double occ = boost::lexical_cast<double>(oc);
719 // These occupations can only be 1 or 0
720 if (occ == 1) {
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;
726 } else {
727 throw runtime_error(
728 "Encountered spin down orbital with occupancy != 0 or 1");
729 }
730 std::string e = results[2];
731 boost::trim(e);
732 energies_beta[i] = boost::lexical_cast<double>(e);
733 }
734 }
735 }
736
737 std::string::size_type success =
738 line.find("* SUCCESS *");
739 if (success != std::string::npos) {
740 found_success = true;
741 }
742 }
743
745 if (options_.exists("ecp")) {
746 orbitals.setECPName(options_.get("ecp").as<std::string>());
747 }
748
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;
757
758 /************************************************************/
759
760 // copying information to the orbitals object
761 orbitals.setNumberOfAlphaElectrons(number_of_electrons);
762 orbitals.setNumberOfOccupiedLevels(occupied_levels);
763
764 // copying energies to a vector
765 orbitals.MOs().eigenvalues().resize(levels);
766 // level_ = 1;
767 for (Index i = 0; i < number_of_electrons + number_of_virtuals; i++) {
768 orbitals.MOs().eigenvalues()[i] = energies[i];
769 }
770
771 if (orbitals.isOpenShell()) {
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
778 << flush;
779
780 orbitals.setNumberOfBetaElectrons(number_of_electrons_beta);
781 orbitals.setNumberOfOccupiedLevelsBeta(number_of_electrons_beta);
782 orbitals.MOs_beta().eigenvalues().resize(levels);
783 for (Index i = 0; i < number_of_electrons_beta + number_of_virtuals_beta;
784 i++) {
785 orbitals.MOs_beta().eigenvalues()[i] = energies_beta[i];
786 }
787 }
788
789 XTP_LOG(Log::error, *pLog_) << "Done reading Log file" << flush;
790
791 return found_success;
792}
793template <class T>
794void Orca::GetCoordinates(T& mol, string& line, ifstream& input_file) const {
795 std::string::size_type coordinates_pos =
796 line.find("CARTESIAN COORDINATES (ANGSTROEM)");
797
798 using Atom = typename T::Atom_Type;
799
800 if (coordinates_pos != std::string::npos) {
801 XTP_LOG(Log::error, *pLog_) << "Getting the coordinates" << flush;
802 bool has_QMAtoms = mol.size() > 0;
803 // three garbage lines
804 tools::getline(input_file, line);
805 // now starts the data in format
806 // id_ type Qnuc x y z
807 vector<string> row = GetLineAndSplit(input_file, "\t ");
808 Index nfields = Index(row.size());
809 Index atom_id = 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));
815 row = GetLineAndSplit(input_file, "\t ");
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));
821 } else {
822 Atom& pAtom = mol.at(atom_id);
823 pAtom.setPos(pos);
824 }
825 atom_id++;
826 }
827 }
828}
829
831 // check if the log file exists
832 ifstream input_file(run_dir_ + "/" + log_file_name_);
833
834 if (input_file.fail()) {
835 XTP_LOG(Log::error, *pLog_) << "Orca LOG is not found" << flush;
836 return false;
837 };
838
839 std::string line;
840 while (input_file) {
841 tools::getline(input_file, line);
842 boost::trim(line);
843 std::string::size_type error = line.find("FATAL ERROR ENCOUNTERED");
844 if (error != std::string::npos) {
845 XTP_LOG(Log::error, *pLog_) << "ORCA encountered a fatal error, maybe a "
846 "look in the log file may help."
847 << flush;
848 return false;
849 }
850 error = line.find(
851 "mpirun detected that one or more processes exited with non-zero "
852 "status");
853 if (error != std::string::npos) {
855 << "ORCA had an mpi problem, maybe your openmpi version is not good."
856 << flush;
857 return false;
858 }
859 }
860 return true;
861}
862
863// Parses the molden file from orca and stores data in the Orbitals object
865 if (!CheckLogFile()) {
866 return false;
867 }
868 std::vector<double> coefficients;
869 Index basis_size = orbitals.getBasisSetSize();
870 if (basis_size == 0) {
871 throw runtime_error(
872 "Basis size not set, calculator does not parse log file first");
873 }
874
875 XTP_LOG(Log::error, *pLog_) << "Reading Molden file" << flush;
876
877 Molden molden(*pLog_);
878
879 if (orbitals.getDFTbasisName() == "") {
880 throw runtime_error(
881 "Basisset names should be set before reading the molden file.");
882 }
884
885 std::string file_name = run_dir_ + "/" +
886 mo_file_name_.substr(0, mo_file_name_.size() - 4) +
887 ".molden.input";
888 XTP_LOG(Log::error, *pLog_) << "Molden file: " << file_name << flush;
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");
898 }
899 molden.parseMoldenFile(file_name, orbitals);
900
901 XTP_LOG(Log::error, *pLog_) << "Done parsing" << flush;
902
903 // ECP charge correction is only applied in Fill() of ECPAOBasis
904 if (orbitals.getECPName() != "") {
905 ECPBasisSet ecpbasisset;
906 ecpbasisset.Load(orbitals.getECPName());
907 ECPAOBasis ecp;
908 ecp.Fill(ecpbasisset, orbitals.QMAtoms());
909 }
910
911 return true;
912}
913
914std::string Orca::indent(const double& number) {
915 std::stringstream ssnumber;
916 if (number >= 0) {
917 ssnumber << " ";
918 } else {
919 ssnumber << " ";
920 }
921 ssnumber << setiosflags(ios::fixed) << setprecision(15) << std::scientific
922 << number;
923 std::string snumber = ssnumber.str();
924 return snumber;
925}
926
927std::string Orca::CreateInputSection(const std::string& key) const {
928 std::stringstream stream;
929 std::string section = key.substr(key.find(".") + 1);
930 stream << "%" << section;
931 if (KeywordIsSingleLine(key)) {
932 stream << " " << options_.get(key).as<std::string>() << "\n";
933 } else {
934 stream << "\n"
935 << options_.get(key).as<std::string>() << "\n"
936 << "end\n";
937 }
938
939 return stream.str();
940}
941
942bool Orca::KeywordIsSingleLine(const std::string& key) const {
943 tools::Tokenizer values(this->options_.get(key).as<std::string>(), " ");
944 std::vector<std::string> words = values.ToVector();
945 return ((words.size() <= 1) ? true : false);
946}
947
948std::string Orca::WriteMethod() const {
949 std::stringstream stream;
950 std::string opt = (options_.get("optimize").as<bool>()) ? "Opt" : "";
951 const tools::Property& orca = options_.get("orca");
952 std::string user_method =
953 (orca.exists("method")) ? orca.get("method").as<std::string>() : "";
954 std::string convergence = "";
955 if (!orca.exists("scf")) {
956 convergence = this->convergence_map_.at(
957 options_.get("convergence_tightness").as<std::string>());
958 }
959 stream << "! DFT " << this->GetOrcaFunctionalName() << " " << convergence
960 << " " << opt
961 << " "
962 // additional properties provided by the user
963 << user_method << "\n";
964 return stream.str();
965}
966
967std::string Orca::GetOrcaFunctionalName() const {
968
969 std::map<std::string, std::string> votca_to_orca;
970
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";
976
977 std::string votca_functional =
978 options_.get("functional").as<std::vector<std::string>>()[0];
979
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)) {
984 orca_functional =
985 options_.get("orca." + votca_functional).as<std::string>();
986 } else {
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.");
991 }
992 return orca_functional;
993}
994
995} // namespace xtp
996} // namespace votca
const Eigen::VectorXd & eigenvalues() const
Definition eigensystem.h:30
information about an element
Definition elements.h:42
std::string getEleFull(std::string eleshort)
Definition elements.cc:137
class to manage program options with xml serialization functionality
Definition property.h:55
Property & get(const std::string &key)
get existing property
Definition property.cc:79
bool exists(const std::string &key) const
check whether property exists
Definition property.cc:122
T as() const
return value as type
Definition property.h:283
break string into words
Definition tokenizer.h:72
std::vector< T > ToVector()
store all words in a vector of type T, does type conversion.
Definition tokenizer.h:109
void push_back(const T &atom)
const T & at(Index index) const
std::vector< std::string > FindUniqueElements() const
void setPos(const Eigen::Vector3d &r)
Definition atom.h:67
const Element & getElement(std::string element_type) const
Definition basisset.cc:212
void Load(const std::string &name)
Definition basisset.cc:149
Container to hold ECPs for all atoms.
Definition ecpaobasis.h:43
std::vector< std::string > Fill(const ECPBasisSet &bs, QMMolecule &atoms)
Definition ecpaobasis.cc:69
void Load(const std::string &name)
const ECPElement & getElement(std::string element_type) const
void parseMoldenFile(const std::string &filename, Orbitals &orbitals) const
Definition molden.cc:245
void setBasissetInfo(const std::string &basisset_name, const std::string &aux_basisset_name="")
Definition molden.h:41
Container for molecular orbitals and derived one-particle data.
Definition orbitals.h:47
void setScaHFX(double ScaHFX)
Store the fraction of exact exchange associated with the functional.
Definition orbitals.h:419
double getDFTTotalEnergy() const
Return the stored total DFT energy.
Definition orbitals.h:289
const tools::EigenSystem & MOs_beta() const
Return read-only access to beta-spin molecular orbitals.
Definition orbitals.h:202
void setNumberOfAlphaElectrons(Index electrons)
Store the total number of alpha electrons.
Definition orbitals.h:139
void setNumberOfBetaElectrons(Index electrons)
Store the total number of beta electrons.
Definition orbitals.h:144
void setECPName(const std::string &ECP)
Store the effective core potential label.
Definition orbitals.h:170
void setXCGrid(std::string grid)
Store the numerical XC grid quality label.
Definition orbitals.h:280
void setNumberOfOccupiedLevels(Index occupied_levels)
Definition orbitals.h:122
Index getBasisSetSize() const
Return the number of AO basis functions in the DFT basis.
Definition orbitals.h:72
void setQMEnergy(double qmenergy)
Store the total DFT energy.
Definition orbitals.h:292
const tools::EigenSystem & MOs() const
Return read-only access to alpha/restricted molecular orbitals.
Definition orbitals.h:192
void setQMpackage(const std::string &qmpackage)
Store the name of the QM package that produced these orbitals.
Definition orbitals.h:181
const QMMolecule & QMAtoms() const
Return read-only access to the molecular geometry.
Definition orbitals.h:259
void setNumberOfOccupiedLevelsBeta(Index occupied_levels_beta)
Store the number of occupied beta-spin orbitals.
Definition orbitals.h:134
void setChargeAndSpin(Index charge, Index spin)
Definition orbitals.h:243
const std::string & getECPName() const
Return the effective core potential label.
Definition orbitals.h:167
const std::string & getDFTbasisName() const
Return the DFT basis-set name.
Definition orbitals.h:301
void SetupDftBasis(std::string basis_name)
Build and attach the DFT AO basis from the stored molecular geometry.
Definition orbitals.cc:99
bool isOpenShell() const
Report whether the stored state corresponds to an open-shell system.
Definition orbitals.h:253
void setXCFunctionalName(std::string functionalname)
Definition orbitals.h:273
void GetCoordinates(T &mol, std::string &line, std::ifstream &input_file) const
Definition orca.cc:794
void WriteCoordinates(std::ofstream &inp_file, const QMMolecule &)
Definition orca.cc:112
void CleanUp() override
Definition orca.cc:408
void WriteECP(std::ofstream &inp_file, const QMMolecule &)
Definition orca.cc:130
bool WriteShellScript()
Definition orca.cc:314
void WriteBasisset(const QMMolecule &qmatoms, std::string &bs_name, std::string &el_file_name)
Definition orca.cc:74
std::string GetOrcaFunctionalName() const
Definition orca.cc:967
std::string indent(const double &number)
Definition orca.cc:914
std::string getPackageName() const override
Definition orca.h:40
Eigen::Matrix3d GetPolarizability() const override
Definition orca.cc:500
void WriteChargeOption() override
Definition orca.cc:182
std::map< std::string, std::string > convergence_map_
Definition orca.h:114
std::string CreateInputSection(const std::string &key) const
Definition orca.cc:927
StaticSegment GetCharges() const override
Definition orca.cc:452
std::string WriteMethod() const
Definition orca.cc:948
void WriteBackgroundCharges()
Definition orca.cc:190
bool CheckLogFile()
Definition orca.cc:830
bool RunDFT() override
Definition orca.cc:377
bool ParseLogFile(Orbitals &orbitals) override
Definition orca.cc:549
void ParseSpecificOptions(const tools::Property &options) final
Definition orca.cc:59
bool ParseMOsFile(Orbitals &orbitals) override
Definition orca.cc:864
bool KeywordIsSingleLine(const std::string &key) const
Definition orca.cc:942
bool WriteInputFile(const Orbitals &orbitals) override
Definition orca.cc:231
container for QM atoms
Definition qmatom.h:37
std::string log_file_name_
Definition qmpackage.h:135
std::vector< std::string > GetLineAndSplit(std::ifstream &input_file, const std::string separators) const
Definition qmpackage.cc:144
std::string run_dir_
Definition qmpackage.h:137
std::vector< std::unique_ptr< StaticSite > > externalsites_
Definition qmpackage.h:144
std::vector< MinimalMMCharge > SplitMultipoles(const StaticSite &site) const
Definition qmpackage.cc:109
std::string mo_file_name_
Definition qmpackage.h:136
std::string basisset_name_
Definition qmpackage.h:132
std::string cleanup_
Definition qmpackage.h:133
tools::Property options_
Definition qmpackage.h:140
std::string shell_file_name_
Definition qmpackage.h:139
std::string scratch_dir_
Definition qmpackage.h:138
std::string input_file_name_
Definition qmpackage.h:134
Class to represent Atom/Site in electrostatic.
Definition staticsite.h:37
const std::string & getElement() const
Definition staticsite.h:79
void setCharge(double q)
Definition staticsite.h:87
#define XTP_LOG(level, log)
Definition logger.h:40
const double ang2bohr
Definition constants.h:48
const double bohr2ang
Definition constants.h:49
std::istream & getline(std::istream &is, std::string &str)
Wrapper for a getline function.
Definition getline.h:35
Index getMaxThreads()
Definition eigen.h:128
Charge transport classes.
Definition ERIs.h:28
std::string EnumToString(L l)
Definition basisset.cc:60
ClassicalSegment< StaticSite > StaticSegment
int run_command_spawn(const std::string &command)
Definition orca.cc:352
Provides a means for comparing floating point numbers.
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26
char ** environ