OOPS
SolveMatrixEquation.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_SOLVEMATRIXEQUATION_H_
9 #define TEST_ASSIMILATION_SOLVEMATRIXEQUATION_H_
10 
11 #include <string>
12 #include <vector>
13 
14 #include "eckit/testing/Test.h"
15 
19 #include "oops/assimilation/IPCG.h"
21 #include "oops/assimilation/PCG.h"
24 #include "oops/runs/Test.h"
25 #include "oops/util/Expect.h"
26 #include "oops/util/FloatCompare.h"
27 
29 
30 namespace test {
31 
33 
34  template <typename MINIMIZER>
35  void test_SolveMatrixEquation(MINIMIZER minimizer)
36  {
37  // Solve Ax = b using the MINIMIZER algorithm
38 
39  // Starting guess
40  Vector3D x(1, 2, 3);
41  const Vector3D b(30, 15, -20);
42  // A is a diagonal matrix
43  const Vector3D diagA(10, 5, -5);
44  const Matrix3D A(diagA);
45  // Deliberately incorrect preconditioning in order to trigger
46  // multiple iterations of the minimisation
47  const Vector3D diagprecond(1, 1, 1);
48  const Matrix3D precond(diagprecond);
49  const int maxiter = 10;
50  const double tolerance = 1e-9;
51 
52  const double normReduction = minimizer(x, b, A, precond, maxiter, tolerance);
53  EXPECT(oops::is_close_absolute(x.x(), 3.0, 1.0e-9));
54  EXPECT(oops::is_close_absolute(x.y(), 3.0, 1.0e-9));
55  EXPECT(oops::is_close_absolute(x.z(), 4.0, 1.0e-9));
56  EXPECT(oops::is_close_absolute(normReduction, 0.0, 1.0e-9));
57  }
58 
59  CASE("assimilation/SolveMatrixEquation/FGMRES") {
60  test_SolveMatrixEquation(oops::FGMRES<Vector3D, Matrix3D, Matrix3D>);
61  }
62 
63  CASE("assimilation/SolveMatrixEquation/GMRESR") {
64  test_SolveMatrixEquation(oops::GMRESR<Vector3D, Matrix3D, Matrix3D>);
65  }
66 
67  CASE("assimilation/SolveMatrixEquation/MINRES") {
68  test_SolveMatrixEquation(oops::MINRES<Vector3D, Matrix3D, Matrix3D>);
69  }
70 
71  CASE("assimilation/SolveMatrixEquation/IPCG") {
72  test_SolveMatrixEquation(oops::IPCG<Vector3D, Matrix3D, Matrix3D>);
73  }
74 
75  CASE("assimilation/SolveMatrixEquation/PCG") {
76  test_SolveMatrixEquation(oops::PCG<Vector3D, Matrix3D, Matrix3D>);
77  }
78 
79  CASE("assimilation/SolveMatrixEquation/PLanczos") {
80  test_SolveMatrixEquation(oops::PLanczos<Vector3D, Matrix3D, Matrix3D>);
81  }
82 
84  private:
85  std::string testid() const override {return "test::SolveMatrixEquation";}
86  void register_tests() const override {}
87  void clear() const override {}
88  };
89 
90 } // namespace test
91 
92 #endif // TEST_ASSIMILATION_SOLVEMATRIXEQUATION_H_
FGMRES solver for Ax=b.
GMRESR solver for Ax=b.
Inexact-Preconditioned Conjugate Gradients solver.
MINRES solver for Ax=b.
Preconditioned Conjugate Gradients solver.
Preconditioned Lanczos solver.
Diagonal matrix.
std::string testid() const override
void register_tests() const override
void clear() const override
double x() const
Definition: Vector3D.h:30
double y() const
Definition: Vector3D.h:31
double z() const
Definition: Vector3D.h:32
void test_SolveMatrixEquation(MINIMIZER minimizer)
oops::DiagonalMatrix< Vector3D > Matrix3D
CASE("test_linearmodelparameterswrapper_valid_name")