votca 2024.2-dev
Loading...
Searching...
No Matches
bfgs_trm.cc
Go to the documentation of this file.
1/*
2 * Copyright 2009-2020 The VOTCA Development Team
3 * (http://www.votca.org)
4 *
5 * Licensed under the Apache License, Version 2.0 (the "License")
6 *
7 * You may not use this file except in compliance with the License.
8 * You may obtain a copy of the License at
9 *
10 * http://www.apache.org/licenses/LICENSE-2.0
11 *
12 * Unless required by applicable law or agreed to in writing, software
13 * distributed under the License is distributed on an "AS IS" BASIS,
14 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
15 * See the License for the specific language governing permissions and
16 * limitations under the License.
17 *
18 */
19
20// Third party includes
21#include <boost/format.hpp>
22
23// Local VOTCA includes
24#include "votca/xtp/atom.h"
25#include "votca/xtp/bfgs_trm.h"
27
28namespace votca {
29namespace xtp {
30
31void BFGSTRM::Optimize(const Eigen::VectorXd& initialparameters) {
32 parameters_ = initialparameters;
34 double lastcost = cost_;
35 Eigen::VectorXd gradient = costfunction_.EvaluateGradient(parameters_);
36 for (auto& func : callbacks_) {
37 func();
38 }
39
40 Eigen::VectorXd delta_p_trial = Eigen::VectorXd::Zero(parameters_.size());
41 Eigen::VectorXd last_gradient = Eigen::VectorXd::Zero(parameters_.size());
42 double delta_cost = 0;
44 for (Index i = 0; i < 100; i++) {
45 TrustRegion subproblem;
46 delta_p_trial =
47 subproblem.CalculateStep(gradient, hessian_, trust_radius_);
48 double trialcost =
49 costfunction_.EvaluateCost(parameters_ + delta_p_trial);
50 delta_cost = trialcost - lastcost;
51 bool step_accepted =
52 AcceptRejectStep(delta_p_trial, gradient, delta_cost);
53 if (step_accepted) {
54 cost_ = trialcost;
55 parameters_ += delta_p_trial;
56 break;
57 }
58 }
60 if (iteration_ > 1) {
61 UpdateHessian(delta_p_trial, gradient - last_gradient);
62 }
63 lastcost = cost_;
64 last_gradient = gradient;
65 for (auto& func : callbacks_) {
66 func();
67 }
68 if (costfunction_.Converged(delta_p_trial, delta_cost, gradient)) {
69 break;
70 } else if (iteration_ == max_iteration_) {
71 success_ = false;
73 << (boost::format("BFGS-TRM @iteration %1$d: not converged after "
74 "%2$d iterations ") %
76 .str()
77 << std::flush;
78 }
79 }
80 return;
81}
82
83/* Accept/reject the new geometry and adjust trust radius, if required */
84bool BFGSTRM::AcceptRejectStep(const Eigen::VectorXd& delta_p,
85 const Eigen::VectorXd& gradient,
86 double cost_delta) {
87 bool step_accepted = false;
88 if (cost_delta > 0.0) {
89 // total energy has unexpectedly increased, half the trust radius
92 << (boost::format("BFGS-TRM @iteration %1$d: DeltaCost %2$2.4e step "
93 "rejected ") %
94 iteration_ % cost_delta)
95 .str()
96 << std::flush;
98 << (boost::format(
99 "BFGS-TRM @iteration %1$d: new trust radius %2$2.4e") %
101 .str()
102 << std::flush;
103
104 } else {
105 // total energy has decreased, we accept the step but might update the trust
106 // radius
107 step_accepted = true;
108 // adjust trust radius, if required
109 double tr_check = cost_delta / QuadraticEnergy(gradient, delta_p);
110 double norm_delta_p = delta_p.squaredNorm();
111 if (tr_check > 0.75 &&
112 1.25 * norm_delta_p > trust_radius_ * trust_radius_) {
114 } else if (tr_check < 0.25) {
116 }
118 << (boost::format(
119 "BFGS-TRM @iteration %1$d: DeltaCost/QuadraticApprox %2$2.4f "
120 "step accepted ") %
121 iteration_ % tr_check)
122 .str()
123 << std::flush;
125 << (boost::format(
126 "BFGS-TRM @iteration %1$d: new trust radius %2$2.4e") %
128 .str()
129 << std::flush;
130 }
131 return step_accepted;
132}
133
134void BFGSTRM::UpdateHessian(const Eigen::VectorXd& delta_pos,
135 const Eigen::VectorXd& delta_gradient) {
136 // second term in BFGS update (needs current Hessian)
137 hessian_ -= hessian_ * delta_pos * delta_pos.transpose() *
138 hessian_.transpose() /
139 (delta_pos.transpose() * hessian_ * delta_pos).value();
140 // first term in BFGS update
141 hessian_ += (delta_gradient * delta_gradient.transpose()) /
142 (delta_gradient.transpose() * delta_pos);
143 // symmetrize Hessian (since d2E/dxidxj should be symmetric)
144 hessian_ = 0.5 * (hessian_ + hessian_.transpose());
145 return;
146}
147
148/* Estimate energy change based on quadratic approximation */
149double BFGSTRM::QuadraticEnergy(const Eigen::VectorXd& gradient,
150 const Eigen::VectorXd& delta_pos) const {
151 return (gradient.transpose() * delta_pos).value() +
152 0.5 * (delta_pos.transpose() * hessian_ * delta_pos).value();
153}
154
155} // namespace xtp
156} // namespace votca
void Optimize(const Eigen::VectorXd &initialparameters)
Definition bfgs_trm.cc:31
Eigen::VectorXd parameters_
Definition bfgs_trm.h:84
double QuadraticEnergy(const Eigen::VectorXd &gradient, const Eigen::VectorXd &delta_pos) const
Definition bfgs_trm.cc:149
bool AcceptRejectStep(const Eigen::VectorXd &delta_pos, const Eigen::VectorXd &gradient, double energy_delta)
Definition bfgs_trm.cc:84
Eigen::MatrixXd hessian_
Definition bfgs_trm.h:83
Optimiser_costfunction & costfunction_
Definition bfgs_trm.h:68
void UpdateHessian(const Eigen::VectorXd &delta_pos, const Eigen::VectorXd &delta_gradient)
Definition bfgs_trm.cc:134
std::vector< std::function< void()> > callbacks_
Definition bfgs_trm.h:81
virtual double EvaluateCost(const Eigen::VectorXd &parameters)=0
virtual Eigen::VectorXd EvaluateGradient(const Eigen::VectorXd &parameters)=0
virtual bool Converged(const Eigen::VectorXd &delta_parameters, double delta_cost, const Eigen::VectorXd &gradient)=0
Eigen::VectorXd CalculateStep(const Eigen::VectorXd &gradient, const Eigen::MatrixXd &Hessian, double delta) const
#define XTP_LOG(level, log)
Definition logger.h:40
base class for all analysis tools
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26