votca 2024.1-dev
Loading...
Searching...
No Matches
trustregion.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// Standard includes
21#include <iostream>
22#include <vector>
23
24// Local VOTCA includes
26
27namespace votca {
28namespace xtp {
29
30Eigen::VectorXd TrustRegion::CalculateStep(const Eigen::VectorXd& gradient,
31 const Eigen::MatrixXd& Hessian,
32 double delta) const {
33 // calculate unrestricted step
34 Eigen::VectorXd freestep = Hessian.colPivHouseholderQr().solve(-gradient);
35 // if inside use the step;
36 if (freestep.norm() < delta) {
37 return freestep;
38 }
39
40 // calculate step on the boundary
41 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> es(Hessian);
42 const Eigen::VectorXd factor =
43 (es.eigenvectors().transpose() * gradient).cwiseAbs2();
44 double lambda = 0;
45 // hard case
46 if (std::abs(factor[0]) < 1e-18) {
47 std::vector<Index> index;
48 Eigen::VectorXd diag = Eigen::VectorXd::Zero(factor.size());
49 // construct pseudo inverse
50 for (Index i = 0; i < factor.size(); i++) {
51 double entry = es.eigenvalues()(i) - es.eigenvalues()(0);
52 if (std::abs(entry) < 1e-14) {
53 index.push_back(i);
54 continue;
55 } else {
56 diag(i) = 1 / entry;
57 }
58 }
59 Eigen::MatrixXd pseudo_inv =
60 es.eigenvectors() * diag.asDiagonal() * es.eigenvectors().transpose();
61 Eigen::VectorXd step = -pseudo_inv * gradient;
62 if (step.norm() < delta) {
63 double tau = std::sqrt(delta * delta - step.squaredNorm());
64 step += tau * es.eigenvectors().col(index[0]);
65 return step;
66 }
67 }
68 // sort out all the factors for the small eigenvalues,
69 Index start_index = 0;
70 for (; start_index < factor.size(); start_index++) {
71 if (factor[start_index] < 1e-18) {
72 continue;
73 }
74 break;
75 }
76
77 if (start_index == factor.size()) {
78 throw std::runtime_error(
79 "trustregion.cc: all the factors are close to zero trust region method "
80 "will not converge further.");
81 }
82
83 // start value for lambda a bit higher than lowest eigenvalue of Hessian
84 lambda =
85 -es.eigenvalues()(start_index) + std::sqrt(factor(start_index)) / delta;
86 TrustRegionFunction TRF = TrustRegionFunction(factor, es, delta, start_index);
87 for (Index iter = 0; iter < 100; iter++) {
88 std::pair<double, double> result = TRF.Evaluate(lambda);
89 double func_value = result.first;
90 double update = result.second;
91
92 if (update < 1e-14 || std::abs(func_value) < 1e-12) {
93 break;
94 }
95 lambda += update;
96 }
97
98 // this is effectively the solution of (H+I\lambda)*\Delta p=-g with \lambda
99 // adjusted so that ||p||=delta
100 Eigen::VectorXd new_delta_pos = Eigen::VectorXd::Zero(gradient.size());
101 for (Index i = start_index; i < gradient.size(); i++) {
102 new_delta_pos -= es.eigenvectors().col(i) *
103 (es.eigenvectors().col(i).transpose() * gradient) /
104 (es.eigenvalues()(i) + lambda);
105 }
106
107 // this is for safety
108 if (new_delta_pos.norm() > delta) {
109 new_delta_pos *= (delta / new_delta_pos.norm());
110 }
111 return new_delta_pos;
112}
113
114} // namespace xtp
115} // namespace votca
std::pair< double, double > Evaluate(double lambda)
Definition trustregion.h:56
Eigen::VectorXd CalculateStep(const Eigen::VectorXd &gradient, const Eigen::MatrixXd &Hessian, double delta) const
base class for all analysis tools
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26