MFEM  v3.3
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
solvers.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
11 
12 #ifndef MFEM_SOLVERS
13 #define MFEM_SOLVERS
14 
15 #include "../config/config.hpp"
16 #include "operator.hpp"
17 
18 #ifdef MFEM_USE_MPI
19 #include <mpi.h>
20 #endif
21 
22 #ifdef MFEM_USE_SUITESPARSE
23 #include "sparsemat.hpp"
24 #include <umfpack.h>
25 #include <klu.h>
26 #endif
27 
28 namespace mfem
29 {
30 
31 /// Abstract base class for iterative solver
32 class IterativeSolver : public Solver
33 {
34 #ifdef MFEM_USE_MPI
35 private:
36  int dot_prod_type; // 0 - local, 1 - global over 'comm'
37  MPI_Comm comm;
38 #endif
39 
40 protected:
41  const Operator *oper;
43 
45  double rel_tol, abs_tol;
46 
47  // stats
48  mutable int final_iter, converged;
49  mutable double final_norm;
50 
51  double Dot(const Vector &x, const Vector &y) const;
52  double Norm(const Vector &x) const { return sqrt(Dot(x, x)); }
53 
54 public:
56 
57 #ifdef MFEM_USE_MPI
58  IterativeSolver(MPI_Comm _comm);
59 #endif
60 
61  void SetRelTol(double rtol) { rel_tol = rtol; }
62  void SetAbsTol(double atol) { abs_tol = atol; }
63  void SetMaxIter(int max_it) { max_iter = max_it; }
64  void SetPrintLevel(int print_lvl);
65 
66  int GetNumIterations() const { return final_iter; }
67  int GetConverged() const { return converged; }
68  double GetFinalNorm() const { return final_norm; }
69 
70  /// This should be called before SetOperator
71  virtual void SetPreconditioner(Solver &pr);
72 
73  /// Also calls SetOperator for the preconditioner if there is one
74  virtual void SetOperator(const Operator &op);
75 };
76 
77 
78 /// Stationary linear iteration: x <- x + B (b - A x)
79 class SLISolver : public IterativeSolver
80 {
81 protected:
82  mutable Vector r, z;
83 
84  void UpdateVectors();
85 
86 public:
87  SLISolver() { }
88 
89 #ifdef MFEM_USE_MPI
90  SLISolver(MPI_Comm _comm) : IterativeSolver(_comm) { }
91 #endif
92 
93  virtual void SetOperator(const Operator &op)
95 
96  virtual void Mult(const Vector &x, Vector &y) const;
97 };
98 
99 /// Stationary linear iteration. (tolerances are squared)
100 void SLI(const Operator &A, const Vector &b, Vector &x,
101  int print_iter = 0, int max_num_iter = 1000,
102  double RTOLERANCE = 1e-12, double ATOLERANCE = 1e-24);
103 
104 /// Preconditioned stationary linear iteration. (tolerances are squared)
105 void SLI(const Operator &A, Solver &B, const Vector &b, Vector &x,
106  int print_iter = 0, int max_num_iter = 1000,
107  double RTOLERANCE = 1e-12, double ATOLERANCE = 1e-24);
108 
109 
110 /// Conjugate gradient method
111 class CGSolver : public IterativeSolver
112 {
113 protected:
114  mutable Vector r, d, z;
115 
116  void UpdateVectors();
117 
118 public:
119  CGSolver() { }
120 
121 #ifdef MFEM_USE_MPI
122  CGSolver(MPI_Comm _comm) : IterativeSolver(_comm) { }
123 #endif
124 
125  virtual void SetOperator(const Operator &op)
127 
128  virtual void Mult(const Vector &x, Vector &y) const;
129 };
130 
131 /// Conjugate gradient method. (tolerances are squared)
132 void CG(const Operator &A, const Vector &b, Vector &x,
133  int print_iter = 0, int max_num_iter = 1000,
134  double RTOLERANCE = 1e-12, double ATOLERANCE = 1e-24);
135 
136 /// Preconditioned conjugate gradient method. (tolerances are squared)
137 void PCG(const Operator &A, Solver &B, const Vector &b, Vector &x,
138  int print_iter = 0, int max_num_iter = 1000,
139  double RTOLERANCE = 1e-12, double ATOLERANCE = 1e-24);
140 
141 
142 /// GMRES method
144 {
145 protected:
146  int m;
147 
148 public:
149  GMRESSolver() { m = 50; }
150 
151 #ifdef MFEM_USE_MPI
152  GMRESSolver(MPI_Comm _comm) : IterativeSolver(_comm) { m = 50; }
153 #endif
154 
155  void SetKDim(int dim) { m = dim; }
156 
157  virtual void Mult(const Vector &x, Vector &y) const;
158 };
159 
160 /// FGMRES method
162 {
163 protected:
164  int m;
165 
166 public:
167  FGMRESSolver() { m = 50; }
168 
169 #ifdef MFEM_USE_MPI
170  FGMRESSolver(MPI_Comm _comm) : IterativeSolver(_comm) { m = 50; }
171 #endif
172 
173  void SetKDim(int dim) { m = dim; }
174 
175  virtual void Mult(const Vector &x, Vector &y) const;
176 };
177 
178 /// GMRES method. (tolerances are squared)
179 int GMRES(const Operator &A, Vector &x, const Vector &b, Solver &M,
180  int &max_iter, int m, double &tol, double atol, int printit);
181 
182 /// GMRES method. (tolerances are squared)
183 void GMRES(const Operator &A, Solver &B, const Vector &b, Vector &x,
184  int print_iter = 0, int max_num_iter = 1000, int m = 50,
185  double rtol = 1e-12, double atol = 1e-24);
186 
187 
188 /// BiCGSTAB method
190 {
191 protected:
192  mutable Vector p, phat, s, shat, t, v, r, rtilde;
193 
194  void UpdateVectors();
195 
196 public:
198 
199 #ifdef MFEM_USE_MPI
200  BiCGSTABSolver(MPI_Comm _comm) : IterativeSolver(_comm) { }
201 #endif
202 
203  virtual void SetOperator(const Operator &op)
205 
206  virtual void Mult(const Vector &x, Vector &y) const;
207 };
208 
209 /// BiCGSTAB method. (tolerances are squared)
210 int BiCGSTAB(const Operator &A, Vector &x, const Vector &b, Solver &M,
211  int &max_iter, double &tol, double atol, int printit);
212 
213 /// BiCGSTAB method. (tolerances are squared)
214 void BiCGSTAB(const Operator &A, Solver &B, const Vector &b, Vector &x,
215  int print_iter = 0, int max_num_iter = 1000,
216  double rtol = 1e-12, double atol = 1e-24);
217 
218 
219 /// MINRES method
221 {
222 protected:
223  mutable Vector v0, v1, w0, w1, q;
224  mutable Vector u1; // used in the preconditioned version
225 
226 public:
228 
229 #ifdef MFEM_USE_MPI
230  MINRESSolver(MPI_Comm _comm) : IterativeSolver(_comm) { }
231 #endif
232 
233  virtual void SetPreconditioner(Solver &pr)
234  {
236  if (oper) { u1.SetSize(width); }
237  }
238 
239  virtual void SetOperator(const Operator &op);
240 
241  virtual void Mult(const Vector &b, Vector &x) const;
242 };
243 
244 /// MINRES method without preconditioner. (tolerances are squared)
245 void MINRES(const Operator &A, const Vector &b, Vector &x, int print_it = 0,
246  int max_it = 1000, double rtol = 1e-12, double atol = 1e-24);
247 
248 /// MINRES method with preconditioner. (tolerances are squared)
249 void MINRES(const Operator &A, Solver &B, const Vector &b, Vector &x,
250  int print_it = 0, int max_it = 1000,
251  double rtol = 1e-12, double atol = 1e-24);
252 
253 
254 /// Newton's method for solving F(x)=b for a given operator F.
255 /** The method GetGradient() must be implemented for the operator F.
256  The preconditioner is used (in non-iterative mode) to evaluate
257  the action of the inverse gradient of the operator. */
259 {
260 protected:
261  mutable Vector r, c;
262 
263 public:
265 
266 #ifdef MFEM_USE_MPI
267  NewtonSolver(MPI_Comm _comm) : IterativeSolver(_comm) { }
268 #endif
269  virtual void SetOperator(const Operator &op);
270 
271  /// Set the linear solver for inverting the Jacobian.
272  /** This method is equivalent to calling SetPreconditioner(). */
273  virtual void SetSolver(Solver &solver) { prec = &solver; }
274 
275  /// Solve the nonlinear system with right-hand side @a b.
276  /** If `b.Size() != Height()`, then @a b is assumed to be zero. */
277  virtual void Mult(const Vector &b, Vector &x) const;
278 };
279 
280 
281 /** Adaptive restarted GMRES.
282  m_max and m_min(=1) are the maximal and minimal restart parameters.
283  m_step(=1) is the step to use for going from m_max and m_min.
284  cf(=0.4) is a desired convergence factor. */
285 int aGMRES(const Operator &A, Vector &x, const Vector &b,
286  const Operator &M, int &max_iter,
287  int m_max, int m_min, int m_step, double cf,
288  double &tol, double &atol, int printit);
289 
290 
291 /** SLBQP: (S)ingle (L)inearly Constrained with (B)ounds (Q)uadratic (P)rogram
292 
293  minimize 1/2 ||x - x_t||^2, subject to:
294  lo_i <= x_i <= hi_i
295  sum_i w_i x_i = a
296 */
298 {
299 protected:
301  double a;
302 
303  /// Solve QP at fixed lambda
304  inline double solve(double l, const Vector &xt, Vector &x, int &nclip) const
305  {
306  add(xt, l, w, x);
307  x.median(lo,hi);
308  nclip++;
309  return Dot(w,x)-a;
310  }
311 
312  inline void print_iteration(int it, double r, double l) const;
313 
314 public:
316 
317 #ifdef MFEM_USE_MPI
318  SLBQPOptimizer(MPI_Comm _comm) : IterativeSolver(_comm) {}
319 #endif
320 
321  void SetBounds(const Vector &_lo, const Vector &_hi);
322  void SetLinearConstraint(const Vector &_w, double _a);
323 
324  // For this problem type, we let the target values play the role of the
325  // initial vector xt, from which the operator generates the optimal vector x.
326  virtual void Mult(const Vector &xt, Vector &x) const;
327 
328  /// These are not currently meaningful for this solver and will error out.
329  virtual void SetPreconditioner(Solver &pr);
330  virtual void SetOperator(const Operator &op);
331 };
332 
333 
334 #ifdef MFEM_USE_SUITESPARSE
335 
336 /// Direct sparse solver using UMFPACK
337 class UMFPackSolver : public Solver
338 {
339 protected:
342  void *Numeric;
343  SuiteSparse_long *AI, *AJ;
344 
345  void Init();
346 
347 public:
348  double Control[UMFPACK_CONTROL];
349  mutable double Info[UMFPACK_INFO];
350 
351  /** @brief For larger matrices, if the solver fails, set the parameter @a
352  _use_long_ints = true. */
353  UMFPackSolver(bool _use_long_ints = false)
354  : use_long_ints(_use_long_ints) { Init(); }
355  /** @brief Factorize the given SparseMatrix using the defaults. For larger
356  matrices, if the solver fails, set the parameter @a _use_long_ints =
357  true. */
358  UMFPackSolver(SparseMatrix &A, bool _use_long_ints = false)
359  : use_long_ints(_use_long_ints) { Init(); SetOperator(A); }
360 
361  /** @brief Factorize the given Operator @a op which must be a SparseMatrix.
362 
363  The factorization uses the parameters set in the #Control data member.
364  @note This method calls SparseMatrix::SortColumnIndices() with @a op,
365  modifying the matrix if the column indices are not already sorted. */
366  virtual void SetOperator(const Operator &op);
367 
368  /// Set the print level field in the #Control data member.
369  void SetPrintLevel(int print_lvl) { Control[UMFPACK_PRL] = print_lvl; }
370 
371  virtual void Mult(const Vector &b, Vector &x) const;
372  virtual void MultTranspose(const Vector &b, Vector &x) const;
373 
374  virtual ~UMFPackSolver();
375 };
376 
377 /// Direct sparse solver using KLU
378 class KLUSolver : public Solver
379 {
380 protected:
382  klu_symbolic *Symbolic;
383  klu_numeric *Numeric;
384 
385  void Init();
386 
387 public:
389  : mat(0),Symbolic(0),Numeric(0)
390  { Init(); }
392  : mat(0),Symbolic(0),Numeric(0)
393  { Init(); SetOperator(A); }
394 
395  // Works on sparse matrices only; calls SparseMatrix::SortColumnIndices().
396  virtual void SetOperator(const Operator &op);
397 
398  virtual void Mult(const Vector &b, Vector &x) const;
399  virtual void MultTranspose(const Vector &b, Vector &x) const;
400 
401  virtual ~KLUSolver();
402 
403  mutable klu_common Common;
404 };
405 
406 #endif // MFEM_USE_SUITESPARSE
407 
408 }
409 
410 #endif // MFEM_SOLVERS
Conjugate gradient method.
Definition: solvers.hpp:111
virtual void MultTranspose(const Vector &b, Vector &x) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
Definition: solvers.cpp:1806
double Info[UMFPACK_INFO]
Definition: solvers.hpp:349
SparseMatrix * mat
Definition: solvers.hpp:341
int GetNumIterations() const
Definition: solvers.hpp:66
const Operator * oper
Definition: solvers.hpp:41
SuiteSparse_long * AI
Definition: solvers.hpp:343
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:310
NewtonSolver(MPI_Comm _comm)
Definition: solvers.hpp:267
FGMRES method.
Definition: solvers.hpp:161
Direct sparse solver using KLU.
Definition: solvers.hpp:378
virtual void SetSolver(Solver &solver)
Set the linear solver for inverting the Jacobian.
Definition: solvers.hpp:273
virtual void MultTranspose(const Vector &b, Vector &x) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
Definition: solvers.cpp:1907
UMFPackSolver(SparseMatrix &A, bool _use_long_ints=false)
Factorize the given SparseMatrix using the defaults. For larger matrices, if the solver fails...
Definition: solvers.hpp:358
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: solvers.cpp:294
MINRES method.
Definition: solvers.hpp:220
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: solvers.cpp:116
void add(const Vector &v1, const Vector &v2, Vector &v)
Definition: vector.cpp:264
double Norm(const Vector &x) const
Definition: solvers.hpp:52
double GetFinalNorm() const
Definition: solvers.hpp:68
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: solvers.cpp:879
int BiCGSTAB(const Operator &A, Vector &x, const Vector &b, Solver &M, int &max_iter, double &tol, double atol, int printit)
BiCGSTAB method. (tolerances are squared)
Definition: solvers.cpp:1008
void SetKDim(int dim)
Definition: solvers.hpp:173
FGMRESSolver(MPI_Comm _comm)
Definition: solvers.hpp:170
Direct sparse solver using UMFPACK.
Definition: solvers.hpp:337
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: solvers.cpp:521
UMFPackSolver(bool _use_long_ints=false)
For larger matrices, if the solver fails, set the parameter _use_long_ints = true.
Definition: solvers.hpp:353
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:1045
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: solvers.cpp:701
void SetLinearConstraint(const Vector &_w, double _a)
Definition: solvers.cpp:1462
int dim
Definition: ex3.cpp:47
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.hpp:233
void SetPrintLevel(int print_lvl)
Definition: solvers.cpp:71
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.cpp:98
double solve(double l, const Vector &xt, Vector &x, int &nclip) const
Solve QP at fixed lambda.
Definition: solvers.hpp:304
Data type sparse matrix.
Definition: sparsemat.hpp:38
void MINRES(const Operator &A, const Vector &b, Vector &x, int print_it, int max_it, double rtol, double atol)
MINRES method without preconditioner. (tolerances are squared)
Definition: solvers.cpp:1203
virtual ~UMFPackSolver()
Definition: solvers.cpp:1840
void SetMaxIter(int max_it)
Definition: solvers.hpp:63
void median(const Vector &lo, const Vector &hi)
v = median(v,lo,hi) entrywise. Implementation assumes lo &lt;= hi.
Definition: vector.cpp:445
SuiteSparse_long * AJ
Definition: solvers.hpp:343
void SetKDim(int dim)
Definition: solvers.hpp:155
Newton&#39;s method for solving F(x)=b for a given operator F.
Definition: solvers.hpp:258
void CG(const Operator &A, const Vector &b, Vector &x, int print_iter, int max_num_iter, double RTOLERANCE, double ATOLERANCE)
Conjugate gradient method. (tolerances are squared)
Definition: solvers.cpp:443
void print_iteration(int it, double r, double l) const
Definition: solvers.cpp:1480
void PCG(const Operator &A, Solver &B, const Vector &b, Vector &x, int print_iter, int max_num_iter, double RTOLERANCE, double ATOLERANCE)
Preconditioned conjugate gradient method. (tolerances are squared)
Definition: solvers.cpp:456
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: solvers.cpp:1862
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.cpp:1031
double Dot(const Vector &x, const Vector &y) const
Definition: solvers.cpp:50
int GetConverged() const
Definition: solvers.hpp:67
BiCGSTABSolver(MPI_Comm _comm)
Definition: solvers.hpp:200
klu_symbolic * Symbolic
Definition: solvers.hpp:382
SLBQPOptimizer(MPI_Comm _comm)
Definition: solvers.hpp:318
virtual void Mult(const Vector &b, Vector &x) const
Solve the nonlinear system with right-hand side b.
Definition: solvers.cpp:1240
CGSolver(MPI_Comm _comm)
Definition: solvers.hpp:122
void SetPrintLevel(int print_lvl)
Set the print level field in the Control data member.
Definition: solvers.hpp:369
Stationary linear iteration: x &lt;- x + B (b - A x)
Definition: solvers.hpp:79
void SetAbsTol(double atol)
Definition: solvers.hpp:62
void SetRelTol(double rtol)
Definition: solvers.hpp:61
double Control[UMFPACK_CONTROL]
Definition: solvers.hpp:348
Abstract base class for iterative solver.
Definition: solvers.hpp:32
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:1893
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.cpp:1474
SLISolver(MPI_Comm _comm)
Definition: solvers.hpp:90
GMRES method.
Definition: solvers.hpp:143
GMRESSolver(MPI_Comm _comm)
Definition: solvers.hpp:152
void SetBounds(const Vector &_lo, const Vector &_hi)
Definition: solvers.cpp:1456
void UpdateVectors()
Definition: solvers.cpp:110
klu_common Common
Definition: solvers.hpp:403
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:203
virtual void Mult(const Vector &xt, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:1487
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:125
Vector data type.
Definition: vector.hpp:36
MINRESSolver(MPI_Comm _comm)
Definition: solvers.hpp:230
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:92
KLUSolver(SparseMatrix &A)
Definition: solvers.hpp:391
int aGMRES(const Operator &A, Vector &x, const Vector &b, const Operator &M, int &max_iter, int m_max, int m_min, int m_step, double cf, double &tol, double &atol, int printit)
Definition: solvers.cpp:1304
virtual ~KLUSolver()
Definition: solvers.cpp:1921
void UpdateVectors()
Definition: solvers.cpp:287
BiCGSTAB method.
Definition: solvers.hpp:189
SparseMatrix * mat
Definition: solvers.hpp:381
Base class for solvers.
Definition: operator.hpp:257
virtual void SetPreconditioner(Solver &pr)
These are not currently meaningful for this solver and will error out.
Definition: solvers.cpp:1468
Abstract operator.
Definition: operator.hpp:21
klu_numeric * Numeric
Definition: solvers.hpp:383
int width
Dimension of the input / number of columns in the matrix.
Definition: operator.hpp:25
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:1774
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:93
void SLI(const Operator &A, const Vector &b, Vector &x, int print_iter, int max_num_iter, double RTOLERANCE, double ATOLERANCE)
Stationary linear iteration. (tolerances are squared)
Definition: solvers.cpp:259
int GMRES(const Operator &A, Vector &x, const Vector &b, Solver &M, int &max_iter, int m, double &tol, double atol, int printit)
GMRES method. (tolerances are squared)
Definition: solvers.cpp:843
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.cpp:1229
virtual void SetOperator(const Operator &op)
Factorize the given Operator op which must be a SparseMatrix.
Definition: solvers.cpp:1674