7 #ifndef READUCT_TSOPTIMIZATIONTASK_H_
8 #define READUCT_TSOPTIMIZATIONTASK_H_
24 #include <boost/exception/diagnostic_information.hpp>
25 #include <boost/filesystem.hpp>
44 :
Task(std::move(input), std::move(output), std::move(logger)) {
47 std::string
name()
const override {
48 return "TS Optimization";
55 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);
63 if (calc->getStructure()->size() == 1) {
64 throw std::runtime_error(
"Cannot calculate transition state for monoatomic systems.");
69 std::string optimizertype = taskSettings.extract(
"optimizer", std::string{
"bofill"});
70 std::transform(optimizertype.begin(), optimizertype.end(), optimizertype.begin(), ::tolower);
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.");
79 std::vector<int> relevantAtoms = taskSettings.extract(
"automatic_mode_selection", std::vector<int>{});
81 if (!relevantAtoms.empty()) {
82 for (
const auto& index : relevantAtoms) {
84 throw std::logic_error(
"You gave an atom index smaller than 0 in automatic_mode_selection. "
85 "This does not make sense.");
91 auto cout = _logger->output;
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);
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);
103 else if (optimizertype ==
"dimer") {
104 auto tmp = std::make_shared<Utils::GeometryOptimizer<Utils::Dimer>>(*calc);
107 auto coordinateSystem = tmp->coordinateSystem;
108 if (taskSettings.valueExists(
"geoopt_coordinate_system")) {
109 coordinateSystem = Utils::CoordinateSystemInterpreter::getCoordinateSystemFromString(
110 taskSettings.getString(
"geoopt_coordinate_system"));
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 "
120 if (taskSettings.valueExists(
"dimer_follow_mode") || !relevantAtoms.empty()) {
121 useEigenvectorForFirstStep =
true;
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>{
"",
""});
129 if (discreteGuessesFileNames.size() != 2) {
130 throw std::runtime_error(
"Problem with discrete guesses. You may only give list with two filepaths.");
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"
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. "
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. "
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. "
161 else if (useEigenvectorForFirstStep) {
162 tmp->optimizer.skipFirstRotation =
true;
163 auto modes = getNormalModes(*calc);
164 auto system = calc->getStructure();
165 if (!relevantAtoms.empty()) {
166 followedMode = selectMode(*calc, relevantAtoms);
168 auto mode = modes.getMode(followedMode);
170 if (coordinateSystem == Utils::CoordinateSystem::Internal) {
173 tmp->optimizer.guessVector = std::make_shared<Eigen::VectorXd>(internalCoordinates.
coordinatesToInternal(mode));
175 else if (coordinateSystem == Utils::CoordinateSystem::CartesianWithoutRotTrans) {
178 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> es;
179 auto hessian = calc->results().get<Utils::Property::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();
188 tmp->optimizer.guessVector = std::make_shared<Eigen::VectorXd>(internalCoordinates.
coordinatesToInternal(mode));
190 Eigen::MatrixXd switchedTransformedHessian = internalCoordinates.
hessianToInternal(switchedHessian);
191 tmp->optimizer.invH = switchedTransformedHessian.inverse();
193 else if (coordinateSystem == Utils::CoordinateSystem::Cartesian) {
195 tmp->optimizer.guessVector =
196 std::make_shared<Eigen::VectorXd>(Eigen::Map<const Eigen::VectorXd>(mode.data(), mode.cols() * mode.rows()));
199 throw std::logic_error(
"Unknown coordinate system " +
200 Utils::CoordinateSystemInterpreter::getStringFromCoordinateSystem(coordinateSystem));
204 else if (!guessVectorFileName.empty()) {
206 std::vector<double> container;
208 std::ifstream file(guessVectorFileName, std::ios::in);
209 if (!file.is_open()) {
210 throw std::runtime_error(
"Problem when opening file " + guessVectorFileName);
212 while (file >> num) {
213 container.push_back(num);
216 Eigen::Map<Eigen::VectorXd> guessVector(container.data(), container.size());
217 auto system = calc->getStructure();
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);
225 else if (coordinateSystem == Utils::CoordinateSystem::CartesianWithoutRotTrans) {
226 transformation = std::make_shared<Utils::InternalCoordinates>(*system,
true);
230 Eigen::Map<Utils::PositionCollection>(guessVector.data(), guessVector.size() / 3, 3);
232 tmp->optimizer.guessVector =
233 std::make_shared<Eigen::VectorXd>(transformation->coordinatesToInternal(guessPosition));
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;
241 else if (coordinateSystem == Utils::CoordinateSystem::Cartesian &&
242 static_cast<unsigned long>(guessVector.size()) != 3 * system->getElements().size()) {
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) {
248 transformation = std::make_shared<Utils::InternalCoordinates>(*system,
true);
251 transformation = std::make_shared<Utils::InternalCoordinates>(*system);
255 auto cartPos0 = system->getPositions();
256 auto ircVec0 = transformation->coordinatesToInternal(cartPos0);
258 Eigen::VectorXd ircVec1 = ircVec0 + guessVector;
260 auto cartPos1 = transformation->coordinatesToCartesian(ircVec1);
262 tmp->optimizer.guessVector =
263 std::make_shared<Eigen::VectorXd>(Eigen::Map<Eigen::VectorXd>(posDiff.data(), posDiff.size()));
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;
272 tmp->optimizer.guessVector = std::make_shared<Eigen::VectorXd>(guessVector);
276 else if (!discreteGuessesFileNames.at(0).empty() && !discreteGuessesFileNames.at(1).empty()) {
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);
284 system = XyzHandler::read(file1);
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);
294 system = XyzHandler::read(file2);
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);
304 else if (coordinateSystem == Utils::CoordinateSystem::CartesianWithoutRotTrans) {
305 transformation = std::make_shared<Utils::InternalCoordinates>(system,
true);
307 if (transformation) {
308 tmp->optimizer.guessVector = std::make_shared<Eigen::VectorXd>(transformation->coordinatesToInternal(positionDiff));
311 tmp->optimizer.guessVector = std::make_shared<Eigen::VectorXd>(mode);
314 optimizer = std::move(tmp);
317 throw std::runtime_error(
318 "Unknown Optimizer requested for TS optimization, available are: Bofill, EVF and Dimer!");
322 bool stopOnError = stopOnErrorExtraction(taskSettings);
324 auto settings = optimizer->getSettings();
325 settings.merge(taskSettings);
329 if (!relevantAtoms.empty() && !testRunOnly) {
331 if (optimizertype !=
"dimer" && optimizertype !=
"bofill") {
332 settings.modifyValue(
"ev_follow_mode", selectMode(*calc, relevantAtoms));
334 else if (optimizertype ==
"bofill") {
335 settings.modifyValue(
"bofill_follow_mode", selectMode(*calc, relevantAtoms));
338 if (!settings.valid()) {
339 settings.throwIncorrectSettings();
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");
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()) {
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());
372 auto structure = calc->getStructure();
373 XyzHandler::write(trajectory, *structure, std::to_string(energy));
375 optimizer->addObserver(func);
378 auto structure = calc->getStructure();
381 cycles = optimizer->optimize(*structure, *_logger);
384 XyzHandler::write(trajectory, *calc->getStructure());
386 _logger->error <<
"TS Optimization failed with error!" <<
Core::Log::endl;
390 _logger->error << boost::current_exception_diagnostic_information() <<
Core::Log::endl;
395 int maxiter = settings.getInt(
"convergence_max_iterations");
396 if (cycles < maxiter) {
406 throw std::runtime_error(
"Problem: TS optimization did not converge.");
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()));
417 return cycles < maxiter;
421 int selectMode(
Core::Calculator& calc,
const std::vector<int>& relevantAtoms)
const {
422 _logger->output <<
" Automatically selecting normal mode " <<
Core::Log::endl;
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.");
430 Utils::NormalModesContainer modes = getNormalModes(calc);
431 auto waveNumbers = modes.getWaveNumbers();
432 std::vector<double> contributions;
434 for (
int i = 0; i < static_cast<int>(waveNumbers.size()); ++i) {
435 if (waveNumbers[i] >= 0.0) {
437 throw std::runtime_error(
"Structure has no imaginary normal mode.");
441 auto mode = modes.getMode(i);
442 double contribution = 0.0;
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();
449 contributions.push_back(contribution);
451 int selection = std::distance(contributions.begin(), std::max_element(contributions.begin(), contributions.end()));
452 double tolerance = contributions[selection] * 0.1;
455 for (
unsigned int i = 0; i < contributions.size(); i++) {
456 if ((contributions[selection] - contributions[i]) < tolerance) {
461 _logger->output <<
" Automatically selected mode " << std::to_string(selection) <<
" with wave number "
462 << std::to_string(waveNumbers[selection]) <<
" cm^-1" <<
Core::Log::endl;
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");
471 auto hessian = calc.results().get<Utils::Property::Hessian>();
472 auto system = calc.getStructure();
473 return Utils::NormalModeAnalysis::calculateNormalModes(hessian, system->getElements(), system->getPositions());
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