votca 2024.2-dev
Loading...
Searching...
No Matches
adiis_costfunction.h
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#pragma once
21#ifndef VOTCA_XTP_ADIIS_COSTFUNCTION_H
22#define VOTCA_XTP_ADIIS_COSTFUNCTION_H
23
24// Local VOTCA includes
26
27namespace votca {
28namespace xtp {
29
31 public:
32 ADIIS_costfunction(Eigen::VectorXd DiF, Eigen::MatrixXd DiFj) {
33 DiF_ = DiF;
34 DiFj_ = DiFj;
35 }
36
37 double EvaluateCost(const Eigen::VectorXd& parameters) override {
38 Eigen::VectorXd c = parameters.cwiseAbs2();
39 double xnorm = c.sum();
40 c /= xnorm;
41 return (2 * c.transpose() * DiF_ + c.transpose() * DiFj_ * c).value();
42 }
43
44 Eigen::VectorXd EvaluateGradient(const Eigen::VectorXd& parameters) override {
45 Eigen::VectorXd c = parameters.cwiseAbs2();
46 double xnorm = c.sum();
47 c /= xnorm;
48 Eigen::VectorXd dEdc = 2.0 * DiF_ + DiFj_ * c + DiFj_.transpose() * c;
49 Eigen::MatrixXd jac = Eigen::MatrixXd::Zero(c.size(), c.size());
50 for (Index i = 0; i < jac.rows(); i++) {
51 for (Index j = 0; j < jac.cols(); j++) {
52 jac(i, j) = -c(i) * 2.0 * parameters(j) / xnorm;
53 }
54 // Extra term on diagonal
55 jac(i, i) += 2.0 * parameters(i) / xnorm;
56 }
57 return jac.transpose() * dEdc;
58 }
59
60 Index NumParameters() const override { return Index(DiF_.size()); }
61
62 bool Converged(const Eigen::VectorXd&, double,
63 const Eigen::VectorXd& gradient) override {
64 return gradient.cwiseAbs().maxCoeff() < 1.e-7;
65 }
66
67 private:
68 Eigen::VectorXd DiF_;
69 Eigen::MatrixXd DiFj_;
70};
71
72} // namespace xtp
73} // namespace votca
74#endif // VOTCA_XTP_ADIIS_COSTFUNCTION_H
Index NumParameters() const override
bool Converged(const Eigen::VectorXd &, double, const Eigen::VectorXd &gradient) override
Eigen::VectorXd EvaluateGradient(const Eigen::VectorXd &parameters) override
double EvaluateCost(const Eigen::VectorXd &parameters) override
ADIIS_costfunction(Eigen::VectorXd DiF, Eigen::MatrixXd DiFj)
base class for all analysis tools
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26