OpenMEEG
gmres.h
Go to the documentation of this file.
1 /*
2 Project Name : OpenMEEG
3 
4 © INRIA and ENPC (contributors: Geoffray ADDE, Maureen CLERC, Alexandre
5 GRAMFORT, Renaud KERIVEN, Jan KYBIC, Perrine LANDREAU, Théodore PAPADOPOULO,
6 Emmanuel OLIVI
7 Maureen.Clerc.AT.inria.fr, keriven.AT.certis.enpc.fr,
8 kybic.AT.fel.cvut.cz, papadop.AT.inria.fr)
9 
10 The OpenMEEG software is a C++ package for solving the forward/inverse
11 problems of electroencephalography and magnetoencephalography.
12 
13 This software is governed by the CeCILL-B license under French law and
14 abiding by the rules of distribution of free software. You can use,
15 modify and/ or redistribute the software under the terms of the CeCILL-B
16 license as circulated by CEA, CNRS and INRIA at the following URL
17 "http://www.cecill.info".
18 
19 As a counterpart to the access to the source code and rights to copy,
20 modify and redistribute granted by the license, users are provided only
21 with a limited warranty and the software's authors, the holders of the
22 economic rights, and the successive licensors have only limited
23 liability.
24 
25 In this respect, the user's attention is drawn to the risks associated
26 with loading, using, modifying and/or developing or reproducing the
27 software by the user in light of its specific status of free software,
28 that may mean that it is complicated to manipulate, and that also
29 therefore means that it is reserved for developers and experienced
30 professionals having in-depth computer knowledge. Users are therefore
31 encouraged to load and test the software's suitability as regards their
32 requirements in conditions enabling the security of their systems and/or
33 data to be ensured and, more generally, to use and operate it in the
34 same conditions as regards security.
35 
36 The fact that you are presently reading this means that you have had
37 knowledge of the CeCILL-B license and that you accept its terms.
38 */
39 
40 #pragma once
41 
42 #include "vector.h"
43 #include "matrix.h"
44 #include "sparse_matrix.h"
45 
46 #include <OpenMEEG_Export.h>
47 
48 namespace OpenMEEG {
49 
50  // ===================================
51  // = Define a Jacobi preconditionner =
52  // ===================================
53  template <typename M>
54  class Jacobi {
55  public:
56  Jacobi (const M& m): J(m.nlin()) {
57  for ( unsigned i = 0; i < m.nlin(); ++i) {
58  J(i, i) = 1.0 / m(i,i);
59  }
60  }
61 
62  Vector operator()(const Vector& g) const {
63  return J*g;
64  }
65 
66  ~Jacobi () {};
67  private:
68  SparseMatrix J; // diagonal
69  };
70 
71  // =========================
72  // = Define a GMRes solver =
73  // =========================
74 
75  void GeneratePlaneRotation(double &dx, double &dy, double &cs, double &sn)
76  {
77  if (dy == 0.0) {
78  cs = 1.0;
79  sn = 0.0;
80  } else if (std::abs(dy) > std::abs(dx)) {
81  double temp = dx / dy;
82  sn = 1.0 / sqrt( 1.0 + temp*temp );
83  cs = temp * sn;
84  } else {
85  double temp = dy / dx;
86  cs = 1.0 / sqrt( 1.0 + temp*temp );
87  sn = temp * cs;
88  }
89  }
90 
91  void ApplyPlaneRotation(double &dx, double &dy, double &cs, double &sn)
92  {
93  double temp = cs * dx + sn * dy;
94  dy = -sn * dx + cs * dy;
95  dx = temp;
96  }
97 
98  template<class T>
99  void Update(Vector &x, int k, T &h, Vector &s, Vector v[])
100  {
101  Vector y(s);
102  // Backsolve:
103  for (int i = k; i >= 0; i--) {
104  y(i) /= h(i,i);
105  for (int j = i - 1; j >= 0; j--)
106  y(j) -= h(j,i) * y(i);
107  }
108  for (int j = 0; j <= k; j++) {
109  x += v[j] * y(j);
110  }
111  }
112 
113  // code taken from http://www.netlib.org/templates/cpp/gmres.h and modified
114  template<class T,class P> // T should be a linear operator, and P a preconditionner
115  unsigned GMRes(const T& A, const P& M, Vector &x, const Vector& b, int max_iter, double tol,unsigned m) {
116 
117  // m is the size of the Krylov subspace, if m<A.nlin(), it is a restarted GMRes (for saving memory)
118  Matrix H(m+1,m);
119  x.set(0.0);
120 
121  double resid;
122  int i, j = 1, k;
123  Vector s(m+1), cs(m+1), sn(m+1), w;
124 
125  double normb = (M(b)).norm();//(M*b).norm()
126  Vector r = M(b-A*x);//M.solve(b - A * x);
127  double beta = r.norm();
128 
129  if (normb == 0.0)
130  normb = 1;
131 
132  if ((resid = r.norm() / normb) <= tol) {
133  tol = resid;
134  max_iter = 0;
135  return 0;
136  }
137  Vector *v = new Vector[m+1];
138 
139  while (j <= max_iter) {
140  v[0] = r * (1.0 / beta);
141  s.set(0.0);
142  s(0) = beta;
143 
144  for (i = 0; i < m && j <= max_iter; i++, j++) {
145  w = M(A*v[i]); //M.solve(A * v[i]);
146  for (k = 0; k <= i; k++) {
147  H(k, i) = w*v[k];
148  w -= H(k, i) * v[k];
149  }
150  H(i+1, i) = w.norm();
151  v[i+1] = (w / H(i+1, i));
152 
153  for (k = 0; k < i; k++)
154  ApplyPlaneRotation(H(k,i), H(k+1,i), cs(k), sn(k));
155 
156  GeneratePlaneRotation(H(i,i), H(i+1,i), cs(i), sn(i));
157  ApplyPlaneRotation(H(i,i), H(i+1,i), cs(i), sn(i));
158  ApplyPlaneRotation(s(i), s(i+1), cs(i), sn(i));
159 
160  if ((resid = std::abs(s(i+1)) / normb) < tol) {
161  Update(x, i, H, s, v);
162  tol = resid;
163  max_iter = j;
164  // std::cout<<max_iter <<std::endl;
165  delete [] v;
166  return 0;
167  }
168  }
169  Update(x, i - 1, H, s, v);
170  r = M(b-A*x);//M.solve(b - A * x);
171  beta = r.norm();
172  if ((resid = beta / normb) < tol) {
173  tol = resid;
174  max_iter = j;
175  // std::cout<<max_iter <<std::endl;
176  delete [] v;
177  return 0;
178  }
179  }
180 
181  tol = resid;
182  delete [] v;
183  return 1;
184  }
185 }
sparse_matrix.h
OpenMEEG::Jacobi
Definition: gmres.h:54
OpenMEEG::Jacobi::Jacobi
Jacobi(const M &m)
Definition: gmres.h:56
OpenMEEG::Jacobi::operator()
Vector operator()(const Vector &g) const
Definition: gmres.h:62
OpenMEEG::Update
void Update(Vector &x, int k, T &h, Vector &s, Vector v[])
Definition: gmres.h:99
OpenMEEG::Vector::set
void set(double x)
OpenMEEG
Definition: analytics.h:45
matrix.h
OpenMEEG::SparseMatrix
Definition: sparse_matrix.h:62
OpenMEEG::ApplyPlaneRotation
void ApplyPlaneRotation(double &dx, double &dy, double &cs, double &sn)
Definition: gmres.h:91
OpenMEEG::Vector::norm
double norm() const
Definition: vector.h:215
OpenMEEG::Jacobi::~Jacobi
~Jacobi()
Definition: gmres.h:66
OpenMEEG::Vector
Definition: vector.h:56
OpenMEEG::Matrix
Matrix class.
Definition: matrix.h:61
vector.h
OpenMEEG::GMRes
unsigned GMRes(const T &A, const P &M, Vector &x, const Vector &b, int max_iter, double tol, unsigned m)
Definition: gmres.h:115
OpenMEEG::Jacobi::J
SparseMatrix J
Definition: gmres.h:66
OpenMEEG::GeneratePlaneRotation
void GeneratePlaneRotation(double &dx, double &dy, double &cs, double &sn)
Definition: gmres.h:75