OOPS
oops/assimilation/TriDiagSolve.h
Go to the documentation of this file.
1 /*
2  * (C) Copyright 2009-2016 ECMWF.
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  * In applying this licence, ECMWF does not waive the privileges and immunities
7  * granted to it by virtue of its status as an intergovernmental organisation nor
8  * does it submit to any jurisdiction.
9  */
10 
11 #ifndef OOPS_ASSIMILATION_TRIDIAGSOLVE_H_
12 #define OOPS_ASSIMILATION_TRIDIAGSOLVE_H_
13 
14 #include <Eigen/Dense>
15 #include <Eigen/Eigenvalues>
16 
17 #include <vector>
18 
19 #include "oops/mpi/mpi.h"
20 #include "oops/util/dot_product.h"
21 #include "oops/util/Logger.h"
22 
23 namespace oops {
24 
25 void TriDiagSolve(const std::vector<double> & diag, const std::vector<double> & sub,
26  const std::vector<double> & rhs, std::vector<double> & sol) {
27  const double n = rhs.size();
28  ASSERT(sub.size() == n-1);
29  ASSERT(diag.size() == n);
30  sol.resize(n);
31  std::vector<double> c(sub);
32  std::vector<double> d(rhs);
33 
34  // Forward substitution
35  c[0] = c[0] / diag[0];
36  d[0] = d[0] / diag[0];
37  if (n > 2) {
38  for (int ii = 1; ii <= n-2; ++ii) {
39  double denom = diag[ii]-sub[ii-1]*c[ii-1];
40  c[ii] = c[ii]/denom;
41  d[ii] = (d[ii]-sub[ii-1]*d[ii-1])/denom;
42  }
43  }
44  d[n-1] = (d[n-1]-sub[n-2]*d[n-2])/(diag[n-1]-sub[n-2]*c[n-2]);
45  // Backward substitution
46  sol[n-1] = d[n-1];
47  for (int ii = n-2; ii >= 0; --ii) {
48  sol[ii] = (d[ii]-(c[ii]*sol[ii+1]));
49  }
50 }
51 
52 //-------------------------------------------------------------------------------------------------
53 // Computes T and e1*beta0, gives out ss
54 void blockTriDiagSolve(const std::vector<Eigen::MatrixXd> & alphas,
55  const std::vector<Eigen::MatrixXd> & betas,
56  const Eigen::MatrixXd & beta0, Eigen::MatrixXd & ss,
57  bool & complexValues,
58  const int members) {
59  const int iter = alphas.size();
60  complexValues = false;
61  Eigen::MatrixXd TT = Eigen::MatrixXd::Zero(iter * members, iter * members);
62 
63  for (int ii = 0; ii < iter; ++ii) {
64  TT.block(ii*members, ii*members, members, members) = alphas[ii];
65  if (ii > 0) TT.block((ii-1)*members, ii*members, members, members) = (betas[ii-1]).transpose();
66  if (ii < iter - 1) TT.block((ii+1)*members, ii*members, members, members) = (betas[ii]);
67  }
68 
69  Eigen::MatrixXd e1b0 = Eigen::MatrixXd::Zero(iter * members, members);
70  e1b0.block(0, 0, members, members) = beta0;
71 
72  ss = TT.llt().solve(e1b0);
73 
74  Eigen::VectorXcd eivals = TT.eigenvalues();
75  Eigen::MatrixXd eivalsImag = eivals.imag();
76  complexValues = (eivalsImag.array() != 0.0).any();
77  if (complexValues) {
78  throw eckit::BadValue("T matrix has complex values.");
79  }
80 }
81 
82 
83 } // namespace oops
84 
85 #endif // OOPS_ASSIMILATION_TRIDIAGSOLVE_H_
The namespace for the main oops code.
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)