votca 2026-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
40
41// Restricted ADIIS coefficients. The minimization is carried out for the
42// approximate energy functional
43//
44// E(c) = sum_i c_i D_iF + 1/2 sum_ij c_i c_j D_iF_j,
45//
46// under the simplex constraints c_i >= 0 and sum_i c_i = 1. The latter is
47// enforced through the squared optimizer parameters followed by
48// renormalization.
49Eigen::VectorXd ADIIS::CalcCoeff(const std::vector<Eigen::MatrixXd>& dmathist,
50 const std::vector<Eigen::MatrixXd>& mathist) {
51 success = true;
52 Index size = dmathist.size();
53
54 const Eigen::MatrixXd& dmat = dmathist.back();
55 const Eigen::MatrixXd& H = mathist.back();
56 Eigen::VectorXd DiF = Eigen::VectorXd::Zero(size);
57 Eigen::MatrixXd DiFj = Eigen::MatrixXd::Zero(size, size);
58
59 for (Index i = 0; i < size; i++) {
60 DiF(i) = ((dmathist[i]) - dmat).cwiseProduct(H).sum();
61 }
62
63 for (Index i = 0; i < size; i++) {
64 for (Index j = 0; j < size; j++) {
65 DiFj(i, j) = ((dmathist[i]) - dmat).cwiseProduct((mathist[j]) - H).sum();
66 }
67 }
68
69 ADIIS_costfunction a_cost = ADIIS_costfunction(DiF, DiFj);
70 BFGSTRM optimizer = BFGSTRM(a_cost);
71 Logger log;
72 optimizer.setLog(&log);
73 optimizer.setNumofIterations(1000);
74 optimizer.setTrustRadius(0.01);
75 // Starting point: equal weights on all matrices
76 Eigen::VectorXd coeffs = Eigen::VectorXd::Constant(size, 1.0 / double(size));
77 optimizer.Optimize(coeffs);
78 success = optimizer.Success();
79 coeffs = optimizer.getParameters().cwiseAbs2();
80 double xnorm = coeffs.sum();
81 coeffs /= xnorm;
82
83 if (std::abs(coeffs.tail(1).value()) < 0.001) {
84 success = false;
85 }
86 return coeffs;
87}
88
89// Unrestricted ADIIS variant. The same quadratic surrogate is assembled, but
90// with alpha and beta traces added so that the objective approximates the spin-
91// summed SCF energy change.
92Eigen::VectorXd ADIIS::CalcCoeff(
93 const std::vector<Eigen::MatrixXd>& dmathist_alpha,
94 const std::vector<Eigen::MatrixXd>& dmathist_beta,
95 const std::vector<Eigen::MatrixXd>& mathist_alpha,
96 const std::vector<Eigen::MatrixXd>& mathist_beta) {
97
98 success = true;
99 Index size = dmathist_alpha.size();
100
101 const Eigen::MatrixXd& dmat_alpha = dmathist_alpha.back();
102 const Eigen::MatrixXd& dmat_beta = dmathist_beta.back();
103 const Eigen::MatrixXd& H_alpha = mathist_alpha.back();
104 const Eigen::MatrixXd& H_beta = mathist_beta.back();
105
106 Eigen::VectorXd DiF = Eigen::VectorXd::Zero(size);
107 Eigen::MatrixXd DiFj = Eigen::MatrixXd::Zero(size, size);
108
109 for (Index i = 0; i < size; ++i) {
110 DiF(i) = (dmathist_alpha[i] - dmat_alpha).cwiseProduct(H_alpha).sum() +
111 (dmathist_beta[i] - dmat_beta).cwiseProduct(H_beta).sum();
112 }
113
114 for (Index i = 0; i < size; ++i) {
115 for (Index j = 0; j < size; ++j) {
116 DiFj(i, j) = (dmathist_alpha[i] - dmat_alpha)
117 .cwiseProduct(mathist_alpha[j] - H_alpha)
118 .sum() +
119 (dmathist_beta[i] - dmat_beta)
120 .cwiseProduct(mathist_beta[j] - H_beta)
121 .sum();
122 }
123 }
124
125 ADIIS_costfunction a_cost(DiF, DiFj);
126 BFGSTRM optimizer(a_cost);
127 Logger log;
128 optimizer.setLog(&log);
129 optimizer.setNumofIterations(1000);
130 optimizer.setTrustRadius(0.01);
131
132 Eigen::VectorXd coeffs = Eigen::VectorXd::Constant(size, 1.0 / double(size));
133 optimizer.Optimize(coeffs);
134 success = optimizer.Success();
135
136 coeffs = optimizer.getParameters().cwiseAbs2();
137 double xnorm = coeffs.sum();
138 coeffs /= xnorm;
139
140 if (std::abs(coeffs.tail(1).value()) < 0.001) {
141 success = false;
142 }
143 return coeffs;
144}
145
146} // namespace xtp
147} // namespace votca
Eigen::VectorXd CalcCoeff(const std::vector< Eigen::MatrixXd > &dmathist, const std::vector< Eigen::MatrixXd > &mathist)
Definition adiis.cc:49
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
Charge transport classes.
Definition ERIs.h:28
Provides a means for comparing floating point numbers.
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26