11 #ifndef OOPS_ASSIMILATION_QNEWTONLMP_H_
12 #define OOPS_ASSIMILATION_QNEWTONLMP_H_
18 #include "eckit/config/LocalConfiguration.h"
19 #include "oops/util/dot_product.h"
20 #include "oops/util/Logger.h"
40 template<
typename VECTOR,
typename BMATRIX>
class QNewtonLMP {
45 void push(
const VECTOR &,
const VECTOR &,
const VECTOR &,
const double &);
60 std::vector<VECTOR>
P_;
61 std::vector<VECTOR>
Ph_;
62 std::vector<VECTOR>
AP_;
75 template<
typename VECTOR,
typename BMATRIX>
77 : maxpairs_(0), maxnewpairs_(0), useoldpairs_(false), maxouter_(0), update_(1)
80 Log::info() <<
"QNewtonLMP: maxouter : " <<
maxouter_ << std::endl;
82 if (conf.has(
"preconditioner")) {
83 const eckit::LocalConfiguration precond(conf,
"preconditioner");
84 if (precond.has(
"maxpairs"))
maxpairs_ = precond.getInt(
"maxpairs");
87 if (precond.has(
"useoldpairs")) old = precond.getString(
"useoldpairs");
101 template<
typename VECTOR,
typename BMATRIX>
103 const VECTOR & ap,
const double & rho) {
104 ASSERT(savedP_.size() <= maxnewpairs_);
105 if (maxnewpairs_ > 0 && update_ < maxouter_) {
106 if (savedP_.size() == maxnewpairs_) {
107 savedP_.erase(savedP_.begin());
108 savedPh_.erase(savedPh_.begin());
109 savedAP_.erase(savedAP_.begin());
110 savedrhos_.erase(savedrhos_.begin());
112 savedP_.push_back(p);
113 savedPh_.push_back(ph);
114 savedAP_.push_back(ap);
115 savedrhos_.push_back(1.0/rho);
121 template<
typename VECTOR,
typename BMATRIX>
123 Log::debug() <<
"QNewtonLMP size saved Ph = " << savedPh_.size() << std::endl;
124 Log::debug() <<
"QNewtonLMP size saved P = " << savedP_.size() << std::endl;
125 Log::debug() <<
"QNewtonLMP size saved AP = " << savedAP_.size() << std::endl;
126 Log::debug() <<
"QNewtonLMP size saved rhos = " << savedrhos_.size() << std::endl;
128 const unsigned nvec = savedPh_.size();
129 ASSERT(savedP_.size() == nvec);
130 ASSERT(savedAP_.size() == nvec);
131 ASSERT(savedrhos_.size() == nvec);
133 if (!useoldpairs_ || nvec >= maxpairs_) {
139 for (
unsigned kiter = 0; kiter < usedpairIndx_.size(); ++kiter) {
140 usedpairIndx_[kiter] = 0;
144 if (nvec > 0 && update_ < maxouter_) {
145 Log::info() <<
"QNewtonLMP: update " << update_ <<
", max = " << maxouter_-1 << std::endl;
146 unsigned newpairs = std::min(nvec, maxpairs_);
147 unsigned oldpairs = P_.size();
148 unsigned rmpairs = 0;
150 if (oldpairs + newpairs > maxpairs_) {
151 rmpairs = oldpairs + newpairs - maxpairs_;
152 for (
unsigned jv = 0; jv < rmpairs; ++jv) {
153 Ph_.erase(Ph_.begin());
154 P_.erase(P_.begin());
155 AP_.erase(AP_.begin());
156 BAP_.erase(BAP_.begin());
157 rhos_.erase(rhos_.begin());
161 unsigned removed = rmpairs;
162 for (
unsigned kiter = 0; kiter < usedpairIndx_.size(); ++kiter) {
163 unsigned rmiter = std::min(usedpairIndx_[kiter], removed);
164 usedpairIndx_[kiter] -= rmiter;
167 ASSERT(removed == 0);
169 ASSERT(P_.size() == oldpairs);
170 ASSERT(oldpairs + newpairs <= maxpairs_);
173 for (
unsigned jv = nvec - newpairs; jv < nvec; ++jv) {
174 P_.push_back(savedP_[jv]);
175 Ph_.push_back(savedPh_[jv]);
176 AP_.push_back(savedAP_[jv]);
177 rhos_.push_back(savedrhos_[jv]);
179 VECTOR ww(savedAP_[jv]);
180 Bmat.multiply(savedAP_[jv], ww);
183 ASSERT(P_.size() == oldpairs + newpairs);
185 Log::info() <<
"Number of inner iterations : " << nvec << std::endl;
186 Log::info() <<
"Number of maximum pairs : " << maxpairs_ << std::endl;
187 Log::info() <<
"Number of used total pairs : " << oldpairs + newpairs << std::endl;
188 Log::info() <<
"Number of new pairs : " << newpairs << std::endl;
189 Log::info() <<
"Number of removed old pairs : " << rmpairs << std::endl;
190 for (
unsigned kiter = 0; kiter < usedpairIndx_.size(); ++kiter) {
191 Log::info() <<
"Number of used pairs from outer loop " << kiter + 1
192 <<
" : " << usedpairIndx_[kiter];
205 template<
typename VECTOR,
typename BMATRIX>
208 const unsigned nvec = P_.size();
210 std::vector<double> etas;
212 for (
int iiter = nvec-1; iiter >= 0; iiter--) {
213 double eta = dot_product(b, P_[iiter]);
216 b.axpy(-eta, AP_[iiter]);
218 b *= dot_product(AP_[nvec-1], AP_[nvec-1])/dot_product(AP_[nvec-1], Ph_[nvec-1]);
219 for (
unsigned iiter = 0; iiter < nvec; ++iiter) {
220 double sigma = dot_product(b, BAP_[iiter]);
221 sigma *= rhos_[iiter];
222 sigma -= etas[iiter];
223 b.axpy(-sigma, Ph_[iiter]);
230 template<
typename VECTOR,
typename BMATRIX>
233 const unsigned nvec = P_.size();
235 std::vector<double> etas;
237 for (
int iiter = nvec-1; iiter >= 0; iiter--) {
238 double eta = dot_product(b, Ph_[iiter]);
241 b.axpy(-eta, BAP_[iiter]);
243 b *= dot_product(AP_[nvec-1], AP_[nvec-1])/dot_product(AP_[nvec-1], Ph_[nvec-1]);
244 for (
unsigned iiter = 0; iiter < nvec; ++iiter) {
245 double sigma = dot_product(b, AP_[iiter]);
246 sigma *= rhos_[iiter];
247 sigma -= etas[iiter];
248 b.axpy(-sigma, P_[iiter]);
std::vector< unsigned > usedpairIndx_
QNewtonLMP(const eckit::Configuration &)
std::vector< VECTOR > Ph_
void push(const VECTOR &, const VECTOR &, const VECTOR &, const double &)
void multiply(const VECTOR &, VECTOR &) const
std::vector< double > rhos_
std::vector< VECTOR > savedAP_
std::vector< VECTOR > BAP_
void update(const BMATRIX &B)
std::vector< VECTOR > savedPh_
std::vector< VECTOR > AP_
void tmultiply(const VECTOR &, VECTOR &) const
std::vector< double > savedrhos_
std::vector< VECTOR > savedP_
The namespace for the main oops code.