Scine::Readuct  4.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) const final {
50  // Warn if only one output system was specified
51  if (_output.size() == 1) {
52  _logger->warning
53  << " Warning: To store IRC results in a new location two output systems have to be specified.\n";
54  }
55 
56  bool silentCalculator = taskSettings.extract("silent_stdout_calculator", true);
57  std::shared_ptr<Core::Calculator> calc;
58  if (!testRunOnly) { // leave out in case of task chaining --> attention calc is NULL
59  // Note: _input is guaranteed not to be empty by Task constructor
60  calc = copyCalculator(systems, _input.front(), name());
61  Utils::CalculationRoutines::setLog(*calc, true, true, !silentCalculator);
62 
63  // Check system size
64  if (calc->getStructure()->size() == 1) {
65  throw std::runtime_error("Cannot perform IRC task for monoatomic systems.");
66  }
67  }
68 
69  // Generate optimizer
70  std::string optimizertype = taskSettings.extract("optimizer", std::string{"SD"});
71  int mode = taskSettings.extract("irc_mode", 0);
72  // Read and delete special settings
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);
78  // Default convergence options
79  optimizer = std::move(tmp);
80  }
81  else if (optimizertype == "BFGS") {
82  auto tmp = std::make_shared<Utils::IrcOptimizer<Utils::Bfgs>>(*calc);
83  // Default convergence options
84  optimizer = std::move(tmp);
85  }
86  else if (optimizertype == "SD" || optimizertype == "STEEPESTDESCENT") {
87  auto tmp = std::make_shared<Utils::IrcOptimizer<Utils::SteepestDescent>>(*calc);
88  // Default convergence options
89  optimizer = std::move(tmp);
90  }
91  else {
92  throw std::runtime_error(
93  "Unknown Optimizer requested for an IRC optimization, available are: SD, BFGS and LBFGS!");
94  }
95  // Apply user settings
96  auto settings = optimizer->getSettings();
97  settings.merge(taskSettings);
98  if (!settings.valid()) {
99  settings.throwIncorrectSettings();
100  }
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");
105  }
106 
107  // If no errors encountered until here, the basic settings should be alright
108  if (testRunOnly) {
109  return true;
110  }
111 
112  // Get/find mode
113  if (!calc->results().has<Utils::Property::Hessian>()) {
114  calc->setRequiredProperties(Utils::Property::Hessian);
115  calc->calculate("Hessian Calculation");
116  }
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!");
123  }
124  auto ircMode = modes.getMode(mode);
125  auto ircModeVector = Eigen::Map<const Eigen::VectorXd>(ircMode.data(), ircMode.cols() * ircMode.rows());
126 
127  /*========================*
128  * Forward optimization
129  *========================*/
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());
133  // Add observer
134  double oldEnergy = 0.0;
135  Eigen::VectorXd oldParams;
136  // Trajectory stream
137  using Writer = Utils::XyzStreamHandler;
138 
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()) {
146  oldParams = params;
147  }
148  if (cycle == 1) {
149  cout.printf("%7s %16s %16s %16s %16s\n", "Cycle", "Energy", "Energy Diff.", "Step RMS", "Max. Step");
150  }
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());
154  oldEnergy = energy;
155  oldParams = params;
156  auto structure = calc->getStructure();
157  Writer::write(trajectoryF, *structure, std::to_string(energy));
158  };
159  optimizer->addObserver(forward);
160 
161  // Run optimization
162  auto structure = systems.at(_input[0])->getStructure();
163  int cycles = 0;
164  bool forwardAborted = false;
165  try {
166  cycles = optimizer->optimize(*structure, *_logger, ircModeVector, true);
167  }
168  catch (...) {
169  Writer::write(trajectoryF, *calc->getStructure());
170  _logger->error << "Forward IRC failed with error!" << Core::Log::endl;
171  if (stopOnError) {
172  trajectoryF.close();
173  throw;
174  }
175  _logger->error << boost::current_exception_diagnostic_information() << Core::Log::endl;
176  forwardAborted = true;
177  }
178  trajectoryF.close();
179 
180  int maxiter = settings.getInt("convergence_max_iterations");
181  const bool forwardConverged = cycles < maxiter && !forwardAborted;
182  if (forwardConverged) {
183  cout << Core::Log::endl
184  << " Converged after " << cycles << " iterations." << Core::Log::endl
185  << Core::Log::endl;
186  }
187  else if (!forwardAborted) {
188  cout << Core::Log::endl
189  << " Stopped after " << maxiter << " iterations." << Core::Log::endl
190  << Core::Log::endl;
191  if (stopOnError) {
192  throw std::runtime_error("Problem: IRC optimization did not converge.");
193  }
194  }
195 
196  // only write separate file with single last frame if no error occurred
197  if (!forwardAborted) {
198  // Print/Store results
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()));
204  xyz.close();
205  }
206  else {
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()));
210  xyz.close();
211  }
212  }
213 
214  /*=========================*
215  * Backward optimization
216  *=========================*/
217  cout << "\n IRC Optimization: Backward\n\n";
218  // Add observer
219  oldEnergy = 0.0;
220  oldParams.resize(0);
221  // Reset optimizer
222  optimizer->setSettings(settings);
223  optimizer->reset();
224  // Trajectory stream
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()) {
231  oldParams = params;
232  }
233  if (cycle == 1) {
234  cout.printf("%7s %16s %16s %16s %16s\n", "Cycle", "Energy", "Energy Diff.", "Step RMS", "Max. Step");
235  }
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());
239  oldEnergy = energy;
240  oldParams = params;
241  auto structure = calc->getStructure();
242  Writer::write(trajectoryB, *structure, std::to_string(energy));
243  };
244  optimizer->clearObservers();
245  optimizer->addObserver(backward);
246 
247  // Run optimization
248  structure = systems.at(_input[0])->getStructure();
249  cycles = 0;
250  try {
251  cycles = optimizer->optimize(*structure, *_logger, ircModeVector, false);
252  }
253  catch (...) {
254  Writer::write(trajectoryB, *calc->getStructure());
255  trajectoryB.close();
256  _logger->error << "Calculation in IRC backward optimization failed!" << Core::Log::endl;
257  if (stopOnError) {
258  throw;
259  }
260  _logger->error << boost::current_exception_diagnostic_information() << Core::Log::endl;
261  return false;
262  }
263  trajectoryB.close();
264 
265  maxiter = settings.getInt("convergence_max_iterations");
266  if (cycles < maxiter) {
267  cout << Core::Log::endl
268  << " Converged after " << cycles << " iterations." << Core::Log::endl
269  << Core::Log::endl;
270  }
271  else {
272  cout << Core::Log::endl
273  << " Stopped after " << maxiter << " iterations." << Core::Log::endl
274  << Core::Log::endl;
275  if (stopOnError) {
276  throw std::runtime_error("Problem: IRC optimization did not converge.");
277  }
278  }
279 
280  // Print/Store results
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()));
286  xyz.close();
287  }
288  else {
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()));
292  xyz.close();
293  }
294 
295  return cycles < maxiter && forwardConverged;
296  }
297 };
298 
299 } // namespace Readuct
300 } // namespace Scine
301 
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