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