Scine::Readuct  4.0.0
This is the SCINE module Readuct.
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Pages
TsOptimizationTask.h
Go to the documentation of this file.
1 
7 #ifndef READUCT_TSOPTIMIZATIONTASK_H_
8 #define READUCT_TSOPTIMIZATIONTASK_H_
9 
10 /* Readuct */
11 #include "Tasks/Task.h"
12 /* Scine */
23 /* External */
24 #include <boost/exception/diagnostic_information.hpp>
25 #include <boost/filesystem.hpp>
26 /* std c++ */
27 #include <cstdio>
28 #include <fstream>
29 #include <iostream>
30 #include <vector>
31 
32 namespace Scine {
33 namespace Readuct {
34 
35 class TsOptimizationTask : public Task {
36  public:
43  TsOptimizationTask(std::vector<std::string> input, std::vector<std::string> output, std::shared_ptr<Core::Log> logger = nullptr)
44  : Task(std::move(input), std::move(output), std::move(logger)) {
45  }
46 
47  std::string name() const override {
48  return "TS Optimization";
49  }
50 
51  bool run(SystemsMap& systems, Utils::UniversalSettings::ValueCollection taskSettings, bool testRunOnly = false) const final {
54 
55  bool silentCalculator = taskSettings.extract("silent_stdout_calculator", true);
56  std::shared_ptr<Core::Calculator> calc;
57  if (!testRunOnly) { // leave out in case of task chaining --> attention calc is NULL
58  // Note: _input is guaranteed not to be empty by Task constructor
59  calc = copyCalculator(systems, _input.front(), name());
60  Utils::CalculationRoutines::setLog(*calc, true, true, !silentCalculator);
61 
62  // Check system size
63  if (calc->getStructure()->size() == 1) {
64  throw std::runtime_error("Cannot calculate transition state for monoatomic systems.");
65  }
66  }
67 
68  // Generate optimizer
69  std::string optimizertype = taskSettings.extract("optimizer", std::string{"bofill"});
70  std::transform(optimizertype.begin(), optimizertype.end(), optimizertype.begin(), ::tolower);
71 
72  // Check for automatic mode selection
73  // sanity check
74  if (taskSettings.valueExists("automatic_mode_selection") &&
75  (taskSettings.valueExists(optimizertype + "_follow_mode") || taskSettings.valueExists("ev_follow_mode"))) {
76  throw std::logic_error("You specified an automatic selection of the mode and gave a specific mode yourself. "
77  "This is not possible. Only give one of those options.");
78  }
79  std::vector<int> relevantAtoms = taskSettings.extract("automatic_mode_selection", std::vector<int>{});
80  // another sanity check
81  if (!relevantAtoms.empty()) {
82  for (const auto& index : relevantAtoms) {
83  if (index < 0) {
84  throw std::logic_error("You gave an atom index smaller than 0 in automatic_mode_selection. "
85  "This does not make sense.");
86  }
87  }
88  }
89 
90  using XyzHandler = Utils::XyzStreamHandler;
91  auto cout = _logger->output;
92 
93  std::shared_ptr<Utils::GeometryOptimizerBase> optimizer;
94  if (optimizertype == "bofill") {
95  auto tmp = std::make_shared<Utils::GeometryOptimizer<Utils::Bofill>>(*calc);
96  optimizer = std::move(tmp);
97  }
98  else if (optimizertype == "ef" || optimizertype == "ev" || optimizertype == "evf" ||
99  optimizertype == "eigenvectorfollowing" || optimizertype == "eigenvector_following") {
100  auto tmp = std::make_shared<Utils::GeometryOptimizer<Utils::EigenvectorFollowing>>(*calc);
101  optimizer = std::move(tmp);
102  }
103  else if (optimizertype == "dimer") {
104  auto tmp = std::make_shared<Utils::GeometryOptimizer<Utils::Dimer>>(*calc);
105 
106  /* get if coordinates are transformed to know whether a provided guess vector has to be transformed */
107  auto coordinateSystem = tmp->coordinateSystem;
108  if (taskSettings.valueExists("geoopt_coordinate_system")) {
109  coordinateSystem = Utils::CoordinateSystemInterpreter::getCoordinateSystemFromString(
110  taskSettings.getString("geoopt_coordinate_system"));
111  }
112 
113  /* Check for possible guess vectors */
114  bool useEigenvectorForFirstStep = taskSettings.extract("dimer_calculate_hessian_once", false);
115  if (!useEigenvectorForFirstStep && !relevantAtoms.empty()) {
116  _logger->warning << "Specified Dimer optimizer with automatic mode selection, this means that a Hessian has to "
117  "be calculated.\n";
118  }
119  /* set to true if mode is specified or automated selection of mode is set in Settings */
120  if (taskSettings.valueExists("dimer_follow_mode") || !relevantAtoms.empty()) {
121  useEigenvectorForFirstStep = true;
122  }
123  int followedMode = taskSettings.extract("dimer_follow_mode", 0);
124  std::string guessVectorFileName = taskSettings.extract("dimer_guess_vector_file", std::string(""));
125  std::vector<std::string> discreteGuessesFileNames =
126  taskSettings.extract("dimer_discrete_guesses", std::vector<std::string>{"", ""});
127 
128  /* sanity checks for Dimer init options */
129  if (discreteGuessesFileNames.size() != 2) {
130  throw std::runtime_error("Problem with discrete guesses. You may only give list with two filepaths.");
131  }
132  // calc cannot be checked in test run
133  if (!testRunOnly && !useEigenvectorForFirstStep && calc->results().has<Utils::Property::Hessian>()) {
134  _logger->warning << "Did not specify to use eigenvector for initialization of Dimer, although your system "
135  "already includes a Hessian. The Hessian will therefore be ignored. "
136  "If you want to change that, activate the setting 'dimer_calculate_hessian_once': true"
137  << Core::Log::endl;
138  }
139  if (useEigenvectorForFirstStep) {
140  if (!guessVectorFileName.empty()) {
141  throw std::logic_error("You specified to use the eigenvector of the Hessian as Dimer initialization, "
142  "but also specified a file to read a guess vector for the initialization. "
143  "Only select one.");
144  }
145  if (!discreteGuessesFileNames[0].empty()) {
146  throw std::logic_error("You specified to use the eigenvector of the Hessian as Dimer initialization, "
147  "but also specified files to read discrete guesses for the initialization. "
148  "Only select one.");
149  }
150  }
151  else if (!guessVectorFileName.empty() && !discreteGuessesFileNames[0].empty()) {
152  throw std::logic_error("You specified a file to read a guess vector for the Dimer initialization. "
153  "but also specified files to read discrete guesses for the initialization. "
154  "Only select one.");
155  }
156 
157  if (testRunOnly) {
158  ; // skip all dimer if statements here, because they require a Hessian calculation or file look-up
159  }
160  /* Generate and pass initial eigenvector and tell optimizer to not perform the first rotation */
161  else if (useEigenvectorForFirstStep) {
162  tmp->optimizer.skipFirstRotation = true;
163  auto modes = getNormalModes(*calc);
164  auto system = calc->getStructure();
165  if (!relevantAtoms.empty()) { // clash with dimer_follow_mode option is already checked prior
166  followedMode = selectMode(*calc, relevantAtoms);
167  }
168  auto mode = modes.getMode(followedMode); // range check for followedMode is included here
169  /* Currently Hessian not possible in internal coordinates */
170  if (coordinateSystem == Utils::CoordinateSystem::Internal) {
171  /* Transform vector of normal mode into internal coordinates */
172  Utils::InternalCoordinates internalCoordinates(*system);
173  tmp->optimizer.guessVector = std::make_shared<Eigen::VectorXd>(internalCoordinates.coordinatesToInternal(mode));
174  }
175  else if (coordinateSystem == Utils::CoordinateSystem::CartesianWithoutRotTrans) {
176  /* Can give proper invH for all but internal */
177  /* Calculate inverse Hessian with one mode flipped for BFGS within Dimer optimizer */
178  Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> es;
179  auto hessian = calc->results().get<Utils::Property::Hessian>(); // guaranteed to have it from getNormalModes
180  es.compute(hessian);
181  int nEigenValues = es.eigenvalues().size();
182  Eigen::MatrixXd eigenValueMatrix = Eigen::MatrixXd::Identity(nEigenValues, nEigenValues);
183  eigenValueMatrix.diagonal() = es.eigenvalues();
184  eigenValueMatrix(followedMode, followedMode) = -es.eigenvalues()[followedMode];
185  Eigen::MatrixXd switchedHessian = es.eigenvectors() * eigenValueMatrix * es.eigenvectors().transpose();
186  /* Transform vector of normal mode into internal coordinates */
187  Utils::InternalCoordinates internalCoordinates(*system, true);
188  tmp->optimizer.guessVector = std::make_shared<Eigen::VectorXd>(internalCoordinates.coordinatesToInternal(mode));
189  /* Get rid of translation and rotation in Hessian */
190  Eigen::MatrixXd switchedTransformedHessian = internalCoordinates.hessianToInternal(switchedHessian);
191  tmp->optimizer.invH = switchedTransformedHessian.inverse();
192  }
193  else if (coordinateSystem == Utils::CoordinateSystem::Cartesian) {
194  /* no not give invH, because this leads to a lot of translations */
195  tmp->optimizer.guessVector =
196  std::make_shared<Eigen::VectorXd>(Eigen::Map<const Eigen::VectorXd>(mode.data(), mode.cols() * mode.rows()));
197  }
198  else {
199  throw std::logic_error("Unknown coordinate system " +
200  Utils::CoordinateSystemInterpreter::getStringFromCoordinateSystem(coordinateSystem));
201  }
202  }
203  /* Read vector from file to give as guess vector */
204  else if (!guessVectorFileName.empty()) {
205  /* push data into std::vector */
206  std::vector<double> container;
207  double num = 0.0;
208  std::ifstream file(guessVectorFileName, std::ios::in);
209  if (!file.is_open()) {
210  throw std::runtime_error("Problem when opening file " + guessVectorFileName);
211  }
212  while (file >> num) {
213  container.push_back(num);
214  }
215  file.close();
216  Eigen::Map<Eigen::VectorXd> guessVector(container.data(), container.size());
217  auto system = calc->getStructure();
218  /* Provided vector is Cartesian and internal shall be used -> transform to internal */
219  if (coordinateSystem != Utils::CoordinateSystem::Cartesian &&
220  static_cast<unsigned long>(guessVector.size()) == 3 * system->getElements().size()) {
221  std::shared_ptr<Utils::InternalCoordinates> transformation = nullptr;
222  if (coordinateSystem == Utils::CoordinateSystem::Internal) {
223  transformation = std::make_shared<Utils::InternalCoordinates>(*system);
224  }
225  else if (coordinateSystem == Utils::CoordinateSystem::CartesianWithoutRotTrans) {
226  transformation = std::make_shared<Utils::InternalCoordinates>(*system, true);
227  }
228  // vector to PositionCollection
229  Utils::PositionCollection guessPosition =
230  Eigen::Map<Utils::PositionCollection>(guessVector.data(), guessVector.size() / 3, 3);
231  try {
232  tmp->optimizer.guessVector =
233  std::make_shared<Eigen::VectorXd>(transformation->coordinatesToInternal(guessPosition));
234  }
235  catch (const Utils::InternalCoordinatesException& e) {
236  _logger->warning << "Could not transform given Cartesian coordinates to Internal coordinates.\n"
237  << "Omitting guess and initializing dimer with random vector. " << Core::Log::endl;
238  }
239  }
240  /* Provided vector is not Cartesian and Cartesian shall be used -> transform to Cartesian */
241  else if (coordinateSystem == Utils::CoordinateSystem::Cartesian &&
242  static_cast<unsigned long>(guessVector.size()) != 3 * system->getElements().size()) {
243  try {
244  std::shared_ptr<Utils::InternalCoordinates> transformation = nullptr;
245  if (static_cast<unsigned long>(guessVector.size()) == 3 * system->getElements().size() - 6 ||
246  static_cast<unsigned long>(guessVector.size()) == 3 * system->getElements().size() - 5) {
247  // cartesianWithoutRotTrans because 5 or 6 missing
248  transformation = std::make_shared<Utils::InternalCoordinates>(*system, true);
249  }
250  else {
251  transformation = std::make_shared<Utils::InternalCoordinates>(*system);
252  }
253  /* single guess vector cannot get backtransformed directly to Cartesian coords */
254  /* first get current position and transform to internal */
255  auto cartPos0 = system->getPositions();
256  auto ircVec0 = transformation->coordinatesToInternal(cartPos0);
257  /* add guessvector */
258  Eigen::VectorXd ircVec1 = ircVec0 + guessVector;
259  /* transform this new vector back to Cartesian and diff of 2 Cartesian PositionCollection is guessvector */
260  auto cartPos1 = transformation->coordinatesToCartesian(ircVec1);
261  Utils::PositionCollection posDiff = cartPos1 - cartPos0;
262  tmp->optimizer.guessVector =
263  std::make_shared<Eigen::VectorXd>(Eigen::Map<Eigen::VectorXd>(posDiff.data(), posDiff.size()));
264  }
265  catch (const Utils::InternalCoordinatesException& e) {
266  _logger->warning << "Could not transform given non-Cartesian coordinates to Cartesian.\n"
267  << "Omitting guess and initializing dimer with random vector. " << Core::Log::endl;
268  }
269  }
270  /* Provided vector should fit with wanted coordinates, optimizer later checks the correct length again */
271  else {
272  tmp->optimizer.guessVector = std::make_shared<Eigen::VectorXd>(guessVector);
273  }
274  }
275  /* Read 2 structures from files to form guess vector from discrete difference of the two */
276  else if (!discreteGuessesFileNames.at(0).empty() && !discreteGuessesFileNames.at(1).empty()) {
277  Utils::AtomCollection system;
278  /* read first */
279  std::string filepath1 = discreteGuessesFileNames.at(0);
280  std::ifstream file1(filepath1, std::ios::in);
281  if (!file1.is_open()) {
282  throw std::runtime_error("Problem when opening file " + filepath1);
283  }
284  system = XyzHandler::read(file1);
285  auto position1 = system.getPositions();
286  file1.close();
287 
288  /* read second */
289  std::string filepath2 = discreteGuessesFileNames.at(1);
290  std::ifstream file2(filepath2, std::ios::in);
291  if (!file2.is_open()) {
292  throw std::runtime_error("Problem when opening file " + filepath2);
293  }
294  system = XyzHandler::read(file2);
295  auto position2 = system.getPositions();
296  file2.close();
297 
298  Utils::PositionCollection positionDiff = position2 - position1;
299  Eigen::Map<Eigen::VectorXd> mode(positionDiff.data(), positionDiff.size());
300  std::shared_ptr<Utils::InternalCoordinates> transformation = nullptr;
301  if (coordinateSystem == Utils::CoordinateSystem::Internal) {
302  transformation = std::make_shared<Utils::InternalCoordinates>(system);
303  }
304  else if (coordinateSystem == Utils::CoordinateSystem::CartesianWithoutRotTrans) {
305  transformation = std::make_shared<Utils::InternalCoordinates>(system, true);
306  }
307  if (transformation) {
308  tmp->optimizer.guessVector = std::make_shared<Eigen::VectorXd>(transformation->coordinatesToInternal(positionDiff));
309  }
310  else {
311  tmp->optimizer.guessVector = std::make_shared<Eigen::VectorXd>(mode);
312  }
313  }
314  optimizer = std::move(tmp);
315  } // end of Dimer handling
316  else {
317  throw std::runtime_error(
318  "Unknown Optimizer requested for TS optimization, available are: Bofill, EVF and Dimer!");
319  }
320 
321  // Read and delete special settings
322  bool stopOnError = stopOnErrorExtraction(taskSettings);
323  // Apply settings
324  auto settings = optimizer->getSettings();
325  settings.merge(taskSettings);
326  // set automatically selected mode, if requested
327  // for Dimer this happened already, because it is GradientBased and does not know anything about modes / Hessians
328  // do not do in testRun because calc would be empty
329  if (!relevantAtoms.empty() && !testRunOnly) {
330  // necessary because multiple strings are accepted for eigenvectorfollowing optimizer
331  if (optimizertype != "dimer" && optimizertype != "bofill") {
332  settings.modifyValue("ev_follow_mode", selectMode(*calc, relevantAtoms));
333  }
334  else if (optimizertype == "bofill") {
335  settings.modifyValue("bofill_follow_mode", selectMode(*calc, relevantAtoms));
336  }
337  }
338  if (!settings.valid()) {
339  settings.throwIncorrectSettings();
340  }
341  optimizer->setSettings(settings);
342  if (!testRunOnly && !Utils::GeometryOptimization::settingsMakeSense(*optimizer)) {
343  throw std::logic_error("The given calculator settings are too inaccurate for the given convergence criteria of "
344  "this optimization Task");
345  }
346 
347  // If no errors encountered until here, the basic settings should be alright
348  if (testRunOnly) {
349  return true;
350  }
351 
352  // Add observer
353  // Trajectory stream
355  const std::string& partialOutput = (!_output.empty() ? _output[0] : _input[0]);
356  boost::filesystem::path dir(partialOutput);
357  boost::filesystem::create_directory(dir);
358  boost::filesystem::path trjfile(partialOutput + ".tsopt.trj.xyz");
359  std::ofstream trajectory((dir / trjfile).string(), std::ofstream::out);
360  cout.printf("%7s %16s %16s %16s %16s\n", "Cycle", "Energy", "Energy Diff.", "Step RMS", "Max. Step");
361  double oldEnergy = 0.0;
362  Eigen::VectorXd oldParams;
363  auto func = [&](const int& cycle, const double& energy, const Eigen::VectorXd& params) {
364  if (oldParams.size() != params.size()) {
365  oldParams = params;
366  }
367  auto diff = (params - oldParams).eval();
368  cout.printf("%7d %+16.9f %+16.9f %+16.9f %+16.9f\n", cycle, energy, energy - oldEnergy,
369  sqrt(diff.squaredNorm() / diff.size()), diff.cwiseAbs().maxCoeff());
370  oldEnergy = energy;
371  oldParams = params;
372  auto structure = calc->getStructure();
373  XyzHandler::write(trajectory, *structure, std::to_string(energy));
374  };
375  optimizer->addObserver(func);
376 
377  // Run optimization
378  auto structure = calc->getStructure();
379  int cycles = 0;
380  try {
381  cycles = optimizer->optimize(*structure, *_logger);
382  }
383  catch (...) {
384  XyzHandler::write(trajectory, *calc->getStructure());
385  trajectory.close();
386  _logger->error << "TS Optimization failed with error!" << Core::Log::endl;
387  if (stopOnError) {
388  throw;
389  }
390  _logger->error << boost::current_exception_diagnostic_information() << Core::Log::endl;
391  return false;
392  }
393  trajectory.close();
394 
395  int maxiter = settings.getInt("convergence_max_iterations");
396  if (cycles < maxiter) {
397  cout << Core::Log::endl
398  << " Converged after " << cycles << " iterations." << Core::Log::endl
399  << Core::Log::endl;
400  }
401  else {
402  cout << Core::Log::endl
403  << " Stopped after " << maxiter << " iterations." << Core::Log::endl
404  << Core::Log::endl;
405  if (stopOnError) {
406  throw std::runtime_error("Problem: TS optimization did not converge.");
407  }
408  }
409 
410  // Print/Store results
411  systems[partialOutput] = calc;
412  boost::filesystem::path xyzfile(partialOutput + ".xyz");
413  std::ofstream xyz((dir / xyzfile).string(), std::ofstream::out);
414  XyzHandler::write(xyz, *(calc->getStructure()));
415  xyz.close();
416 
417  return cycles < maxiter;
418  }
419 
420  private:
421  int selectMode(Core::Calculator& calc, const std::vector<int>& relevantAtoms) const {
422  _logger->output << " Automatically selecting normal mode " << Core::Log::endl;
423  int nAtoms = calc.getStructure()->size();
424  for (const auto& index : relevantAtoms) {
425  if (index >= nAtoms) {
426  throw std::logic_error("You gave an atom index larger than the number of atoms in automatic_mode_selection. "
427  "This does not make sense.");
428  }
429  }
430  Utils::NormalModesContainer modes = getNormalModes(calc);
431  auto waveNumbers = modes.getWaveNumbers();
432  std::vector<double> contributions;
433  // cycle all imag. frequencies
434  for (int i = 0; i < static_cast<int>(waveNumbers.size()); ++i) {
435  if (waveNumbers[i] >= 0.0) {
436  if (i == 0) {
437  throw std::runtime_error("Structure has no imaginary normal mode.");
438  }
439  break;
440  }
441  auto mode = modes.getMode(i);
442  double contribution = 0.0;
443  // cycle all atoms
444  for (int j = 0; j < mode.size(); ++j) {
445  if (std::find(relevantAtoms.begin(), relevantAtoms.end(), j) != relevantAtoms.end()) {
446  contribution += mode.row(j).norm(); // add norm of vector of relevant atom contributing to mode
447  }
448  }
449  contributions.push_back(contribution);
450  }
451  int selection = std::distance(contributions.begin(), std::max_element(contributions.begin(), contributions.end()));
452  double tolerance = contributions[selection] * 0.1;
453  // Contributions are sorted by wave number
454  // -> pick the first one that fits into the tolerance.
455  for (unsigned int i = 0; i < contributions.size(); i++) {
456  if ((contributions[selection] - contributions[i]) < tolerance) {
457  selection = i;
458  break;
459  }
460  }
461  _logger->output << " Automatically selected mode " << std::to_string(selection) << " with wave number "
462  << std::to_string(waveNumbers[selection]) << " cm^-1" << Core::Log::endl;
463  return selection;
464  }
465 
466  static Utils::NormalModesContainer getNormalModes(Core::Calculator& calc) {
467  if (!calc.results().has<Utils::Property::Hessian>()) {
468  calc.setRequiredProperties(Utils::Property::Hessian);
469  calc.calculate("Hessian Calculation");
470  }
471  auto hessian = calc.results().get<Utils::Property::Hessian>();
472  auto system = calc.getStructure();
473  return Utils::NormalModeAnalysis::calculateNormalModes(hessian, system->getElements(), system->getPositions());
474  }
475 };
476 
477 } // namespace Readuct
478 } // namespace Scine
479 
480 #endif // READUCT_TSOPTIMIZATIONTASK_H_
Eigen::MatrixXd hessianToInternal(const HessianMatrix &cartesian) const
void warningIfMultipleOutputsGiven() const
Warn if more than one output system was specified.
Definition: Task.h:99
virtual std::unique_ptr< Utils::AtomCollection > getStructure() const =0
static std::ostream & endl(std::ostream &os)
const PositionCollection & getPositions() const
const std::vector< std::string > & input() const
Getter for the expected names of the input systems.
Definition: Task.h:77
Eigen::VectorXd coordinatesToInternal(const PositionCollection &cartesian) const
std::string name() const override
Getter for the tasks name.
Definition: TsOptimizationTask.h:47
TsOptimizationTask(std::vector< std::string > input, std::vector< std::string > output, std::shared_ptr< Core::Log > logger=nullptr)
Construct a new TsOptimizationTask.
Definition: TsOptimizationTask.h:43
bool run(SystemsMap &systems, Utils::UniversalSettings::ValueCollection taskSettings, bool testRunOnly=false) const final
Executes the actual task represented by this class.
Definition: TsOptimizationTask.h:51
void warningIfMultipleInputsGiven() const
Warn if more than one input system was specified.
Definition: Task.h:90
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
Definition: TsOptimizationTask.h:35