7 #ifndef READUCT_GEOMETRYOPTIMIZATIONTASK_H_
8 #define READUCT_GEOMETRYOPTIMIZATIONTASK_H_
24 #include <boost/exception/diagnostic_information.hpp>
25 #include <boost/filesystem.hpp>
41 std::shared_ptr<Core::Log> logger =
nullptr)
42 :
Task(std::move(input), std::move(output), std::move(logger)) {
45 std::string
name()
const override {
46 return "Geometry Optimization";
51 observers = {})
const final {
54 bool silentCalculator = taskSettings.extract(
"silent_stdout_calculator",
true);
56 std::shared_ptr<Core::Calculator> calc;
59 calc = copyCalculator(systems, _input.front(),
name());
60 Utils::CalculationRoutines::setLog(*calc,
true,
true, !silentCalculator);
64 auto optimizertype = taskSettings.extract(
"optimizer", std::string{
"BFGS"});
65 auto unitcelloptimizertype = taskSettings.extract(
"unitcelloptimizer", std::string{
""});
66 auto optimizer = constructOptimizer(*calc, optimizertype, unitcelloptimizertype);
69 bool stopOnError = stopOnErrorExtraction(taskSettings);
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) {
80 _logger->warning <<
" Cannot treat a monoatomic system in internal coordinates. Switching to Cartesians."
83 settings.modifyString(
"geoopt_coordinate_system", Utils::CoordinateSystemInterpreter::getStringFromCoordinateSystem(
84 Utils::CoordinateSystem::Cartesian));
86 if (!settings.valid()) {
87 settings.throwIncorrectSettings();
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");
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()) {
116 cout.printf(
"%7s %16s %16s %16s %16s\n",
"Cycle",
"Energy",
"Energy Diff.",
"Step RMS",
"Max. Step");
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());
123 auto structure = calc->getStructure();
124 Writer::write(trajectory, *structure, std::to_string(energy));
126 optimizer->addObserver(func);
129 auto customObservers = [&calc, &observers](
const int& cycle,
const double& ,
const Eigen::VectorXd& ) {
130 for (
auto& observer : observers) {
131 auto atoms = calc->getStructure();
132 Utils::Results& results = calc->results();
133 observer(cycle, *atoms, results,
"geometry_optimization");
136 optimizer->addObserver(customObservers);
139 auto structure = calc->getStructure();
142 cycles = optimizer->optimize(*structure, *_logger);
145 Writer::write(trajectory, *calc->getStructure());
147 _logger->error <<
"Optimization failed with error!" <<
Core::Log::endl;
151 _logger->error << boost::current_exception_diagnostic_information() <<
Core::Log::endl;
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
163 cout << Core::Log::endl
164 <<
" Stopped after " << maxiter <<
" iterations." << Core::Log::endl
167 throw std::runtime_error(
"Problem: Structure optimization did not converge.");
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()));
178 return cycles < maxiter;
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);
190 if (cellType ==
"LBFGS") {
191 return std::make_shared<Utils::UnitCellGeometryOptimizer<Utils::Lbfgs, Utils::Lbfgs>>(calc);
193 else if (cellType ==
"BFGS") {
194 return std::make_shared<Utils::UnitCellGeometryOptimizer<Utils::Lbfgs, Utils::Bfgs>>(calc);
196 else if (cellType ==
"SD" || cellType ==
"STEEPESTDESCENT") {
197 return std::make_shared<Utils::UnitCellGeometryOptimizer<Utils::Lbfgs, Utils::SteepestDescent>>(calc);
199 throw std::runtime_error(
200 "Unknown CellOptimizer requested for a geometry optimization, available are: SD, BFGS and LBFGS!");
203 else if (type ==
"BFGS") {
204 if (cellType.empty()) {
205 return std::make_shared<Utils::GeometryOptimizer<Utils::Bfgs>>(calc);
208 if (cellType ==
"LBFGS") {
209 return std::make_shared<Utils::UnitCellGeometryOptimizer<Utils::Bfgs, Utils::Lbfgs>>(calc);
211 else if (cellType ==
"BFGS") {
212 return std::make_shared<Utils::UnitCellGeometryOptimizer<Utils::Bfgs, Utils::Bfgs>>(calc);
214 else if (cellType ==
"SD" || cellType ==
"STEEPESTDESCENT") {
215 return std::make_shared<Utils::UnitCellGeometryOptimizer<Utils::Bfgs, Utils::SteepestDescent>>(calc);
217 throw std::runtime_error(
218 "Unknown CellOptimizer requested for a geometry optimization, available are: SD, BFGS and LBFGS!");
221 else if (type ==
"SD" || type ==
"STEEPESTDESCENT") {
222 if (cellType.empty()) {
223 return std::make_shared<Utils::GeometryOptimizer<Utils::SteepestDescent>>(calc);
226 if (cellType ==
"LBFGS") {
227 return std::make_shared<Utils::UnitCellGeometryOptimizer<Utils::SteepestDescent, Utils::Lbfgs>>(calc);
229 else if (cellType ==
"BFGS") {
230 return std::make_shared<Utils::UnitCellGeometryOptimizer<Utils::SteepestDescent, Utils::Bfgs>>(calc);
232 else if (cellType ==
"SD" || cellType ==
"STEEPESTDESCENT") {
233 return std::make_shared<Utils::UnitCellGeometryOptimizer<Utils::SteepestDescent, Utils::SteepestDescent>>(calc);
235 throw std::runtime_error(
236 "Unknown CellOptimizer requested for a geometry optimization, available are: SD, BFGS and LBFGS!");
239 else if (type ==
"NR" || type ==
"NEWTONRAPHSON") {
240 if (cellType.empty()) {
241 return std::make_shared<Utils::GeometryOptimizer<Utils::NewtonRaphson>>(calc);
244 if (cellType ==
"LBFGS") {
245 return std::make_shared<Utils::UnitCellGeometryOptimizer<Utils::NewtonRaphson, Utils::Lbfgs>>(calc);
247 else if (cellType ==
"BFGS") {
248 return std::make_shared<Utils::UnitCellGeometryOptimizer<Utils::NewtonRaphson, Utils::Bfgs>>(calc);
250 else if (cellType ==
"SD" || cellType ==
"STEEPESTDESCENT") {
251 return std::make_shared<Utils::UnitCellGeometryOptimizer<Utils::NewtonRaphson, Utils::SteepestDescent>>(calc);
253 throw std::runtime_error(
254 "Unknown CellOptimizer requested for a geometry optimization, available are: SD, BFGS and LBFGS!");
257 throw std::runtime_error(
258 "Unknown Optimizer requested for a geometry optimization, available are: SD, NR, BFGS and LBFGS!");
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