OOPS
MinimizerUtils.h
Go to the documentation of this file.
1 /*
2  * (C) Copyright 2020 UCAR.
3  *
4  * This software is licensed under the terms of the Apache Licence Version 2.0
5  * which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
6  */
7 
8 #ifndef OOPS_ASSIMILATION_MINIMIZERUTILS_H_
9 #define OOPS_ASSIMILATION_MINIMIZERUTILS_H_
10 
11 #include <memory>
12 #include <vector>
13 
14 #include "eckit/config/Configuration.h"
15 
18 #include "oops/util/Logger.h"
19 
20 namespace oops {
21 
22 /// Prints to Log::info gradient reduction \p grad and normalized gradient reduction \p norm
23 /// for iteration \p iteration
24 void printNormReduction(int iteration, const double & grad, const double & norm);
25 /// Prints to Log::info cost function values for \p costJ, \p costJb, \p costJoJc for
26 /// iteration \p iteration
27 void printQuadraticCostFunction(int iteration, const double & costJ, const double & costJb,
28  const double & costJoJc);
29 
30 // -----------------------------------------------------------------------------
31 
32 template<typename MODEL, typename OBS>
33 void writeIncrement(const eckit::Configuration & config,
34  const ControlIncrement<MODEL, OBS> & dx, const int & loop) {
35 // Write out the increment
36 
37  if (config.has("online diagnostics")) {
38  const eckit::LocalConfiguration onlineDiag(config, "online diagnostics");
39  bool writeinc = onlineDiag.getBool("write increment", false);
40 
41  if (writeinc) {
42  // print log
43  Log::info() << "Write Increment - starting: " << loop << std::endl << std::endl;
44 
45  const eckit::LocalConfiguration incConf(onlineDiag, "increment");
46 
47  // write increment
48  dx.write(incConf);
49 
50  // print log
51  Log::info() << std::endl << "Write Increment: done." << std::endl;
52  }
53  }
54 }
55 
56 // -----------------------------------------------------------------------------
57 
58 template<typename MODEL, typename OBS>
59 void writeKrylovBasis(const eckit::Configuration & config,
61  const int & loop) {
62 // Write out the increment
63  if (config.has("online diagnostics")) {
64  eckit::LocalConfiguration diagConf(config, "online diagnostics");
65  bool writeinc = diagConf.getBool("write basis", false);
66 
67  if (writeinc) {
68  // print log
69  Log::info() << "Write Krylov Basis: starting: " << loop << std::endl;
70 
71  eckit::LocalConfiguration basisConf(diagConf, "krylov basis");
72  basisConf.set("iteration", loop);
73 
74  // write increment
75  dx.write(basisConf);
76 
77  // print log
78  Log::info() << "Write Krylov Basis: done." << std::endl;
79  }
80  }
81 }
82 
83 //-------------------------------------------------------------------------------------------------
84 template<typename MODEL, typename OBS>
85 void writeEigenvectors(const eckit::Configuration & diagConf,
86  const std::vector<double> & diag,
87  const std::vector<double> & sub,
88  const std::vector<double> & rhs,
89  std::vector<std::unique_ptr<ControlIncrement<MODEL, OBS>>> & zvecs,
90  std::vector<std::unique_ptr<ControlIncrement<MODEL, OBS>>> & hvecs,
91  const HtRinvHMatrix<MODEL, OBS> & HtRinvH,
95  if (diagConf.has("online diagnostics")) {
96  int maxEigen = diagConf.getInt("online diagnostics.max eigenvectors", 0);
97 
98  const double nn = rhs.size();
99  Eigen::MatrixXd TT = Eigen::MatrixXd::Zero(nn, nn);
100 
101  for (int ii = 0; ii < nn; ++ii) {
102  TT(ii, ii) = diag[ii];
103  if (ii > 0) TT(ii-1, ii) = sub[ii-1];
104  if (ii < nn - 1) TT(ii+1, ii) = sub[ii];
105  }
106 
107  Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> soluce(TT);
108  Eigen::MatrixXd eigenvecT = soluce.eigenvectors().real();
109  Eigen::VectorXd eigenvalT = soluce.eigenvalues().real();
110 
111  // Compute the eigenvectors (y = Zx)
112  for (int ii = 0; ii < maxEigen && ii < nn; ++ii) {
113  eigenz.zero();
114  eigenv.zero();
115  for (unsigned int jj = 0; jj < nn; ++jj) {
116  temp.zero();
117  temp = *zvecs[jj];
118  temp *= eigenvecT(jj, nn - 1 - ii);
119  eigenz += temp;
120  temp.zero();
121  temp = *hvecs[jj];
122  temp *= eigenvecT(jj, nn - 1 - ii);
123  eigenv += temp;
124  }
125  // Save the eigenvector
126  eckit::LocalConfiguration basisConf(diagConf, "online diagnostics.eigenvector");
127  basisConf.set("iteration", ii);
128  eigenz.write(basisConf);
129 
130  // Verification that eigenz is an eigenvector:
131  // A.eigenv = eigenv + HtRinvH.Beigenv = eigenv + HtRinvH.eigenz = lambda eigenv
132  temp.zero();
133  HtRinvH.multiply(eigenz, temp);
134  temp += eigenv;
135  eigenv *= eigenvalT(nn - 1 - ii);
136  temp -= eigenv;
137  Log::info() << "Eigenvalue " << ii+1 << " : " << eigenvalT(nn - 1 - ii) << std::endl;
138  Log::info() << "Norm A*y-lambda*y = " << dot_product(temp, temp) << std::endl;
139 
140  Log::test() << "Eigenvalue " << ii+1 << " : " << eigenvalT(nn - 1 - ii) << std::endl;
141  Log::test() << "Norm eigenvector = " << dot_product(eigenz, eigenz) << std::endl;
142  } // end for()
143  } // end if()
144 }
145 
146 // -----------------------------------------------------------------------------
147 
148 } // namespace oops
149 
150 #endif // OOPS_ASSIMILATION_MINIMIZERUTILS_H_
void zero()
Linear algebra operators.
void write(const eckit::Configuration &) const
void multiply(const CtrlInc_ &dx, CtrlInc_ &dz) const
Definition: HtRinvHMatrix.h:63
The namespace for the main oops code.
void writeKrylovBasis(const eckit::Configuration &config, const ControlIncrement< MODEL, OBS > &dx, const int &loop)
void printQuadraticCostFunction(int iteration, const double &costJ, const double &costJb, const double &costJoJc)
void writeIncrement(const eckit::Configuration &config, const ControlIncrement< MODEL, OBS > &dx, const int &loop)
void writeEigenvectors(const eckit::Configuration &diagConf, const std::vector< double > &diag, const std::vector< double > &sub, const std::vector< double > &rhs, std::vector< std::unique_ptr< ControlIncrement< MODEL, OBS >>> &zvecs, std::vector< std::unique_ptr< ControlIncrement< MODEL, OBS >>> &hvecs, const HtRinvHMatrix< MODEL, OBS > &HtRinvH, ControlIncrement< MODEL, OBS > &temp, ControlIncrement< MODEL, OBS > &eigenv, ControlIncrement< MODEL, OBS > &eigenz)
void printNormReduction(int iteration, const double &grad, const double &norm)