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