OOPS
test/generic/VerticalLocEV.h
Go to the documentation of this file.
1 /*
2  * (C) Copyright 2020-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 TEST_GENERIC_VERTICALLOCEV_H_
9 #define TEST_GENERIC_VERTICALLOCEV_H_
10 
11 #include <Eigen/Dense>
12 
13 #include <cfloat>
14 #include <cmath>
15 #include <iostream>
16 #include <memory>
17 #include <string>
18 #include <vector>
19 
20 #define ECKIT_TESTING_SELF_REGISTER_CASES 0
21 
22 #include <boost/noncopyable.hpp>
23 
24 #include "eckit/config/LocalConfiguration.h"
25 #include "eckit/testing/Test.h"
26 #include "oops/base/Geometry.h"
27 #include "oops/base/Increment.h"
28 #include "oops/base/Increment4D.h"
30 #include "oops/base/Variables.h"
32 #include "oops/runs/Test.h"
33 #include "oops/util/DateTime.h"
34 #include "oops/util/dot_product.h"
36 #include "test/TestEnvironment.h"
37 
38 namespace test {
39 
40 // =============================================================================
41 
42 template <typename MODEL> void testVerticalLocEV() {
43  typedef IncrementFixture<MODEL> Test_;
44  typedef oops::Geometry<MODEL> Geometry_;
45  typedef oops::VerticalLocEV<MODEL> VerticalLocEV_;
46  typedef oops::Increment4D<MODEL> Increment4D_;
47  typedef oops::IncrementEnsemble4D<MODEL> IncrementEnsemble_;
48  typedef oops::State<MODEL> State_;
49 
50  const Geometry_ & geometry = Test_::resol();
51  eckit::LocalConfiguration vertlocconf(TestEnvironment::config(), "vertical localization");
52 
53  // make an empty state vector to be used to intitialize VerticalLocEV_
54  State_ x(Test_::resol(), Test_::ctlvars(), Test_::time());
55 
56  VerticalLocEV_ vertloc(vertlocconf, x);
57  oops::Log::test() << "Number of eigenvalues used in VerticalLoc: " << vertloc.neig() << std::endl;
58 
59  //--- check for expected number of eigen modes
60  int nEigExpected = TestEnvironment::config().getInt("expected neig");
61  int neig = vertloc.neig();
62  oops::Log::debug() << "Expected number of eigen modes: " << nEigExpected << std::endl;
63  oops::Log::debug() << "Actual number of eigen modes: " << neig << std::endl;
64  EXPECT(nEigExpected == neig);
65 
66  //--- check that truncation and rescaling was done correctly
67  EXPECT(vertloc.testTruncateEvecs(geometry));
68 
69  //--- check modulation
70  // if the increment dx=1, then modulated ensemble will be the same as (scaled) eigen vectors
71  // then the dot products will satisfy orthogonality condition
72 
73  // need at least 2 eigen vectors for the following check to work.
74  EXPECT(neig > 1);
75 
76  // construct instances of Increment4D_ and IncrementEnsemble_
77  std::vector<util::DateTime> times;
78  times.push_back(Test_::time());
79  Increment4D_ dx1(Test_::resol(), Test_::ctlvars(), times);
80  Increment4D_ dx2(Test_::resol(), Test_::ctlvars(), times);
81  IncrementEnsemble_ incEns(Test_::resol(), Test_::ctlvars(), times, neig);
82  // set incEns to zero
83  for (int i = 1; i < neig; ++i) {incEns[i][0].zero();}
84  for (int i = 1; i < neig; ++i) {
85  double n = incEns[0][0].dot_product_with(incEns[i][0]);
86  EXPECT(n < 100*DBL_EPSILON);
87  }
88 
89  // make increment of ones and check that it is indeed full of ones
90  dx1.ones();
91  dx2.random();
92  double normRandom = dx2.dot_product_with(dx2);
93  dx2.schur_product_with(dx1);
94  EXPECT(std::abs(dx2.dot_product_with(dx2)-normRandom) < normRandom*DBL_EPSILON);
95  oops::Log::debug() << "Increment ones()" << dx1 << std::endl;
96 
97  // modulate increments
98  vertloc.modulateIncrement(dx1, incEns);
99 
100  // check the orthogonality condition
101  double n0 = incEns[0][0].dot_product_with(incEns[0][0]);
102  oops::Log::debug() << "dot product 0: " << n0 << std::endl;
103  // check that eigen vectors are not zeros
104  EXPECT(n0 > 0);
105  double tol = 2*n0*DBL_EPSILON;
106  oops::Log::debug() << "tolerance :" << tol << std::endl;
107 
108  // check that eig[ieig>0] are orthogonal to eig[0]
109  for (int i = 1; i < neig; ++i) {
110  double n = incEns[0][0].dot_product_with(incEns[i][0]);
111  oops::Log::debug() << "dot product " << i << ": " << n << std::endl;
112  EXPECT(n < tol);
113  }
114 
115  // try the second interface for modulateIncrement
116  // this checks that modulation of a single column @geom.begin() works as well
117  IncrementEnsemble_ incEns2(Test_::resol(), Test_::ctlvars(), times, 1);
118  incEns2[0] = dx1;
119 
120  Eigen::MatrixXd modInc = vertloc.modulateIncrement(incEns2, geometry.begin(), 0);
121  Eigen::MatrixXd modIncInner = modInc.transpose()*modInc;
122  oops::Log::debug() << "modInc'*modInc" << modIncInner << std::endl;
123  // modIncInner should be a diagonal matrix
124  for (int i = 1; i < neig; ++i) {
125  EXPECT(modIncInner(0, i) < modIncInner(0, 0)*DBL_EPSILON);
126  }
127 }
128 
129 // =============================================================================
130 
131 template <typename MODEL>
132 class VerticalLocEV : public oops::Test {
133  public:
134  VerticalLocEV() = default;
135  virtual ~VerticalLocEV() = default;
136 
137  private:
138  std::string testid() const override {return "test::VerticalLocEV<" + MODEL::name() + ">";}
139 
140  void register_tests() const override {
141  std::vector<eckit::testing::Test>& ts = eckit::testing::specification();
142 
143  ts.emplace_back(CASE("generic/VerticalLocEV/testVerticalLocEV")
144  { testVerticalLocEV<MODEL>(); });
145  }
146 
147  void clear() const override {}
148 };
149 
150 // =============================================================================
151 
152 } // namespace test
153 
154 #endif // TEST_GENERIC_VERTICALLOCEV_H_
Geometry class used in oops; subclass of interface class interface::Geometry.
4D model state Increment (vector of 3D Increments)
Definition: Increment4D.h:32
Ensemble of 4D increments.
State class used in oops; subclass of interface class interface::State.
static const eckit::Configuration & config()
virtual ~VerticalLocEV()=default
std::string testid() const override
void register_tests() const override
VerticalLocEV()=default
void clear() const override
void testVerticalLocEV()
CASE("test_linearmodelparameterswrapper_valid_name")