OOPS
test/assimilation/TriDiagSolve.h
Go to the documentation of this file.
1 /*
2  * (C) Copyright 2020 Met Office UK
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_ASSIMILATION_TRIDIAGSOLVE_H_
9 #define TEST_ASSIMILATION_TRIDIAGSOLVE_H_
10 
11 #include <Eigen/Dense>
12 
13 #include <string>
14 #include <vector>
15 
16 #include "eckit/testing/Test.h"
17 
20 #include "oops/runs/Test.h"
21 #include "oops/util/Expect.h"
22 #include "oops/util/FloatCompare.h"
23 
24 namespace test {
25 
27  {
28  const std::vector<double> diag{1.0, 1.0, 1.0};
29  const std::vector<double> sub{0.5, 0.5};
30  const std::vector<double> rhs{1.0, 0.5, 0.0};
31  std::vector<double> sol;
32 
33  oops::TriDiagSolve(diag, sub, rhs, sol);
34 
35  EXPECT_EQUAL(sol[0], 1.0);
36  EXPECT_EQUAL(sol[1], 0.0);
37  EXPECT_EQUAL(sol[2], 0.0);
38  }
39 
41  {
42  const int members = 3;
43  Eigen::MatrixXd A0 = Eigen::MatrixXd::Ones(members, members);
44  Eigen::MatrixXd A1 = Eigen::MatrixXd::Ones(members, members);
45  Eigen::MatrixXd B0 = Eigen::MatrixXd::Ones(members, members);
46  Eigen::MatrixXd B1 = Eigen::MatrixXd::Ones(members, members);
47  B0(1, 0) = 0;
48  B0(2, 0) = 0;
49  B0(2, 1) = 0;
50  B1(1, 0) = 0;
51  B1(2, 0) = 0;
52  B1(2, 1) = 0;
53 
54  Eigen::MatrixXd beta0 = Eigen::MatrixXd::Ones(members, members);
55  beta0(1, 0) = 0;
56  beta0(2, 0) = 0;
57  beta0(2, 1) = 0;
58  Eigen::MatrixXd ss;
59  bool complexValues = false;
60 
61  oops::blockTriDiagSolve({A0, A1}, {B0, B1}, beta0, ss, complexValues, members);
62 
63  EXPECT(oops::is_close_absolute(ss(0, 0), 4.0, 1e-9));
64  EXPECT(oops::is_close_absolute(ss(0, 1), 2.0, 1e-9));
65  EXPECT(oops::is_close_absolute(ss(0, 2), 2.0, 1e-9));
66  EXPECT(oops::is_close_absolute(ss(1, 0), -2.0, 1e-9));
67  EXPECT(oops::is_close_absolute(ss(1, 1), 1.0, 1e-9));
68  EXPECT(oops::is_close_absolute(ss(1, 2), 0.0, 1e-9));
69  EXPECT(oops::is_close_absolute(ss(2, 0), 0.0, 1e-9));
70  EXPECT(oops::is_close_absolute(ss(2, 1), -1.0, 1e-9));
71  EXPECT(oops::is_close_absolute(ss(2, 2), 1.0, 1e-9));
72  EXPECT(oops::is_close_absolute(ss(3, 0), -1.0, 1e-9));
73  EXPECT(oops::is_close_absolute(ss(3, 1), -1.0, 1e-9));
74  EXPECT(oops::is_close_absolute(ss(3, 2), -2.0, 1e-9));
75  EXPECT(oops::is_close_absolute(ss(4, 0), 2.0, 1e-9));
76  EXPECT(oops::is_close_absolute(ss(4, 1), 1.0, 1e-9));
77  EXPECT(oops::is_close_absolute(ss(4, 2), 1.0, 1e-9));
78  EXPECT(oops::is_close_absolute(ss(5, 0), -1.0, 1e-9));
79  EXPECT(oops::is_close_absolute(ss(5, 1), 0.0, 1e-9));
80  EXPECT(oops::is_close_absolute(ss(5, 2), 0.0, 1e-9));
81  EXPECT_EQUAL(complexValues, false);
82  }
83 
84  CASE("assimilation/TriDiagSolve/TriDiagSolve") {
86  }
87 
88  CASE("assimilation/TriDiagSolve/blockTriDiagSolve") {
90  }
91 
92  class TriDiagSolve : public oops::Test {
93  private:
94  std::string testid() const override {return "test::TriDiagSolve";}
95  void register_tests() const override {}
96  void clear() const override {}
97  };
98 
99 } // namespace test
100 
101 #endif // TEST_ASSIMILATION_TRIDIAGSOLVE_H_
std::string testid() const override
void register_tests() const override
void TriDiagSolve(const std::vector< double > &diag, const std::vector< double > &sub, const std::vector< double > &rhs, std::vector< double > &sol)
void blockTriDiagSolve(const std::vector< Eigen::MatrixXd > &alphas, const std::vector< Eigen::MatrixXd > &betas, const Eigen::MatrixXd &beta0, Eigen::MatrixXd &ss, bool &complexValues, const int members)
CASE("test_linearmodelparameterswrapper_valid_name")