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