Scine::Readuct  4.0.0
This is the SCINE module Readuct.
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Pages
GeometryOptimizationTask.h
Go to the documentation of this file.
1 
7 #ifndef READUCT_GEOMETRYOPTIMIZATIONTASK_H_
8 #define READUCT_GEOMETRYOPTIMIZATIONTASK_H_
9 
10 /* Readuct */
11 #include "Tasks/Task.h"
12 /* Scine */
23 /* External includes */
24 #include <boost/exception/diagnostic_information.hpp>
25 #include <boost/filesystem.hpp>
26 #include <cstdio>
27 #include <fstream>
28 
29 namespace Scine {
30 namespace Readuct {
31 
33  public:
40  GeometryOptimizationTask(std::vector<std::string> input, std::vector<std::string> output,
41  std::shared_ptr<Core::Log> logger = nullptr)
42  : Task(std::move(input), std::move(output), std::move(logger)) {
43  }
44 
45  std::string name() const override {
46  return "Geometry Optimization";
47  }
48 
49  bool run(SystemsMap& systems, Utils::UniversalSettings::ValueCollection taskSettings, bool testRunOnly = false) const final {
52  bool silentCalculator = taskSettings.extract("silent_stdout_calculator", true);
53  // Get/Copy Calculator
54  std::shared_ptr<Core::Calculator> calc;
55  if (!testRunOnly) { // leave out in case of task chaining --> attention calc is NULL
56  // Note: _input is guaranteed not to be empty by Task constructor
57  calc = copyCalculator(systems, _input.front(), name());
58  Utils::CalculationRoutines::setLog(*calc, true, true, !silentCalculator);
59  }
60 
61  // Generate optimizer
62  auto optimizertype = taskSettings.extract("optimizer", std::string{"BFGS"});
63  auto unitcelloptimizertype = taskSettings.extract("unitcelloptimizer", std::string{""});
64  auto optimizer = constructOptimizer(*calc, optimizertype, unitcelloptimizertype);
65 
66  // Read and delete special settings
67  bool stopOnError = stopOnErrorExtraction(taskSettings);
68  // Apply settings
69  auto settings = optimizer->getSettings();
70  settings.merge(taskSettings);
71  if (!testRunOnly && systems.at(_input[0])->getStructure()->size() == 1 &&
72  Utils::CoordinateSystemInterpreter::getCoordinateSystemFromString(
73  settings.getString("geoopt_coordinate_system")) != Utils::CoordinateSystem::Cartesian) {
74  /* For monoatomic systems avoid transformation to internals resulting into
75  * empty optimization parameters. Necessary because observer and
76  * convergence check use transformed parameters
77  */
78  _logger->warning << " Cannot treat a monoatomic system in internal coordinates. Switching to Cartesians."
80  << Core::Log::endl;
81  settings.modifyString("geoopt_coordinate_system", Utils::CoordinateSystemInterpreter::getStringFromCoordinateSystem(
82  Utils::CoordinateSystem::Cartesian));
83  }
84  if (!settings.valid()) {
85  settings.throwIncorrectSettings();
86  }
87  optimizer->setSettings(settings);
88  if (!testRunOnly && !Utils::GeometryOptimization::settingsMakeSense(*optimizer)) {
89  throw std::logic_error("The given calculator settings are too inaccurate for the given convergence criteria of "
90  "this optimization Task");
91  }
92 
93  // If no errors encountered until here, the basic settings should be alright
94  if (testRunOnly) {
95  return true;
96  }
97 
98  // Add observer
99  // Trajectory stream
100  const std::string& outputSystem = (!_output.empty() ? _output.front() : _input.front());
101  using Writer = Utils::XyzStreamHandler;
102  boost::filesystem::path dir(outputSystem);
103  boost::filesystem::create_directory(dir);
104  boost::filesystem::path trjfile(outputSystem + ".opt.trj.xyz");
105  std::ofstream trajectory((dir / trjfile).string(), std::ofstream::out);
106  double oldEnergy = 0.0;
107  Eigen::VectorXd oldParams;
108  auto cout = _logger->output;
109  auto func = [&](const int& cycle, const double& energy, const Eigen::VectorXd& params) {
110  if (oldParams.size() != params.size()) {
111  oldParams = params;
112  }
113  if (cycle == 1) {
114  cout.printf("%7s %16s %16s %16s %16s\n", "Cycle", "Energy", "Energy Diff.", "Step RMS", "Max. Step");
115  }
116  auto diff = (params - oldParams).eval();
117  cout.printf("%7d %+16.9f %+16.9f %+16.9f %+16.9f\n", cycle, energy, energy - oldEnergy,
118  sqrt(diff.squaredNorm() / diff.size()), diff.cwiseAbs().maxCoeff());
119  oldEnergy = energy;
120  oldParams = params;
121  auto structure = calc->getStructure();
122  Writer::write(trajectory, *structure, std::to_string(energy));
123  };
124  optimizer->addObserver(func);
125 
126  // Run optimization
127  auto structure = calc->getStructure();
128  int cycles = 0;
129  try {
130  cycles = optimizer->optimize(*structure, *_logger);
131  }
132  catch (...) {
133  Writer::write(trajectory, *calc->getStructure());
134  trajectory.close();
135  _logger->error << "Optimization failed with error!" << Core::Log::endl;
136  if (stopOnError) {
137  throw;
138  }
139  _logger->error << boost::current_exception_diagnostic_information() << Core::Log::endl;
140  return false;
141  }
142  trajectory.close();
143 
144  int maxiter = settings.getInt("convergence_max_iterations");
145  if (cycles < maxiter) {
146  cout << Core::Log::endl
147  << " Converged after " << cycles << " iterations." << Core::Log::endl
148  << Core::Log::endl;
149  }
150  else {
151  cout << Core::Log::endl
152  << " Stopped after " << maxiter << " iterations." << Core::Log::endl
153  << Core::Log::endl;
154  if (stopOnError) {
155  throw std::runtime_error("Problem: Structure optimization did not converge.");
156  }
157  }
158 
159  // Print/Store result
160  systems[outputSystem] = calc;
161  boost::filesystem::path xyzfile(outputSystem + ".xyz");
162  std::ofstream xyz((dir / xyzfile).string(), std::ofstream::out);
163  Writer::write(xyz, *(calc->getStructure()));
164  xyz.close();
165 
166  return cycles < maxiter;
167  }
168 
169  inline std::shared_ptr<Utils::GeometryOptimizerBase> constructOptimizer(Core::Calculator& calc, std::string type,
170  std::string cellType) const {
171  std::transform(type.begin(), type.end(), type.begin(), ::toupper);
172  std::transform(cellType.begin(), cellType.end(), cellType.begin(), ::toupper);
173  if (type == "LBFGS") {
174  if (cellType.empty()) {
175  return std::make_shared<Utils::GeometryOptimizer<Utils::Lbfgs>>(calc);
176  }
177  else {
178  if (cellType == "LBFGS") {
179  return std::make_shared<Utils::UnitCellGeometryOptimizer<Utils::Lbfgs, Utils::Lbfgs>>(calc);
180  }
181  else if (cellType == "BFGS") {
182  return std::make_shared<Utils::UnitCellGeometryOptimizer<Utils::Lbfgs, Utils::Bfgs>>(calc);
183  }
184  else if (cellType == "SD" || cellType == "STEEPESTDESCENT") {
185  return std::make_shared<Utils::UnitCellGeometryOptimizer<Utils::Lbfgs, Utils::SteepestDescent>>(calc);
186  }
187  throw std::runtime_error(
188  "Unknown CellOptimizer requested for a geometry optimization, available are: SD, BFGS and LBFGS!");
189  }
190  }
191  else if (type == "BFGS") {
192  if (cellType.empty()) {
193  return std::make_shared<Utils::GeometryOptimizer<Utils::Bfgs>>(calc);
194  }
195  else {
196  if (cellType == "LBFGS") {
197  return std::make_shared<Utils::UnitCellGeometryOptimizer<Utils::Bfgs, Utils::Lbfgs>>(calc);
198  }
199  else if (cellType == "BFGS") {
200  return std::make_shared<Utils::UnitCellGeometryOptimizer<Utils::Bfgs, Utils::Bfgs>>(calc);
201  }
202  else if (cellType == "SD" || cellType == "STEEPESTDESCENT") {
203  return std::make_shared<Utils::UnitCellGeometryOptimizer<Utils::Bfgs, Utils::SteepestDescent>>(calc);
204  }
205  throw std::runtime_error(
206  "Unknown CellOptimizer requested for a geometry optimization, available are: SD, BFGS and LBFGS!");
207  }
208  }
209  else if (type == "SD" || type == "STEEPESTDESCENT") {
210  if (cellType.empty()) {
211  return std::make_shared<Utils::GeometryOptimizer<Utils::SteepestDescent>>(calc);
212  }
213  else {
214  if (cellType == "LBFGS") {
215  return std::make_shared<Utils::UnitCellGeometryOptimizer<Utils::SteepestDescent, Utils::Lbfgs>>(calc);
216  }
217  else if (cellType == "BFGS") {
218  return std::make_shared<Utils::UnitCellGeometryOptimizer<Utils::SteepestDescent, Utils::Bfgs>>(calc);
219  }
220  else if (cellType == "SD" || cellType == "STEEPESTDESCENT") {
221  return std::make_shared<Utils::UnitCellGeometryOptimizer<Utils::SteepestDescent, Utils::SteepestDescent>>(calc);
222  }
223  throw std::runtime_error(
224  "Unknown CellOptimizer requested for a geometry optimization, available are: SD, BFGS and LBFGS!");
225  }
226  }
227  else if (type == "NR" || type == "NEWTONRAPHSON") {
228  if (cellType.empty()) {
229  return std::make_shared<Utils::GeometryOptimizer<Utils::NewtonRaphson>>(calc);
230  }
231  else {
232  if (cellType == "LBFGS") {
233  return std::make_shared<Utils::UnitCellGeometryOptimizer<Utils::NewtonRaphson, Utils::Lbfgs>>(calc);
234  }
235  else if (cellType == "BFGS") {
236  return std::make_shared<Utils::UnitCellGeometryOptimizer<Utils::NewtonRaphson, Utils::Bfgs>>(calc);
237  }
238  else if (cellType == "SD" || cellType == "STEEPESTDESCENT") {
239  return std::make_shared<Utils::UnitCellGeometryOptimizer<Utils::NewtonRaphson, Utils::SteepestDescent>>(calc);
240  }
241  throw std::runtime_error(
242  "Unknown CellOptimizer requested for a geometry optimization, available are: SD, BFGS and LBFGS!");
243  }
244  }
245  throw std::runtime_error(
246  "Unknown Optimizer requested for a geometry optimization, available are: SD, NR, BFGS and LBFGS!");
247  };
248 };
249 
250 } // namespace Readuct
251 } // namespace Scine
252 
253 #endif // READUCT_GEOMETRYOPTIMIZATIONTASK_H_
void warningIfMultipleOutputsGiven() const
Warn if more than one output system was specified.
Definition: Task.h:99
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
std::string name() const override
Getter for the tasks name.
Definition: GeometryOptimizationTask.h:45
bool run(SystemsMap &systems, Utils::UniversalSettings::ValueCollection taskSettings, bool testRunOnly=false) const final
Executes the actual task represented by this class.
Definition: GeometryOptimizationTask.h:49
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
GeometryOptimizationTask(std::vector< std::string > input, std::vector< std::string > output, std::shared_ptr< Core::Log > logger=nullptr)
Construct a new GeometryOptimizationTask.
Definition: GeometryOptimizationTask.h:40
The base class for all tasks in Readuct.
Definition: Task.h:28
Definition: GeometryOptimizationTask.h:32