7 #ifndef READUCT_IRCTASK_H_
8 #define READUCT_IRCTASK_H_
21 #include <boost/exception/diagnostic_information.hpp>
22 #include <boost/filesystem.hpp>
40 IrcTask(std::vector<std::string>
input, std::vector<std::string>
output, std::shared_ptr<Core::Log> logger =
nullptr)
41 :
Task(std::move(input), std::move(output), std::move(logger)) {
44 std::string
name()
const override {
45 return "IRC Optimizations";
51 if (_output.size() == 1) {
53 <<
" Warning: To store IRC results in a new location two output systems have to be specified.\n";
56 bool silentCalculator = taskSettings.extract(
"silent_stdout_calculator",
true);
57 std::shared_ptr<Core::Calculator> calc;
60 calc = copyCalculator(systems, _input.front(),
name());
61 Utils::CalculationRoutines::setLog(*calc,
true,
true, !silentCalculator);
64 if (calc->getStructure()->size() == 1) {
65 throw std::runtime_error(
"Cannot perform IRC task for monoatomic systems.");
70 std::string optimizertype = taskSettings.extract(
"optimizer", std::string{
"SD"});
71 int mode = taskSettings.extract(
"irc_mode", 0);
73 bool stopOnError = stopOnErrorExtraction(taskSettings);
74 std::transform(optimizertype.begin(), optimizertype.end(), optimizertype.begin(), ::toupper);
75 std::shared_ptr<Utils::IrcOptimizerBase> optimizer;
76 if (optimizertype ==
"LBFGS") {
77 auto tmp = std::make_shared<Utils::IrcOptimizer<Utils::Lbfgs>>(*calc);
79 optimizer = std::move(tmp);
81 else if (optimizertype ==
"BFGS") {
82 auto tmp = std::make_shared<Utils::IrcOptimizer<Utils::Bfgs>>(*calc);
84 optimizer = std::move(tmp);
86 else if (optimizertype ==
"SD" || optimizertype ==
"STEEPESTDESCENT") {
87 auto tmp = std::make_shared<Utils::IrcOptimizer<Utils::SteepestDescent>>(*calc);
89 optimizer = std::move(tmp);
92 throw std::runtime_error(
93 "Unknown Optimizer requested for an IRC optimization, available are: SD, BFGS and LBFGS!");
96 auto settings = optimizer->getSettings();
97 settings.merge(taskSettings);
98 if (!settings.valid()) {
99 settings.throwIncorrectSettings();
101 optimizer->setSettings(settings);
102 if (!testRunOnly && !Utils::GeometryOptimization::settingsMakeSense(*optimizer)) {
103 throw std::logic_error(
"The given calculator settings are too inaccurate for the given convergence criteria of "
104 "this optimization Task");
113 if (!calc->results().has<Utils::Property::Hessian>()) {
114 calc->setRequiredProperties(Utils::Property::Hessian);
115 calc->calculate(
"Hessian Calculation");
117 auto hessian = calc->results().get<Utils::Property::Hessian>();
118 auto system = calc->getStructure();
119 auto modes = Utils::NormalModeAnalysis::calculateNormalModes(hessian, system->getElements(), system->getPositions());
120 if (mode < 0 || mode > modes.size() - 1) {
121 throw std::runtime_error(
122 "The chosen normal mode number is smaller than 0 or larger than the number of modes present!");
124 auto ircMode = modes.getMode(mode);
125 auto ircModeVector = Eigen::Map<const Eigen::VectorXd>(ircMode.data(), ircMode.cols() * ircMode.rows());
130 auto cout = _logger->output;
131 cout <<
"\n IRC Optimization: Forward\n\n";
132 std::shared_ptr<Scine::Core::Calculator> forwardCalc = std::move(calc->clone());
134 double oldEnergy = 0.0;
135 Eigen::VectorXd oldParams;
139 const std::string& partialOutput = ((_output.size() > 1) ? _output[0] : _input[0]);
140 boost::filesystem::path dirF(partialOutput);
141 boost::filesystem::create_directory(dirF);
142 boost::filesystem::path trjfileF(partialOutput +
".irc.forward.trj.xyz");
143 std::ofstream trajectoryF((dirF / trjfileF).
string(), std::ofstream::out);
144 auto forward = [&](
const int& cycle,
const double& energy,
const Eigen::VectorXd& params) {
145 if (oldParams.size() != params.size()) {
149 cout.printf(
"%7s %16s %16s %16s %16s\n",
"Cycle",
"Energy",
"Energy Diff.",
"Step RMS",
"Max. Step");
151 auto diff = (params - oldParams).eval();
152 cout.printf(
"%7d %+16.9f %+16.9f %+16.9f %+16.9f\n", cycle, energy, energy - oldEnergy,
153 sqrt(diff.squaredNorm() / diff.size()), diff.cwiseAbs().maxCoeff());
156 auto structure = calc->getStructure();
157 Writer::write(trajectoryF, *structure, std::to_string(energy));
159 optimizer->addObserver(forward);
162 auto structure = systems.at(_input[0])->getStructure();
164 bool forwardAborted =
false;
166 cycles = optimizer->optimize(*structure, *_logger, ircModeVector,
true);
169 Writer::write(trajectoryF, *calc->getStructure());
175 _logger->error << boost::current_exception_diagnostic_information() <<
Core::Log::endl;
176 forwardAborted =
true;
180 int maxiter = settings.getInt(
"convergence_max_iterations");
181 const bool forwardConverged = cycles < maxiter && !forwardAborted;
182 if (forwardConverged) {
187 else if (!forwardAborted) {
192 throw std::runtime_error(
"Problem: IRC optimization did not converge.");
197 if (!forwardAborted) {
199 if (_output.size() > 1) {
200 systems[_output[0]] = std::shared_ptr<Core::Calculator>(calc->clone().release());
201 boost::filesystem::path xyzfile(_output[0] +
".xyz");
202 std::ofstream xyz((dirF / xyzfile).
string(), std::ofstream::out);
203 Writer::write(xyz, *(calc->getStructure()));
207 boost::filesystem::path xyzfile(_input[0] +
".irc.forward.xyz");
208 std::ofstream xyz((dirF / xyzfile).
string(), std::ofstream::out);
209 Writer::write(xyz, *(calc->getStructure()));
217 cout <<
"\n IRC Optimization: Backward\n\n";
222 optimizer->setSettings(settings);
225 boost::filesystem::path dirB(((_output.size() > 1) ? _output[1] : _input[0]));
226 boost::filesystem::create_directory(dirB);
227 boost::filesystem::path trjfileB(((_output.size() > 1) ? _output[1] : _input[0]) +
".irc.backward.trj.xyz");
228 std::ofstream trajectoryB((dirB / trjfileB).
string(), std::ofstream::out);
229 auto backward = [&](
const int& cycle,
const double& energy,
const Eigen::VectorXd& params) {
230 if (oldParams.size() != params.size()) {
234 cout.printf(
"%7s %16s %16s %16s %16s\n",
"Cycle",
"Energy",
"Energy Diff.",
"Step RMS",
"Max. Step");
236 auto diff = (params - oldParams).eval();
237 cout.printf(
"%7d %+16.9f %+16.9f %+16.9f %+16.9f\n", cycle, energy, energy - oldEnergy,
238 sqrt(diff.squaredNorm() / diff.size()), diff.cwiseAbs().maxCoeff());
241 auto structure = calc->getStructure();
242 Writer::write(trajectoryB, *structure, std::to_string(energy));
244 optimizer->clearObservers();
245 optimizer->addObserver(backward);
248 structure = systems.at(_input[0])->getStructure();
251 cycles = optimizer->optimize(*structure, *_logger, ircModeVector,
false);
254 Writer::write(trajectoryB, *calc->getStructure());
256 _logger->error <<
"Calculation in IRC backward optimization failed!" <<
Core::Log::endl;
260 _logger->error << boost::current_exception_diagnostic_information() <<
Core::Log::endl;
265 maxiter = settings.getInt(
"convergence_max_iterations");
266 if (cycles < maxiter) {
276 throw std::runtime_error(
"Problem: IRC optimization did not converge.");
281 if (_output.size() > 1) {
282 systems[_output[1]] = std::shared_ptr<Core::Calculator>(calc->clone().release());
283 boost::filesystem::path xyzfile(_output[1] +
".xyz");
284 std::ofstream xyz((dirB / xyzfile).
string(), std::ofstream::out);
285 Writer::write(xyz, *(calc->getStructure()));
289 boost::filesystem::path xyzfile(_input[0] +
".irc.backward.xyz");
290 std::ofstream xyz((dirB / xyzfile).
string(), std::ofstream::out);
291 Writer::write(xyz, *(calc->getStructure()));
295 return cycles < maxiter && forwardConverged;
302 #endif // READUCT_IRCTASK_H_
static std::ostream & endl(std::ostream &os)
const std::vector< std::string > & input() const
Getter for the expected names of the input systems.
Definition: Task.h:77
bool run(SystemsMap &systems, Utils::UniversalSettings::ValueCollection taskSettings, bool testRunOnly=false) const final
Executes the actual task represented by this class.
Definition: IrcTask.h:48
IrcTask(std::vector< std::string > input, std::vector< std::string > output, std::shared_ptr< Core::Log > logger=nullptr)
Construct a new IrcTask.
Definition: IrcTask.h:40
void warningIfMultipleInputsGiven() const
Warn if more than one input system was specified.
Definition: Task.h:90
std::string name() const override
Getter for the tasks name.
Definition: IrcTask.h:44
const std::vector< std::string > & output() const
Getter for the names of the output systems generated by this task.
Definition: Task.h:84
The base class for all tasks in Readuct.
Definition: Task.h:28
A task to run an IRC optimization using a given normal mode.
Definition: IrcTask.h:32