Scine::Readuct  5.0.0
This is the SCINE module ReaDuct.
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Pages
IrcTask.h
Go to the documentation of this file.
1 
7 #ifndef READUCT_IRCTASK_H_
8 #define READUCT_IRCTASK_H_
9 
10 /* Readuct */
11 #include "Tasks/Task.h"
12 /* Scine */
20 /* std c++ */
21 #include <boost/exception/diagnostic_information.hpp>
22 #include <boost/filesystem.hpp>
23 #include <cstdio>
24 #include <fstream>
25 
26 namespace Scine {
27 namespace Readuct {
28 
32 class IrcTask : public Task {
33  public:
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)) {
42  }
43 
44  std::string name() const override {
45  return "IRC Optimizations";
46  }
47 
48  bool run(SystemsMap& systems, Utils::UniversalSettings::ValueCollection taskSettings, bool testRunOnly = false,
49  std::vector<std::function<void(const int&, const Utils::AtomCollection&, const Utils::Results&, const std::string&)>>
50  observers = {}) const final {
52  // Warn if only one output system was specified
53  if (_output.size() == 1) {
54  _logger->warning
55  << " Warning: To store IRC results in a new location two output systems have to be specified.\n";
56  }
57 
58  bool silentCalculator = taskSettings.extract("silent_stdout_calculator", true);
59  std::shared_ptr<Core::Calculator> calc;
60  if (!testRunOnly) { // leave out in case of task chaining --> attention calc is NULL
61  // Note: _input is guaranteed not to be empty by Task constructor
62  calc = copyCalculator(systems, _input.front(), name());
63  Utils::CalculationRoutines::setLog(*calc, true, true, !silentCalculator);
64 
65  // Check system size
66  if (calc->getStructure()->size() == 1) {
67  throw std::runtime_error("Cannot perform IRC task for monoatomic systems.");
68  }
69  }
70 
71  // Generate optimizer
72  std::string optimizertype = taskSettings.extract("optimizer", std::string{"SD"});
73  int mode = taskSettings.extract("irc_mode", 0);
74  // Read and delete special settings
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);
80  // Default convergence options
81  optimizer = std::move(tmp);
82  }
83  else if (optimizertype == "BFGS") {
84  auto tmp = std::make_shared<Utils::IrcOptimizer<Utils::Bfgs>>(*calc);
85  // Default convergence options
86  optimizer = std::move(tmp);
87  }
88  else if (optimizertype == "SD" || optimizertype == "STEEPESTDESCENT") {
89  auto tmp = std::make_shared<Utils::IrcOptimizer<Utils::SteepestDescent>>(*calc);
90  // Default convergence options
91  optimizer = std::move(tmp);
92  }
93  else {
94  throw std::runtime_error(
95  "Unknown Optimizer requested for an IRC optimization, available are: SD, BFGS and LBFGS!");
96  }
97  // Apply user settings
98  auto settings = optimizer->getSettings();
99  settings.merge(taskSettings);
100  if (!settings.valid()) {
101  settings.throwIncorrectSettings();
102  }
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");
107  }
108 
109  // If no errors encountered until here, the basic settings should be alright
110  if (testRunOnly) {
111  return true;
112  }
113 
114  // Get/find mode
115  if (!calc->results().has<Utils::Property::Hessian>()) {
116  calc->setRequiredProperties(Utils::Property::Hessian);
117  calc->calculate("Hessian Calculation");
118  }
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!");
125  }
126  auto ircMode = modes.getMode(mode);
127  auto ircModeVector = Eigen::Map<const Eigen::VectorXd>(ircMode.data(), ircMode.cols() * ircMode.rows());
128 
129  /*========================*
130  * Forward optimization
131  *========================*/
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());
135  // Add observer
136  double oldEnergy = 0.0;
137  Eigen::VectorXd oldParams;
138  // Trajectory stream
139  using Writer = Utils::XyzStreamHandler;
140 
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()) {
148  oldParams = params;
149  }
150  if (cycle == 1) {
151  cout.printf("%7s %16s %16s %16s %16s\n", "Cycle", "Energy", "Energy Diff.", "Step RMS", "Max. Step");
152  }
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());
156  oldEnergy = energy;
157  oldParams = params;
158  auto structure = calc->getStructure();
159  Writer::write(trajectoryF, *structure, std::to_string(energy));
160  };
161  optimizer->addObserver(forward);
162 
163  // Add custom observers forward
164  auto customObserversForward = [&calc, &observers](const int& cycle, const double& /*energy*/,
165  const Eigen::VectorXd& /*params*/) {
166  for (auto& observer : observers) {
167  auto atoms = calc->getStructure();
168  Utils::Results& results = calc->results();
169  observer(cycle, *atoms, results, "irc_forward");
170  }
171  };
172  optimizer->addObserver(customObserversForward);
173 
174  // Run optimization
175  auto structure = systems.at(_input[0])->getStructure();
176  int cycles = 0;
177  bool forwardAborted = false;
178  try {
179  cycles = optimizer->optimize(*structure, *_logger, ircModeVector, true);
180  }
181  catch (...) {
182  Writer::write(trajectoryF, *calc->getStructure());
183  _logger->error << "Forward IRC failed with error!" << Core::Log::endl;
184  if (stopOnError) {
185  trajectoryF.close();
186  throw;
187  }
188  _logger->error << boost::current_exception_diagnostic_information() << Core::Log::endl;
189  forwardAborted = true;
190  }
191  trajectoryF.close();
192 
193  int maxiter = settings.getInt("convergence_max_iterations");
194  const bool forwardConverged = cycles < maxiter && !forwardAborted;
195  if (forwardConverged) {
196  cout << Core::Log::endl
197  << " Converged after " << cycles << " iterations." << Core::Log::endl
198  << Core::Log::endl;
199  }
200  else if (!forwardAborted) {
201  cout << Core::Log::endl
202  << " Stopped after " << maxiter << " iterations." << Core::Log::endl
203  << Core::Log::endl;
204  if (stopOnError) {
205  throw std::runtime_error("Problem: IRC optimization did not converge.");
206  }
207  }
208 
209  // only write separate file with single last frame if no error occurred
210  if (!forwardAborted) {
211  // Print/Store results
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()));
217  xyz.close();
218  }
219  else {
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()));
223  xyz.close();
224  }
225  }
226 
227  /*=========================*
228  * Backward optimization
229  *=========================*/
230  cout << "\n IRC Optimization: Backward\n\n";
231  // Add observer
232  oldEnergy = 0.0;
233  oldParams.resize(0);
234  // Reset optimizer
235  optimizer->setSettings(settings);
236  optimizer->reset();
237  // Trajectory stream
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()) {
244  oldParams = params;
245  }
246  if (cycle == 1) {
247  cout.printf("%7s %16s %16s %16s %16s\n", "Cycle", "Energy", "Energy Diff.", "Step RMS", "Max. Step");
248  }
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());
252  oldEnergy = energy;
253  oldParams = params;
254  auto structure = calc->getStructure();
255  Writer::write(trajectoryB, *structure, std::to_string(energy));
256  };
257  optimizer->clearObservers();
258  optimizer->addObserver(backward);
259 
260  // Add custom observers backward
261  auto customObserversBackward = [&calc, &observers](const int& cycle, const double& /*energy*/,
262  const Eigen::VectorXd& /*params*/) {
263  for (auto& observer : observers) {
264  auto atoms = calc->getStructure();
265  Utils::Results& results = calc->results();
266  observer(cycle, *atoms, results, "irc_backward");
267  }
268  };
269  optimizer->addObserver(customObserversBackward);
270 
271  // Run optimization
272  structure = systems.at(_input[0])->getStructure();
273  cycles = 0;
274  try {
275  cycles = optimizer->optimize(*structure, *_logger, ircModeVector, false);
276  }
277  catch (...) {
278  Writer::write(trajectoryB, *calc->getStructure());
279  trajectoryB.close();
280  _logger->error << "Calculation in IRC backward optimization failed!" << Core::Log::endl;
281  if (stopOnError) {
282  throw;
283  }
284  _logger->error << boost::current_exception_diagnostic_information() << Core::Log::endl;
285  return false;
286  }
287  trajectoryB.close();
288 
289  maxiter = settings.getInt("convergence_max_iterations");
290  if (cycles < maxiter) {
291  cout << Core::Log::endl
292  << " Converged after " << cycles << " iterations." << Core::Log::endl
293  << Core::Log::endl;
294  }
295  else {
296  cout << Core::Log::endl
297  << " Stopped after " << maxiter << " iterations." << Core::Log::endl
298  << Core::Log::endl;
299  if (stopOnError) {
300  throw std::runtime_error("Problem: IRC optimization did not converge.");
301  }
302  }
303 
304  // Print/Store results
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()));
310  xyz.close();
311  }
312  else {
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()));
316  xyz.close();
317  }
318 
319  return cycles < maxiter && forwardConverged;
320  }
321 };
322 
323 } // namespace Readuct
324 } // namespace Scine
325 
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