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";
50 observers = {})
const final {
53 if (_output.size() == 1) {
55 <<
" Warning: To store IRC results in a new location two output systems have to be specified.\n";
58 bool silentCalculator = taskSettings.extract(
"silent_stdout_calculator",
true);
59 std::shared_ptr<Core::Calculator> calc;
62 calc = copyCalculator(systems, _input.front(),
name());
63 Utils::CalculationRoutines::setLog(*calc,
true,
true, !silentCalculator);
66 if (calc->getStructure()->size() == 1) {
67 throw std::runtime_error(
"Cannot perform IRC task for monoatomic systems.");
72 std::string optimizertype = taskSettings.extract(
"optimizer", std::string{
"SD"});
73 int mode = taskSettings.extract(
"irc_mode", 0);
75 bool stopOnError = stopOnErrorExtraction(taskSettings);
76 std::transform(optimizertype.begin(), optimizertype.end(), optimizertype.begin(), ::toupper);
77 std::shared_ptr<Utils::IrcOptimizerBase> optimizer;
78 if (optimizertype ==
"LBFGS") {
79 auto tmp = std::make_shared<Utils::IrcOptimizer<Utils::Lbfgs>>(*calc);
81 optimizer = std::move(tmp);
83 else if (optimizertype ==
"BFGS") {
84 auto tmp = std::make_shared<Utils::IrcOptimizer<Utils::Bfgs>>(*calc);
86 optimizer = std::move(tmp);
88 else if (optimizertype ==
"SD" || optimizertype ==
"STEEPESTDESCENT") {
89 auto tmp = std::make_shared<Utils::IrcOptimizer<Utils::SteepestDescent>>(*calc);
91 optimizer = std::move(tmp);
94 throw std::runtime_error(
95 "Unknown Optimizer requested for an IRC optimization, available are: SD, BFGS and LBFGS!");
98 auto settings = optimizer->getSettings();
99 settings.merge(taskSettings);
100 if (!settings.valid()) {
101 settings.throwIncorrectSettings();
103 optimizer->setSettings(settings);
104 if (!testRunOnly && !Utils::GeometryOptimization::settingsMakeSense(*optimizer)) {
105 throw std::logic_error(
"The given calculator settings are too inaccurate for the given convergence criteria of "
106 "this optimization Task");
115 if (!calc->results().has<Utils::Property::Hessian>()) {
116 calc->setRequiredProperties(Utils::Property::Hessian);
117 calc->calculate(
"Hessian Calculation");
119 auto hessian = calc->results().get<Utils::Property::Hessian>();
120 auto system = calc->getStructure();
121 auto modes = Utils::NormalModeAnalysis::calculateNormalModes(hessian, system->getElements(), system->getPositions());
122 if (mode < 0 || mode > modes.size() - 1) {
123 throw std::runtime_error(
124 "The chosen normal mode number is smaller than 0 or larger than the number of modes present!");
126 auto ircMode = modes.getMode(mode);
127 auto ircModeVector = Eigen::Map<const Eigen::VectorXd>(ircMode.data(), ircMode.cols() * ircMode.rows());
132 auto cout = _logger->output;
133 cout <<
"\n IRC Optimization: Forward\n\n";
134 std::shared_ptr<Scine::Core::Calculator> forwardCalc = std::move(calc->clone());
136 double oldEnergy = 0.0;
137 Eigen::VectorXd oldParams;
139 using Writer = Utils::XyzStreamHandler;
141 const std::string& partialOutput = ((_output.size() > 1) ? _output[0] : _input[0]);
142 boost::filesystem::path dirF(partialOutput);
143 boost::filesystem::create_directory(dirF);
144 boost::filesystem::path trjfileF(partialOutput +
".irc.forward.trj.xyz");
145 std::ofstream trajectoryF((dirF / trjfileF).
string(), std::ofstream::out);
146 auto forward = [&](
const int& cycle,
const double& energy,
const Eigen::VectorXd& params) {
147 if (oldParams.size() != params.size()) {
151 cout.printf(
"%7s %16s %16s %16s %16s\n",
"Cycle",
"Energy",
"Energy Diff.",
"Step RMS",
"Max. Step");
153 auto diff = (params - oldParams).eval();
154 cout.printf(
"%7d %+16.9f %+16.9f %+16.9f %+16.9f\n", cycle, energy, energy - oldEnergy,
155 sqrt(diff.squaredNorm() / diff.size()), diff.cwiseAbs().maxCoeff());
158 auto structure = calc->getStructure();
159 Writer::write(trajectoryF, *structure, std::to_string(energy));
161 optimizer->addObserver(forward);
164 auto customObserversForward = [&calc, &observers](
const int& cycle,
const double& ,
165 const Eigen::VectorXd& ) {
166 for (
auto& observer : observers) {
167 auto atoms = calc->getStructure();
168 Utils::Results& results = calc->results();
169 observer(cycle, *atoms, results,
"irc_forward");
172 optimizer->addObserver(customObserversForward);
175 auto structure = systems.at(_input[0])->getStructure();
177 bool forwardAborted =
false;
179 cycles = optimizer->optimize(*structure, *_logger, ircModeVector,
true);
182 Writer::write(trajectoryF, *calc->getStructure());
188 _logger->error << boost::current_exception_diagnostic_information() <<
Core::Log::endl;
189 forwardAborted =
true;
193 int maxiter = settings.getInt(
"convergence_max_iterations");
194 const bool forwardConverged = cycles < maxiter && !forwardAborted;
195 if (forwardConverged) {
200 else if (!forwardAborted) {
205 throw std::runtime_error(
"Problem: IRC optimization did not converge.");
210 if (!forwardAborted) {
212 if (_output.size() > 1) {
213 systems[_output[0]] = calc->clone();
214 boost::filesystem::path xyzfile(_output[0] +
".xyz");
215 std::ofstream xyz((dirF / xyzfile).
string(), std::ofstream::out);
216 Writer::write(xyz, *(calc->getStructure()));
220 boost::filesystem::path xyzfile(_input[0] +
".irc.forward.xyz");
221 std::ofstream xyz((dirF / xyzfile).
string(), std::ofstream::out);
222 Writer::write(xyz, *(calc->getStructure()));
230 cout <<
"\n IRC Optimization: Backward\n\n";
235 optimizer->setSettings(settings);
238 boost::filesystem::path dirB(((_output.size() > 1) ? _output[1] : _input[0]));
239 boost::filesystem::create_directory(dirB);
240 boost::filesystem::path trjfileB(((_output.size() > 1) ? _output[1] : _input[0]) +
".irc.backward.trj.xyz");
241 std::ofstream trajectoryB((dirB / trjfileB).
string(), std::ofstream::out);
242 auto backward = [&](
const int& cycle,
const double& energy,
const Eigen::VectorXd& params) {
243 if (oldParams.size() != params.size()) {
247 cout.printf(
"%7s %16s %16s %16s %16s\n",
"Cycle",
"Energy",
"Energy Diff.",
"Step RMS",
"Max. Step");
249 auto diff = (params - oldParams).eval();
250 cout.printf(
"%7d %+16.9f %+16.9f %+16.9f %+16.9f\n", cycle, energy, energy - oldEnergy,
251 sqrt(diff.squaredNorm() / diff.size()), diff.cwiseAbs().maxCoeff());
254 auto structure = calc->getStructure();
255 Writer::write(trajectoryB, *structure, std::to_string(energy));
257 optimizer->clearObservers();
258 optimizer->addObserver(backward);
261 auto customObserversBackward = [&calc, &observers](
const int& cycle,
const double& ,
262 const Eigen::VectorXd& ) {
263 for (
auto& observer : observers) {
264 auto atoms = calc->getStructure();
265 Utils::Results& results = calc->results();
266 observer(cycle, *atoms, results,
"irc_backward");
269 optimizer->addObserver(customObserversBackward);
272 structure = systems.at(_input[0])->getStructure();
275 cycles = optimizer->optimize(*structure, *_logger, ircModeVector,
false);
278 Writer::write(trajectoryB, *calc->getStructure());
280 _logger->error <<
"Calculation in IRC backward optimization failed!" <<
Core::Log::endl;
284 _logger->error << boost::current_exception_diagnostic_information() <<
Core::Log::endl;
289 maxiter = settings.getInt(
"convergence_max_iterations");
290 if (cycles < maxiter) {
291 cout << Core::Log::endl
292 <<
" Converged after " << cycles <<
" iterations." << Core::Log::endl
296 cout << Core::Log::endl
297 <<
" Stopped after " << maxiter <<
" iterations." << Core::Log::endl
300 throw std::runtime_error(
"Problem: IRC optimization did not converge.");
305 if (_output.size() > 1) {
306 systems[_output[1]] = calc->clone();
307 boost::filesystem::path xyzfile(_output[1] +
".xyz");
308 std::ofstream xyz((dirB / xyzfile).
string(), std::ofstream::out);
309 Writer::write(xyz, *(calc->getStructure()));
313 boost::filesystem::path xyzfile(_input[0] +
".irc.backward.xyz");
314 std::ofstream xyz((dirB / xyzfile).
string(), std::ofstream::out);
315 Writer::write(xyz, *(calc->getStructure()));
319 return cycles < maxiter && forwardConverged;
326 #endif // READUCT_IRCTASK_H_
bool run(SystemsMap &systems, Utils::UniversalSettings::ValueCollection taskSettings, bool testRunOnly=false, std::vector< std::function< void(const int &, const Utils::AtomCollection &, const Utils::Results &, const std::string &)>> observers={}) const final
Executes the actual task represented by this class.
Definition: IrcTask.h:48
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:82
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:95
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:89
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