7 #ifndef READUCT_AFIROPTIMIZATIONTASK_H_
8 #define READUCT_AFIROPTIMIZATIONTASK_H_
21 #include <boost/exception/diagnostic_information.hpp>
22 #include <boost/filesystem.hpp>
40 :
Task(std::move(input), std::move(output), std::move(logger)) {
43 std::string
name()
const override {
44 return "AFIR Optimization";
49 observers = {})
const final {
53 bool silentCalculator = taskSettings.extract(
"silent_stdout_calculator",
true);
54 std::shared_ptr<Core::Calculator> calc;
57 calc = copyCalculator(systems, _input.front(),
name());
58 Utils::CalculationRoutines::setLog(*calc,
true,
true, !silentCalculator);
61 if (calc->getStructure()->size() == 1) {
62 throw std::runtime_error(
"Cannot calculate AFIR for monoatomic systems.");
67 auto optimizertype = taskSettings.extract(
"optimizer", std::string{
"BFGS"});
68 std::transform(optimizertype.begin(), optimizertype.end(), optimizertype.begin(), ::toupper);
69 std::shared_ptr<Utils::AfirOptimizerBase> optimizer;
70 if (optimizertype ==
"LBFGS") {
71 auto tmp = std::make_shared<Utils::AfirOptimizer<Utils::Lbfgs>>(*calc);
72 tmp->optimizer.useTrustRadius =
true;
73 tmp->optimizer.trustRadius = 0.1;
74 optimizer = std::move(tmp);
76 else if (optimizertype ==
"BFGS") {
77 auto tmp = std::make_shared<Utils::AfirOptimizer<Utils::Bfgs>>(*calc);
78 tmp->optimizer.useTrustRadius =
true;
79 tmp->optimizer.trustRadius = 0.1;
80 optimizer = std::move(tmp);
82 else if (optimizertype ==
"SD" || optimizertype ==
"STEEPESTDESCENT") {
83 auto tmp = std::make_shared<Utils::AfirOptimizer<Utils::SteepestDescent>>(*calc);
84 optimizer = std::move(tmp);
87 throw std::runtime_error(
"Unknown Optimizer requested for AFIR, available are: SD, BFGS and LBFGS!");
91 bool stopOnError = stopOnErrorExtraction(taskSettings);
94 auto settings = optimizer->getSettings();
95 settings.merge(taskSettings);
96 if (!settings.valid()) {
97 settings.throwIncorrectSettings();
99 optimizer->setSettings(settings);
100 if (!testRunOnly && !Utils::GeometryOptimization::settingsMakeSense(*optimizer)) {
101 throw std::logic_error(
"The given calculator settings are too inaccurate for the given convergence criteria of "
102 "this optimization Task");
112 using Writer = Utils::XyzStreamHandler;
113 const std::string& outputSystem = (!_output.empty() ? _output.front() : _input.front());
114 boost::filesystem::path dir(outputSystem);
115 boost::filesystem::create_directory(dir);
116 boost::filesystem::path trjfile(outputSystem +
".afir.trj.xyz");
117 std::ofstream trajectory((dir / trjfile).
string(), std::ofstream::out);
118 double oldEnergy = 0.0;
119 Eigen::VectorXd oldParams;
120 auto cout = _logger->output;
121 auto func = [&](
const int& cycle,
const double& energy,
const Eigen::VectorXd& params) {
122 if (oldParams.size() != params.size()) {
126 cout.printf(
"%7s %16s %16s %16s %16s\n",
"Cycle",
"Energy",
"Energy Diff.",
"Step RMS",
"Max. Step");
128 auto diff = (params - oldParams).eval();
129 cout.printf(
"%7d %+16.9f %+16.9f %+16.9f %+16.9f\n", cycle, energy, energy - oldEnergy,
130 sqrt(diff.squaredNorm() / diff.size()), diff.cwiseAbs().maxCoeff());
133 auto structure = calc->getStructure();
134 Writer::write(trajectory, *structure, std::to_string(energy));
136 optimizer->addObserver(func);
139 auto customObservers = [&calc, &observers](
const int& cycle,
const double& ,
const Eigen::VectorXd& ) {
140 for (
auto& observer : observers) {
141 auto atoms = calc->getStructure();
142 Utils::Results& results = calc->results();
143 observer(cycle, *atoms, results,
"afir_scan");
146 optimizer->addObserver(customObservers);
149 auto structure = calc->getStructure();
152 cycles = optimizer->optimize(*structure, *_logger);
155 Writer::write(trajectory, *calc->getStructure());
157 _logger->error <<
"AFIR Optimization failed with error!" <<
Core::Log::endl;
161 _logger->error << boost::current_exception_diagnostic_information() <<
Core::Log::endl;
166 int maxiter = settings.getInt(
"convergence_max_iterations");
167 if (maxiter > cycles) {
177 throw std::runtime_error(
"Problem: AFIR optimization did not converge.");
182 boost::filesystem::path xyzfile(outputSystem +
".xyz");
183 std::ofstream xyz((dir / xyzfile).
string(), std::ofstream::out);
184 Writer::write(xyz, *(calc->getStructure()));
186 systems[outputSystem] = std::move(calc);
188 return cycles < maxiter;
195 #endif // READUCT_AFIROPTIMIZATIONTASK_H_
std::string name() const override
Getter for the tasks name.
Definition: AfirOptimizationTask.h:43
void warningIfMultipleOutputsGiven() const
Warn if more than one output system was specified.
Definition: Task.h:107
AfirOptimizationTask(std::vector< std::string > input, std::vector< std::string > output, std::shared_ptr< Core::Log > logger=nullptr)
Construct a new AfirOptimizationTask.
Definition: AfirOptimizationTask.h:39
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:85
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: AfirOptimizationTask.h:47
Definition: AfirOptimizationTask.h:31
void warningIfMultipleInputsGiven() const
Warn if more than one input system was specified.
Definition: Task.h:98
const std::vector< std::string > & output() const
Getter for the names of the output systems generated by this task.
Definition: Task.h:92
The base class for all tasks in Readuct.
Definition: Task.h:29