votca 2025.1-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 <sys/wait.h>
47#include <unistd.h>
48#include <stdexcept>
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
352
353
354
355int run_command_spawn(const std::string& command) {
356 pid_t pid;
357
358 // system() implicitly does: /bin/sh -c "command"
359 std::vector<char*> argv = {
360 const_cast<char*>("/bin/sh"),
361 const_cast<char*>("-c"),
362 const_cast<char*>(command.c_str()),
363 nullptr
364 };
365
366 int rc = posix_spawn(
367 &pid,
368 "/bin/sh",
369 nullptr,
370 nullptr,
371 argv.data(),
372 environ
373 );
374
375 if (rc != 0) {
376 throw std::runtime_error("posix_spawn failed: " + std::to_string(rc));
377 }
378
379 int status;
380 waitpid(pid, &status, 0);
381
382 if (WIFEXITED(status))
383 return WEXITSTATUS(status);
384
385 return -1; // abnormal termination
386}
387
388
389
394
395 XTP_LOG(Log::error, *pLog_) << "Running Orca job\n" << flush;
396
397 if (std::system(nullptr)) {
398
399 std::string command = "cd " + run_dir_ + "; sh " + shell_file_name_;
400 Index check = run_command_spawn(command);
401 if (check == -1) {
403 << input_file_name_ << " failed to start" << flush;
404 return false;
405 }
406 if (CheckLogFile()) {
407 XTP_LOG(Log::error, *pLog_) << "Finished Orca job" << flush;
408 return true;
409 } else {
410 XTP_LOG(Log::error, *pLog_) << "Orca job failed" << flush;
411 }
412 } else {
414 << input_file_name_ << " failed to start" << flush;
415 return false;
416 }
417
418 return true;
419}
420
425
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());
430 }
431 // cleaning up the generated files
432 if (!cleanup_.empty()) {
433 tools::Tokenizer tok_cleanup(cleanup_, ",");
434 for (const std::string& substring : tok_cleanup) {
435 if (substring == "inp") {
436 std::string file_name = run_dir_ + "/" + input_file_name_;
437 remove(file_name.c_str());
438 }
439
440 if (substring == "bas") {
441 std::string file_name = run_dir_ + "/system.bas";
442 remove(file_name.c_str());
443 }
444
445 if (substring == "log") {
446 std::string file_name = run_dir_ + "/" + log_file_name_;
447 remove(file_name.c_str());
448 }
449
450 if (substring == "gbw") {
451 std::string file_name = run_dir_ + "/" + mo_file_name_;
452 remove(file_name.c_str());
453 }
454
455 if (substring == "ges") {
456 std::string file_name = run_dir_ + "/system.ges";
457 remove(file_name.c_str());
458 }
459 if (substring == "prop") {
460 std::string file_name = run_dir_ + "/system.prop";
461 remove(file_name.c_str());
462 }
463 }
464 }
465 return;
466}
467
469
470 StaticSegment result("charges", 0);
471
472 XTP_LOG(Log::error, *pLog_) << "Parsing " << log_file_name_ << flush;
473 std::string log_file_name_full = run_dir_ + "/" + log_file_name_;
474 std::string line;
475
476 std::ifstream input_file(log_file_name_full);
477 while (input_file) {
478 tools::getline(input_file, line);
479 boost::trim(line);
480 GetCoordinates(result, line, input_file);
481
482 std::string::size_type charge_pos = line.find("CHELPG Charges");
483
484 if (charge_pos != std::string::npos) {
485 XTP_LOG(Log::error, *pLog_) << "Getting charges" << flush;
486 tools::getline(input_file, line);
487 std::vector<std::string> row = GetLineAndSplit(input_file, "\t ");
488 Index nfields = Index(row.size());
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));
494 row = GetLineAndSplit(input_file, "\t ");
495 nfields = Index(row.size());
496 if (hasAtoms) {
497 StaticSite& temp = result.at(atom_id);
498 if (temp.getElement() != atom_type) {
499 throw std::runtime_error(
500 "Getting charges failed. Mismatch in elemts:" +
501 temp.getElement() + " vs " + atom_type);
502 }
503 temp.setCharge(atom_charge);
504 } else {
505 StaticSite temp =
506 StaticSite(atom_id, atom_type, Eigen::Vector3d::Zero());
507 temp.setCharge(atom_charge);
508 result.push_back(temp);
509 }
510 }
511 }
512 }
513 return result;
514}
515
516Eigen::Matrix3d Orca::GetPolarizability() const {
517 std::string line;
518 ifstream input_file((run_dir_ + "/" + log_file_name_));
519 bool has_pol = false;
520
521 Eigen::Matrix3d pol = Eigen::Matrix3d::Zero();
522 while (input_file) {
523 tools::getline(input_file, line);
524 boost::trim(line);
525
526 std::string::size_type pol_pos = line.find("POLARIZABILITY TENSOR");
527
528 if (pol_pos != std::string::npos) {
529
530 XTP_LOG(Log::error, *pLog_) << "Getting polarizability" << flush;
531 for (Index i = 0; i < 10; i++){
532 tools::getline(input_file, line);
533 if (line.find("The raw cartesian tensor (atomic units)") !=
534 std::string::npos) {
535 break;
536 }
537 if (i == 9) {
538 throw std::runtime_error(
539 "Could not find cartesian polarization tensor");
540 }
541 }
542
543 for (Index i = 0; i < 3; i++) {
544 tools::getline(input_file, line);
545 std::vector<double> values =
546 tools::Tokenizer(line, " ").ToVector<double>();
547 if (values.size() != 3) {
548 throw std::runtime_error("polarization line " + line +
549 " cannot be parsed");
550 }
551 Eigen::Vector3d row;
552 row << values[0], values[1], values[2];
553 pol.row(i) = row;
554 }
555
556 has_pol = true;
557 }
558 }
559 if (!has_pol) {
560 throw std::runtime_error("Could not find polarization in logfile");
561 }
562 return pol;
563}
564
566 bool found_success = false;
567 orbitals.setQMpackage(getPackageName());
568
569 orbitals.setXCGrid("xfine"); // TODO find a better approximation for the
570 // orca grid.
571 orbitals.setXCFunctionalName(options_.get("functional").as<std::string>());
572
573 XTP_LOG(Log::error, *pLog_) << "Parsing " << log_file_name_ << flush;
574 std::string log_file_name_full = run_dir_ + "/" + log_file_name_;
575 // check if LOG file is complete
576 if (!CheckLogFile()) {
577 return false;
578 }
579 std::map<Index, double> energies;
580 std::map<Index, double> energies_beta;
581
582 std::map<Index, double> occupancy;
583 std::map<Index, double> occupancy_beta;
584
585 std::string line;
586 std::string orca_version;
587 Index orca_major_version = 0;
588
589 Index levels = 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;
594
595 std::vector<std::string> results;
596
597 std::ifstream input_file(log_file_name_full);
598
599 if (input_file.fail()) {
601 << "File " << log_file_name_full << " not found " << flush;
602 return false;
603 } else {
605 << "Reading basic ORCA output from " << log_file_name_full << flush;
606 }
607 // Coordinates of the final configuration depending on whether it is an
608 // optimization or not
609
610 QMMolecule& mol = orbitals.QMAtoms();
611 orbitals.setChargeAndSpin(charge_, spin_);
612 while (input_file) {
613 tools::getline(input_file, line);
614 boost::trim(line);
615
616 GetCoordinates(mol, line, input_file);
617
618 std::string::size_type version_pos = line.find("Program Version");
619 if (version_pos != std::string::npos) {
620 results = tools::Tokenizer(line, " ").ToVector();
621 orca_version = results[2];
622 boost::trim(orca_version);
623 XTP_LOG(Log::error, *pLog_) << "ORCA Version " << orca_version << flush;
624 results = tools::Tokenizer(orca_version, ".").ToVector();
625 orca_major_version = boost::lexical_cast<Index>(results[0]);
627 << "ORCA Major Version " << orca_major_version << flush;
628 }
629
630 std::string::size_type energy_pos = line.find("FINAL SINGLE");
631 if (energy_pos != std::string::npos) {
632 results = tools::Tokenizer(line, " ").ToVector();
633 std::string energy = results[4];
634 boost::trim(energy);
635 orbitals.setQMEnergy(boost::lexical_cast<double>(energy));
636 XTP_LOG(Log::error, *pLog_) << (boost::format("QM energy[Hrt]: %4.6f ") %
637 orbitals.getDFTTotalEnergy())
638 .str()
639 << flush;
640 }
641
642 std::string::size_type HFX_pos = line.find("Fraction HF Exchange ScalHFX");
643
644 if (HFX_pos != std::string::npos) {
645 results = tools::Tokenizer(line, " ").ToVector();
646 double ScaHFX = boost::lexical_cast<double>(results.back());
647
648 orbitals.setScaHFX(ScaHFX);
650 << "DFT with " << ScaHFX << " of HF exchange!" << flush;
651 }
652
653 std::string::size_type dim_pos = line.find("Basis Dimension");
654 if (dim_pos != std::string::npos) {
655 results = tools::Tokenizer(line, " ").ToVector();
656 std::string dim =
657 results[4]; // The 4th element of results vector is the Basis Dim
658 boost::trim(dim);
659 levels = boost::lexical_cast<Index>(dim);
660 XTP_LOG(Log::info, *pLog_) << "Basis Dimension: " << levels << flush;
661 XTP_LOG(Log::info, *pLog_) << "Energy levels: " << levels << flush;
662 }
663
664 std::string::size_type OE_pos = line.find("ORBITAL ENERGIES");
665 if (OE_pos != std::string::npos) {
666
667 number_of_electrons = 0;
668 tools::getline(input_file, line);
669 tools::getline(input_file, line); // for open shell systems, this line
670 // will have "SPIN UP ORBITALS"
671 if (orbitals.isOpenShell() &&
672 line.find("SPIN UP ORBITALS") == std::string::npos) {
673 throw runtime_error(
674 "Expected to read an open-shell system but found no spin orbitals");
675 }
676 tools::getline(input_file, line);
677 if (line.find("E(Eh)") == std::string::npos) {
679 << "Warning: Orbital Energies not found in log file" << flush;
680 }
681 for (Index i = 0; i < levels; i++) {
682 if (number_of_virtuals == 10 && orca_major_version > 5) {
683 break;
684 }
685 results = GetLineAndSplit(input_file, " ");
686 std::string no = results[0];
687 boost::trim(no);
688 Index levelnumber = boost::lexical_cast<Index>(no);
689 if (levelnumber != i) {
690 XTP_LOG(Log::error, *pLog_) << "Have a look at the orbital energies "
691 "something weird is going on"
692 << flush;
693 }
694 std::string oc = results[1];
695 boost::trim(oc);
696 double occ = boost::lexical_cast<double>(oc);
697 // We only count alpha electrons, each orbital must be empty or doubly
698 // occupied
699 if (occ == 2 || occ == 1) {
700 number_of_electrons++;
701 occupancy[i] = occ;
702 } else if (occ == 0) {
703 number_of_virtuals++;
704 occupancy[i] = occ;
705 }
706 std::string e = results[2];
707 boost::trim(e);
708 energies[i] = boost::lexical_cast<double>(e);
709 }
710
711 // now read spin down energies, if needed
712 if (orbitals.isOpenShell()) {
713 number_of_electrons_beta = 0;
714 number_of_virtuals_beta = 0;
715 tools::getline(input_file, line);
716 tools::getline(input_file, line);
717 tools::getline(input_file, line);
718 for (Index i = 0; i < levels; i++) {
719 if (number_of_virtuals == 10 && orca_major_version > 5) {
720 break;
721 }
722 results = GetLineAndSplit(input_file, " ");
723 std::string no = results[0];
724 boost::trim(no);
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"
730 << flush;
731 }
732 std::string oc = results[1];
733 boost::trim(oc);
734 double occ = boost::lexical_cast<double>(oc);
735 // These occupations can only be 1 or 0
736 if (occ == 1) {
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;
742 } else {
743 throw runtime_error(
744 "Encountered spin down orbital with occupancy != 0 or 1");
745 }
746 std::string e = results[2];
747 boost::trim(e);
748 energies_beta[i] = boost::lexical_cast<double>(e);
749 }
750 }
751 }
752
753 std::string::size_type success =
754 line.find("* SUCCESS *");
755 if (success != std::string::npos) {
756 found_success = true;
757 }
758 }
759
761 if (options_.exists("ecp")) {
762 orbitals.setECPName(options_.get("ecp").as<std::string>());
763 }
764
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;
773
774 /************************************************************/
775
776 // copying information to the orbitals object
777 orbitals.setNumberOfAlphaElectrons(number_of_electrons);
778 orbitals.setNumberOfOccupiedLevels(occupied_levels);
779
780 // copying energies to a vector
781 orbitals.MOs().eigenvalues().resize(levels);
782 // level_ = 1;
783 for (Index i = 0; i < number_of_electrons + number_of_virtuals; i++) {
784 orbitals.MOs().eigenvalues()[i] = energies[i];
785 }
786
787 if (orbitals.isOpenShell()) {
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
794 << flush;
795
796 orbitals.setNumberOfBetaElectrons(number_of_electrons_beta);
797 orbitals.setNumberOfOccupiedLevelsBeta(number_of_electrons_beta);
798 orbitals.MOs_beta().eigenvalues().resize(levels);
799 for (Index i = 0; i < number_of_electrons_beta + number_of_virtuals_beta;
800 i++) {
801 orbitals.MOs_beta().eigenvalues()[i] = energies_beta[i];
802 }
803 }
804
805 XTP_LOG(Log::error, *pLog_) << "Done reading Log file" << flush;
806
807 return found_success;
808}
809template <class T>
810void Orca::GetCoordinates(T& mol, string& line, ifstream& input_file) const {
811 std::string::size_type coordinates_pos =
812 line.find("CARTESIAN COORDINATES (ANGSTROEM)");
813
814 using Atom = typename T::Atom_Type;
815
816 if (coordinates_pos != std::string::npos) {
817 XTP_LOG(Log::error, *pLog_) << "Getting the coordinates" << flush;
818 bool has_QMAtoms = mol.size() > 0;
819 // three garbage lines
820 tools::getline(input_file, line);
821 // now starts the data in format
822 // id_ type Qnuc x y z
823 vector<string> row = GetLineAndSplit(input_file, "\t ");
824 Index nfields = Index(row.size());
825 Index atom_id = 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));
831 row = GetLineAndSplit(input_file, "\t ");
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));
837 } else {
838 Atom& pAtom = mol.at(atom_id);
839 pAtom.setPos(pos);
840 }
841 atom_id++;
842 }
843 }
844}
845
847 // check if the log file exists
848 ifstream input_file(run_dir_ + "/" + log_file_name_);
849
850 if (input_file.fail()) {
851 XTP_LOG(Log::error, *pLog_) << "Orca LOG is not found" << flush;
852 return false;
853 };
854
855 std::string line;
856 while (input_file) {
857 tools::getline(input_file, line);
858 boost::trim(line);
859 std::string::size_type error = line.find("FATAL ERROR ENCOUNTERED");
860 if (error != std::string::npos) {
861 XTP_LOG(Log::error, *pLog_) << "ORCA encountered a fatal error, maybe a "
862 "look in the log file may help."
863 << flush;
864 return false;
865 }
866 error = line.find(
867 "mpirun detected that one or more processes exited with non-zero "
868 "status");
869 if (error != std::string::npos) {
871 << "ORCA had an mpi problem, maybe your openmpi version is not good."
872 << flush;
873 return false;
874 }
875 }
876 return true;
877}
878
879// Parses the molden file from orca and stores data in the Orbitals object
881 if (!CheckLogFile()) {
882 return false;
883 }
884 std::vector<double> coefficients;
885 Index basis_size = orbitals.getBasisSetSize();
886 if (basis_size == 0) {
887 throw runtime_error(
888 "Basis size not set, calculator does not parse log file first");
889 }
890
891 XTP_LOG(Log::error, *pLog_) << "Reading Molden file" << flush;
892
893 Molden molden(*pLog_);
894
895 if (orbitals.getDFTbasisName() == "") {
896 throw runtime_error(
897 "Basisset names should be set before reading the molden file.");
898 }
900
901 std::string file_name = run_dir_ + "/" +
902 mo_file_name_.substr(0, mo_file_name_.size() - 4) +
903 ".molden.input";
904 XTP_LOG(Log::error, *pLog_) << "Molden file: " << file_name << flush;
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");
914 }
915 molden.parseMoldenFile(file_name, orbitals);
916
917 XTP_LOG(Log::error, *pLog_) << "Done parsing" << flush;
918
919 // ECP charge correction is only applied in Fill() of ECPAOBasis
920 if (orbitals.getECPName() != "") {
921 ECPBasisSet ecpbasisset;
922 ecpbasisset.Load(orbitals.getECPName());
923 ECPAOBasis ecp;
924 ecp.Fill(ecpbasisset, orbitals.QMAtoms());
925 }
926
927 return true;
928}
929
930std::string Orca::indent(const double& number) {
931 std::stringstream ssnumber;
932 if (number >= 0) {
933 ssnumber << " ";
934 } else {
935 ssnumber << " ";
936 }
937 ssnumber << setiosflags(ios::fixed) << setprecision(15) << std::scientific
938 << number;
939 std::string snumber = ssnumber.str();
940 return snumber;
941}
942
943std::string Orca::CreateInputSection(const std::string& key) const {
944 std::stringstream stream;
945 std::string section = key.substr(key.find(".") + 1);
946 stream << "%" << section;
947 if (KeywordIsSingleLine(key)) {
948 stream << " " << options_.get(key).as<std::string>() << "\n";
949 } else {
950 stream << "\n"
951 << options_.get(key).as<std::string>() << "\n"
952 << "end\n";
953 }
954
955 return stream.str();
956}
957
958bool Orca::KeywordIsSingleLine(const std::string& key) const {
959 tools::Tokenizer values(this->options_.get(key).as<std::string>(), " ");
960 std::vector<std::string> words = values.ToVector();
961 return ((words.size() <= 1) ? true : false);
962}
963
964std::string Orca::WriteMethod() const {
965 std::stringstream stream;
966 std::string opt = (options_.get("optimize").as<bool>()) ? "Opt" : "";
967 const tools::Property& orca = options_.get("orca");
968 std::string user_method =
969 (orca.exists("method")) ? orca.get("method").as<std::string>() : "";
970 std::string convergence = "";
971 if (!orca.exists("scf")) {
972 convergence = this->convergence_map_.at(
973 options_.get("convergence_tightness").as<std::string>());
974 }
975 stream << "! DFT " << this->GetOrcaFunctionalName() << " " << convergence
976 << " " << opt
977 << " "
978 // additional properties provided by the user
979 << user_method << "\n";
980 return stream.str();
981}
982
983std::string Orca::GetOrcaFunctionalName() const {
984
985 std::map<std::string, std::string> votca_to_orca;
986
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";
992
993 std::string votca_functional =
994 options_.get("functional").as<std::vector<std::string>>()[0];
995
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)) {
1000 orca_functional =
1001 options_.get("orca." + votca_functional).as<std::string>();
1002 } else {
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.");
1007 }
1008 return orca_functional;
1009}
1010
1011} // namespace xtp
1012} // 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
Definition orbitals.h:46
void setScaHFX(double ScaHFX)
Definition orbitals.h:294
double getDFTTotalEnergy() const
Definition orbitals.h:193
const tools::EigenSystem & MOs_beta() const
Definition orbitals.h:128
void setNumberOfAlphaElectrons(Index electrons)
Definition orbitals.h:96
void setNumberOfBetaElectrons(Index electrons)
Definition orbitals.h:99
void setECPName(const std::string &ECP)
Definition orbitals.h:106
void setXCGrid(std::string grid)
Definition orbitals.h:187
void setNumberOfOccupiedLevels(Index occupied_levels)
Definition orbitals.h:78
Index getBasisSetSize() const
Definition orbitals.h:64
void setQMEnergy(double qmenergy)
Definition orbitals.h:195
const tools::EigenSystem & MOs() const
Definition orbitals.h:122
void setQMpackage(const std::string &qmpackage)
Definition orbitals.h:114
const QMMolecule & QMAtoms() const
Definition orbitals.h:172
void setNumberOfOccupiedLevelsBeta(Index occupied_levels_beta)
Definition orbitals.h:81
void setChargeAndSpin(Index charge, Index spin)
Definition orbitals.h:161
const std::string & getECPName() const
Definition orbitals.h:104
const std::string & getDFTbasisName() const
Definition orbitals.h:202
void SetupDftBasis(std::string basis_name)
Definition orbitals.cc:90
bool isOpenShell() const
Definition orbitals.h:168
void setXCFunctionalName(std::string functionalname)
Definition orbitals.h:182
void GetCoordinates(T &mol, std::string &line, std::ifstream &input_file) const
Definition orca.cc:810
void WriteCoordinates(std::ofstream &inp_file, const QMMolecule &)
Definition orca.cc:112
void CleanUp() override
Definition orca.cc:424
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:983
std::string indent(const double &number)
Definition orca.cc:930
std::string getPackageName() const override
Definition orca.h:40
Eigen::Matrix3d GetPolarizability() const override
Definition orca.cc:516
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:943
StaticSegment GetCharges() const override
Definition orca.cc:468
std::string WriteMethod() const
Definition orca.cc:964
void WriteBackgroundCharges()
Definition orca.cc:190
bool CheckLogFile()
Definition orca.cc:846
bool RunDFT() override
Definition orca.cc:393
bool ParseLogFile(Orbitals &orbitals) override
Definition orca.cc:565
void ParseSpecificOptions(const tools::Property &options) final
Definition orca.cc:59
bool ParseMOsFile(Orbitals &orbitals) override
Definition orca.cc:880
bool KeywordIsSingleLine(const std::string &key) const
Definition orca.cc:958
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:355
Provides a means for comparing floating point numbers.
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26
char ** environ