votca 2024.2-dev
Loading...
Searching...
No Matches
adiis.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/adiis.h"
26#include "votca/xtp/bfgs_trm.h"
27#include "votca/xtp/logger.h"
28
29namespace votca {
30namespace xtp {
31
32Eigen::VectorXd ADIIS::CalcCoeff(const std::vector<Eigen::MatrixXd>& dmathist,
33 const std::vector<Eigen::MatrixXd>& mathist) {
34 success = true;
35 Index size = dmathist.size();
36
37 const Eigen::MatrixXd& dmat = dmathist.back();
38 const Eigen::MatrixXd& H = mathist.back();
39 Eigen::VectorXd DiF = Eigen::VectorXd::Zero(size);
40 Eigen::MatrixXd DiFj = Eigen::MatrixXd::Zero(size, size);
41
42 for (Index i = 0; i < size; i++) {
43 DiF(i) = ((dmathist[i]) - dmat).cwiseProduct(H).sum();
44 }
45
46 for (Index i = 0; i < size; i++) {
47 for (Index j = 0; j < size; j++) {
48 DiFj(i, j) = ((dmathist[i]) - dmat).cwiseProduct((mathist[j]) - H).sum();
49 }
50 }
51
52 ADIIS_costfunction a_cost = ADIIS_costfunction(DiF, DiFj);
53 BFGSTRM optimizer = BFGSTRM(a_cost);
54 Logger log;
55 optimizer.setLog(&log);
56 optimizer.setNumofIterations(1000);
57 optimizer.setTrustRadius(0.01);
58 // Starting point: equal weights on all matrices
59 Eigen::VectorXd coeffs = Eigen::VectorXd::Constant(size, 1.0 / double(size));
60 optimizer.Optimize(coeffs);
61 success = optimizer.Success();
62 coeffs = optimizer.getParameters().cwiseAbs2();
63 double xnorm = coeffs.sum();
64 coeffs /= xnorm;
65
66 if (std::abs(coeffs.tail(1).value()) < 0.001) {
67 success = false;
68 }
69 return coeffs;
70}
71
72} // namespace xtp
73} // namespace votca
Eigen::VectorXd CalcCoeff(const std::vector< Eigen::MatrixXd > &dmathist, const std::vector< Eigen::MatrixXd > &mathist)
Definition adiis.cc:32
void Optimize(const Eigen::VectorXd &initialparameters)
Definition bfgs_trm.cc:31
const Eigen::VectorXd getParameters() const
Definition bfgs_trm.h:63
bool Success() const
Definition bfgs_trm.h:56
void setTrustRadius(double trust_radius)
Definition bfgs_trm.h:44
void setNumofIterations(Index iterations)
Definition bfgs_trm.h:52
void setLog(Logger *pLog)
Definition bfgs_trm.h:42
Logger is used for thread-safe output of messages.
Definition logger.h:164
base class for all analysis tools
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26