MFEM  v4.1.0
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-2020, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #ifndef MFEM_SOLVERS
13 #define MFEM_SOLVERS
14 
15 #include "../config/config.hpp"
16 #include "densemat.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 class BilinearForm;
32 
33 /// Abstract base class for iterative solver
34 class IterativeSolver : public Solver
35 {
36 #ifdef MFEM_USE_MPI
37 private:
38  int dot_prod_type; // 0 - local, 1 - global over 'comm'
39  MPI_Comm comm;
40 #endif
41 
42 protected:
43  const Operator *oper;
45 
47  double rel_tol, abs_tol;
48 
49  // stats
50  mutable int final_iter, converged;
51  mutable double final_norm;
52 
53  double Dot(const Vector &x, const Vector &y) const;
54  double Norm(const Vector &x) const { return sqrt(Dot(x, x)); }
55 
56 public:
58 
59 #ifdef MFEM_USE_MPI
60  IterativeSolver(MPI_Comm _comm);
61 #endif
62 
63  void SetRelTol(double rtol) { rel_tol = rtol; }
64  void SetAbsTol(double atol) { abs_tol = atol; }
65  void SetMaxIter(int max_it) { max_iter = max_it; }
66  void SetPrintLevel(int print_lvl);
67 
68  int GetNumIterations() const { return final_iter; }
69  int GetConverged() const { return converged; }
70  double GetFinalNorm() const { return final_norm; }
71 
72  /// This should be called before SetOperator
73  virtual void SetPreconditioner(Solver &pr);
74 
75  /// Also calls SetOperator for the preconditioner if there is one
76  virtual void SetOperator(const Operator &op);
77 };
78 
79 
80 /// Jacobi smoothing for a given bilinear form (no matrix necessary).
81 /** Useful with tensorized, partially assembled operators. Can also be defined
82  by given diagonal vector. This is basic Jacobi iteration; for tolerances,
83  iteration control, etc. wrap with SLISolver. */
85 {
86 public:
87  /** Setup a Jacobi smoother with the diagonal of @a a obtained by calling
88  a.AssembleDiagonal(). It is assumed that the underlying operator acts as
89  the identity on entries in ess_tdof_list, corresponding to (assembled)
90  DIAG_ONE policy or ConstratinedOperator in the matrix-free setting. */
92  const Array<int> &ess_tdof_list,
93  const double damping=1.0);
94 
95  /** Application is by the *inverse* of the given vector. It is assumed that
96  the underlying operator acts as the identity on entries in ess_tdof_list,
97  corresponding to (assembled) DIAG_ONE policy or ConstratinedOperator in
98  the matrix-free setting. */
100  const Array<int> &ess_tdof_list,
101  const double damping=1.0);
103 
104  void Mult(const Vector &x, Vector &y) const;
105  void SetOperator(const Operator &op) { oper = &op; }
106  void Setup(const Vector &diag);
107 
108 private:
109  const int N;
110  Vector dinv;
111  const double damping;
112  const Array<int> &ess_tdof_list;
113  mutable Vector residual;
114 
115  const Operator *oper;
116 };
117 
118 
119 /// Stationary linear iteration: x <- x + B (b - A x)
121 {
122 protected:
123  mutable Vector r, z;
124 
125  void UpdateVectors();
126 
127 public:
128  SLISolver() { }
129 
130 #ifdef MFEM_USE_MPI
131  SLISolver(MPI_Comm _comm) : IterativeSolver(_comm) { }
132 #endif
133 
134  virtual void SetOperator(const Operator &op)
136 
137  virtual void Mult(const Vector &b, Vector &x) const;
138 };
139 
140 /// Stationary linear iteration. (tolerances are squared)
141 void SLI(const Operator &A, const Vector &b, Vector &x,
142  int print_iter = 0, int max_num_iter = 1000,
143  double RTOLERANCE = 1e-12, double ATOLERANCE = 1e-24);
144 
145 /// Preconditioned stationary linear iteration. (tolerances are squared)
146 void SLI(const Operator &A, Solver &B, const Vector &b, Vector &x,
147  int print_iter = 0, int max_num_iter = 1000,
148  double RTOLERANCE = 1e-12, double ATOLERANCE = 1e-24);
149 
150 
151 /// Conjugate gradient method
152 class CGSolver : public IterativeSolver
153 {
154 protected:
155  mutable Vector r, d, z;
156 
157  void UpdateVectors();
158 
159 public:
160  CGSolver() { }
161 
162 #ifdef MFEM_USE_MPI
163  CGSolver(MPI_Comm _comm) : IterativeSolver(_comm) { }
164 #endif
165 
166  virtual void SetOperator(const Operator &op)
168 
169  virtual void Mult(const Vector &b, Vector &x) const;
170 };
171 
172 /// Conjugate gradient method. (tolerances are squared)
173 void CG(const Operator &A, const Vector &b, Vector &x,
174  int print_iter = 0, int max_num_iter = 1000,
175  double RTOLERANCE = 1e-12, double ATOLERANCE = 1e-24);
176 
177 /// Preconditioned conjugate gradient method. (tolerances are squared)
178 void PCG(const Operator &A, Solver &B, const Vector &b, Vector &x,
179  int print_iter = 0, int max_num_iter = 1000,
180  double RTOLERANCE = 1e-12, double ATOLERANCE = 1e-24);
181 
182 
183 /// GMRES method
185 {
186 protected:
187  int m; // see SetKDim()
188 
189 public:
190  GMRESSolver() { m = 50; }
191 
192 #ifdef MFEM_USE_MPI
193  GMRESSolver(MPI_Comm _comm) : IterativeSolver(_comm) { m = 50; }
194 #endif
195 
196  /// Set the number of iteration to perform between restarts, default is 50.
197  void SetKDim(int dim) { m = dim; }
198 
199  virtual void Mult(const Vector &b, Vector &x) const;
200 };
201 
202 /// FGMRES method
204 {
205 protected:
206  int m;
207 
208 public:
209  FGMRESSolver() { m = 50; }
210 
211 #ifdef MFEM_USE_MPI
212  FGMRESSolver(MPI_Comm _comm) : IterativeSolver(_comm) { m = 50; }
213 #endif
214 
215  void SetKDim(int dim) { m = dim; }
216 
217  virtual void Mult(const Vector &b, Vector &x) const;
218 };
219 
220 /// GMRES method. (tolerances are squared)
221 int GMRES(const Operator &A, Vector &x, const Vector &b, Solver &M,
222  int &max_iter, int m, double &tol, double atol, int printit);
223 
224 /// GMRES method. (tolerances are squared)
225 void GMRES(const Operator &A, Solver &B, const Vector &b, Vector &x,
226  int print_iter = 0, int max_num_iter = 1000, int m = 50,
227  double rtol = 1e-12, double atol = 1e-24);
228 
229 
230 /// BiCGSTAB method
232 {
233 protected:
234  mutable Vector p, phat, s, shat, t, v, r, rtilde;
235 
236  void UpdateVectors();
237 
238 public:
240 
241 #ifdef MFEM_USE_MPI
242  BiCGSTABSolver(MPI_Comm _comm) : IterativeSolver(_comm) { }
243 #endif
244 
245  virtual void SetOperator(const Operator &op)
247 
248  virtual void Mult(const Vector &b, Vector &x) const;
249 };
250 
251 /// BiCGSTAB method. (tolerances are squared)
252 int BiCGSTAB(const Operator &A, Vector &x, const Vector &b, Solver &M,
253  int &max_iter, double &tol, double atol, int printit);
254 
255 /// BiCGSTAB method. (tolerances are squared)
256 void BiCGSTAB(const Operator &A, Solver &B, const Vector &b, Vector &x,
257  int print_iter = 0, int max_num_iter = 1000,
258  double rtol = 1e-12, double atol = 1e-24);
259 
260 
261 /// MINRES method
263 {
264 protected:
265  mutable Vector v0, v1, w0, w1, q;
266  mutable Vector u1; // used in the preconditioned version
267 
268 public:
270 
271 #ifdef MFEM_USE_MPI
272  MINRESSolver(MPI_Comm _comm) : IterativeSolver(_comm) { }
273 #endif
274 
275  virtual void SetPreconditioner(Solver &pr)
276  {
278  if (oper) { u1.SetSize(width); }
279  }
280 
281  virtual void SetOperator(const Operator &op);
282 
283  virtual void Mult(const Vector &b, Vector &x) const;
284 };
285 
286 /// MINRES method without preconditioner. (tolerances are squared)
287 void MINRES(const Operator &A, const Vector &b, Vector &x, int print_it = 0,
288  int max_it = 1000, double rtol = 1e-12, double atol = 1e-24);
289 
290 /// MINRES method with preconditioner. (tolerances are squared)
291 void MINRES(const Operator &A, Solver &B, const Vector &b, Vector &x,
292  int print_it = 0, int max_it = 1000,
293  double rtol = 1e-12, double atol = 1e-24);
294 
295 
296 /// Newton's method for solving F(x)=b for a given operator F.
297 /** The method GetGradient() must be implemented for the operator F.
298  The preconditioner is used (in non-iterative mode) to evaluate
299  the action of the inverse gradient of the operator. */
301 {
302 protected:
303  mutable Vector r, c;
304 
305 public:
307 
308 #ifdef MFEM_USE_MPI
309  NewtonSolver(MPI_Comm _comm) : IterativeSolver(_comm) { }
310 #endif
311  virtual void SetOperator(const Operator &op);
312 
313  /// Set the linear solver for inverting the Jacobian.
314  /** This method is equivalent to calling SetPreconditioner(). */
315  virtual void SetSolver(Solver &solver) { prec = &solver; }
316 
317  /// Solve the nonlinear system with right-hand side @a b.
318  /** If `b.Size() != Height()`, then @a b is assumed to be zero. */
319  virtual void Mult(const Vector &b, Vector &x) const;
320 
321  /** @brief This method can be overloaded in derived classes to implement line
322  search algorithms. */
323  /** The base class implementation (NewtonSolver) simply returns 1. A return
324  value of 0 indicates a failure, interrupting the Newton iteration. */
325  virtual double ComputeScalingFactor(const Vector &x, const Vector &b) const
326  { return 1.0; }
327 
328  /** @brief This method can be overloaded in derived classes to perform
329  computations that need knowledge of the newest Newton state. */
330  virtual void ProcessNewState(const Vector &x) const { }
331 };
332 
333 /** Adaptive restarted GMRES.
334  m_max and m_min(=1) are the maximal and minimal restart parameters.
335  m_step(=1) is the step to use for going from m_max and m_min.
336  cf(=0.4) is a desired convergence factor. */
337 int aGMRES(const Operator &A, Vector &x, const Vector &b,
338  const Operator &M, int &max_iter,
339  int m_max, int m_min, int m_step, double cf,
340  double &tol, double &atol, int printit);
341 
342 
343 /** Defines operators and constraints for the following optimization problem:
344  *
345  * Find x that minimizes the objective function F(x), subject to
346  * C(x) = c_e,
347  * d_lo <= D(x) <= d_hi,
348  * x_lo <= x <= x_hi.
349  *
350  * The operators F, C, D must take input of the same size (same width).
351  * Gradients of F, C, D might be needed, depending on the OptimizationSolver.
352  * When used with Hiop, gradients of C and D must be DenseMatrices.
353  * F always returns a scalar value, see CalcObjective(), CalcObjectiveGrad().
354  * C and D can have arbitrary heights.
355  * C and D can be NULL, meaning that their constraints are not used.
356  *
357  * When used in parallel, all Vectors are assumed to be true dof vectors, and
358  * the operators are expected to be defined for tdof vectors. */
360 {
361 protected:
362  /// Not owned, some can remain unused (NULL).
363  const Operator *C, *D;
364  const Vector *c_e, *d_lo, *d_hi, *x_lo, *x_hi;
365 
366 public:
367  const int input_size;
368 
369  /// In parallel, insize is the number of the local true dofs.
370  OptimizationProblem(int insize, const Operator *C_, const Operator *D_);
371 
372  /// Objective F(x). In parallel, the result should be reduced over tasks.
373  virtual double CalcObjective(const Vector &x) const = 0;
374  /// The result grad is expected to enter with the correct size.
375  virtual void CalcObjectiveGrad(const Vector &x, Vector &grad) const
376  { MFEM_ABORT("The objective gradient is not implemented."); }
377 
378  void SetEqualityConstraint(const Vector &c);
379  void SetInequalityConstraint(const Vector &dl, const Vector &dh);
380  void SetSolutionBounds(const Vector &xl, const Vector &xh);
381 
382  const Operator *GetC() const { return C; }
383  const Operator *GetD() const { return D; }
384  const Vector *GetEqualityVec() const { return c_e; }
385  const Vector *GetInequalityVec_Lo() const { return d_lo; }
386  const Vector *GetInequalityVec_Hi() const { return d_hi; }
387  const Vector *GetBoundsVec_Lo() const { return x_lo; }
388  const Vector *GetBoundsVec_Hi() const { return x_hi; }
389 
390  int GetNumConstraints() const;
391 };
392 
393 /// Abstract solver for OptimizationProblems.
395 {
396 protected:
398 
399 public:
401 #ifdef MFEM_USE_MPI
402  OptimizationSolver(MPI_Comm _comm): IterativeSolver(_comm), problem(NULL) { }
403 #endif
404  virtual ~OptimizationSolver() { }
405 
406  /** This function is virtual as solvers might need to perform some initial
407  * actions (e.g. validation) with the OptimizationProblem. */
408  virtual void SetOptimizationProblem(const OptimizationProblem &prob)
409  { problem = &prob; }
410 
411  virtual void Mult(const Vector &xt, Vector &x) const = 0;
412 
413  virtual void SetPreconditioner(Solver &pr)
414  { MFEM_ABORT("Not meaningful for this solver."); }
415  virtual void SetOperator(const Operator &op)
416  { MFEM_ABORT("Not meaningful for this solver."); }
417 };
418 
419 /** SLBQP optimizer:
420  * (S)ingle (L)inearly Constrained with (B)ounds (Q)uadratic (P)rogram
421  *
422  * Minimize || x-x_t ||, subject to
423  * sum w_i x_i = a,
424  * x_lo <= x <= x_hi.
425  */
427 {
428 protected:
430  double a;
431 
432  /// Solve QP at fixed lambda
433  inline double solve(double l, const Vector &xt, Vector &x, int &nclip) const
434  {
435  add(xt, l, w, x);
436  if (problem == NULL) { x.median(lo,hi); }
437  else
438  {
441  }
442  nclip++;
443  if (problem == NULL) { return Dot(w, x) - a; }
444  else
445  {
446  Vector c(1);
447  // Includes parallel communication.
448  problem->GetC()->Mult(x, c);
449 
450  return c(0) - (*problem->GetEqualityVec())(0);
451  }
452  }
453 
454  inline void print_iteration(int it, double r, double l) const;
455 
456 public:
458 
459 #ifdef MFEM_USE_MPI
460  SLBQPOptimizer(MPI_Comm _comm) : OptimizationSolver(_comm) { }
461 #endif
462 
463  /** Setting an OptimizationProblem will overwrite the Vectors given by
464  * SetBounds and SetLinearConstraint. The objective function remains
465  * unchanged. */
466  virtual void SetOptimizationProblem(const OptimizationProblem &prob);
467 
468  void SetBounds(const Vector &_lo, const Vector &_hi);
469  void SetLinearConstraint(const Vector &_w, double _a);
470 
471  /** We let the target values play the role of the initial vector xt, from
472  * which the operator generates the optimal vector x. */
473  virtual void Mult(const Vector &xt, Vector &x) const;
474 };
475 
476 /** Block ILU solver:
477  * Performs a block ILU(k) approximate factorization with specified block
478  * size. Currently only k=0 is supported. This is useful as a preconditioner
479  * for DG-type discretizations, where the system matrix has a natural
480  * (elemental) block structure.
481  *
482  * In the case of DG discretizations, the block size should usually be set to
483  * either ndofs_per_element or vdim*ndofs_per_element (if the finite element
484  * space has Ordering::byVDIM). The block size must evenly divide the size of
485  * the matrix.
486  *
487  * Renumbering the blocks is also supported by specifying a reordering method.
488  * Currently greedy minimum discarded fill ordering and no reordering are
489  * supported. Renumbering the blocks can lead to a much better approximate
490  * factorization.
491  */
492 class BlockILU : public Solver
493 {
494 public:
495 
496  /// The reordering method used by the BlockILU factorization.
497  enum class Reordering
498  {
500  NONE
501  };
502 
503  /** Create an "empty" BlockILU solver. SetOperator must be called later to
504  * actually form the factorization
505  */
506  BlockILU(int block_size_,
508  int k_fill_ = 0);
509 
510  /** Create a block ILU approximate factorization for the matrix @a op.
511  * @a op should be of type either SparseMatrix or HypreParMatrix. In the
512  * case that @a op is a HypreParMatrix, the ILU factorization is performed
513  * on the diagonal blocks of the parallel decomposition.
514  */
515  BlockILU(Operator &op, int block_size_ = 1,
517  int k_fill_ = 0);
518 
519  /** Perform the block ILU factorization for the matrix @a op.
520  * As in the constructor, @a op must either be a SparseMatrix or
521  * HypreParMatrix
522  */
523  void SetOperator(const Operator &op);
524 
525  /// Solve the system `LUx = b`, where `L` and `U` are the block ILU factors.
526  void Mult(const Vector &b, Vector &x) const;
527 
528  /** Get the I array for the block CSR representation of the factorization.
529  * Similar to SparseMatrix::GetI(). Mostly used for testing.
530  */
531  int *GetBlockI() { return IB.GetData(); }
532 
533  /** Get the J array for the block CSR representation of the factorization.
534  * Similar to SparseMatrix::GetJ(). Mostly used for testing.
535  */
536  int *GetBlockJ() { return JB.GetData(); }
537 
538  /** Get the data array for the block CSR representation of the factorization.
539  * Similar to SparseMatrix::GetData(). Mostly used for testing.
540  */
541  double *GetBlockData() { return AB.Data(); }
542 
543 private:
544  /// Set up the block CSR structure corresponding to a sparse matrix @a A
545  void CreateBlockPattern(const class SparseMatrix &A);
546 
547  /// Perform the block ILU factorization
548  void Factorize();
549 
550  int block_size;
551 
552  /// Fill level for block ILU(k) factorizations. Only k=0 is supported.
553  int k_fill;
554 
555  Reordering reordering;
556 
557  /// Temporary vector used in the Mult() function.
558  mutable Vector y;
559 
560  /// Permutation and inverse permutation vectors for the block reordering.
561  Array<int> P, Pinv;
562 
563  /** Block CSR storage of the factorization. The block upper triangular part
564  * stores the U factor. The L factor implicitly has identity on the diagonal
565  * blocks, and the rest of L is given by the strictly block lower triangular
566  * part.
567  */
568  Array<int> IB, ID, JB;
569  DenseTensor AB;
570 
571  /// DB(i) stores the LU factorization of the i'th diagonal block
572  mutable DenseTensor DB;
573  /// Pivot arrays for the LU factorizations given by #DB
574  mutable Array<int> ipiv;
575 };
576 
577 #ifdef MFEM_USE_SUITESPARSE
578 
579 /// Direct sparse solver using UMFPACK
580 class UMFPackSolver : public Solver
581 {
582 protected:
585  void *Numeric;
586  SuiteSparse_long *AI, *AJ;
587 
588  void Init();
589 
590 public:
591  double Control[UMFPACK_CONTROL];
592  mutable double Info[UMFPACK_INFO];
593 
594  /** @brief For larger matrices, if the solver fails, set the parameter @a
595  _use_long_ints = true. */
596  UMFPackSolver(bool _use_long_ints = false)
597  : use_long_ints(_use_long_ints) { Init(); }
598  /** @brief Factorize the given SparseMatrix using the defaults. For larger
599  matrices, if the solver fails, set the parameter @a _use_long_ints =
600  true. */
601  UMFPackSolver(SparseMatrix &A, bool _use_long_ints = false)
602  : use_long_ints(_use_long_ints) { Init(); SetOperator(A); }
603 
604  /** @brief Factorize the given Operator @a op which must be a SparseMatrix.
605 
606  The factorization uses the parameters set in the #Control data member.
607  @note This method calls SparseMatrix::SortColumnIndices() with @a op,
608  modifying the matrix if the column indices are not already sorted. */
609  virtual void SetOperator(const Operator &op);
610 
611  /// Set the print level field in the #Control data member.
612  void SetPrintLevel(int print_lvl) { Control[UMFPACK_PRL] = print_lvl; }
613 
614  virtual void Mult(const Vector &b, Vector &x) const;
615  virtual void MultTranspose(const Vector &b, Vector &x) const;
616 
617  virtual ~UMFPackSolver();
618 };
619 
620 /// Direct sparse solver using KLU
621 class KLUSolver : public Solver
622 {
623 protected:
625  klu_symbolic *Symbolic;
626  klu_numeric *Numeric;
627 
628  void Init();
629 
630 public:
632  : mat(0),Symbolic(0),Numeric(0)
633  { Init(); }
635  : mat(0),Symbolic(0),Numeric(0)
636  { Init(); SetOperator(A); }
637 
638  // Works on sparse matrices only; calls SparseMatrix::SortColumnIndices().
639  virtual void SetOperator(const Operator &op);
640 
641  virtual void Mult(const Vector &b, Vector &x) const;
642  virtual void MultTranspose(const Vector &b, Vector &x) const;
643 
644  virtual ~KLUSolver();
645 
646  mutable klu_common Common;
647 };
648 
649 #endif // MFEM_USE_SUITESPARSE
650 
651 }
652 
653 #endif // MFEM_SOLVERS
const Vector * GetInequalityVec_Hi() const
Definition: solvers.hpp:386
Conjugate gradient method.
Definition: solvers.hpp:152
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:2414
double Info[UMFPACK_INFO]
Definition: solvers.hpp:592
void Mult(const Vector &b, Vector &x) const
Solve the system LUx = b, where L and U are the block ILU factors.
Definition: solvers.cpp:2218
SparseMatrix * mat
Definition: solvers.hpp:584
int GetNumIterations() const
Definition: solvers.hpp:68
const Operator * oper
Definition: solvers.hpp:43
SuiteSparse_long * AI
Definition: solvers.hpp:586
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:358
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:415
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:407
virtual void ProcessNewState(const Vector &x) const
This method can be overloaded in derived classes to perform computations that need knowledge of the n...
Definition: solvers.hpp:330
NewtonSolver(MPI_Comm _comm)
Definition: solvers.hpp:309
T * GetData()
Returns the data.
Definition: array.hpp:98
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: solvers.cpp:151
FGMRES method.
Definition: solvers.hpp:203
int * GetBlockI()
Definition: solvers.hpp:531
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:774
Direct sparse solver using KLU.
Definition: solvers.hpp:621
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
virtual void SetSolver(Solver &solver)
Set the linear solver for inverting the Jacobian.
Definition: solvers.hpp:315
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:2515
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:601
void SetOperator(const Operator &op)
Definition: solvers.cpp:1987
MINRES method.
Definition: solvers.hpp:262
const Operator * GetC() const
Definition: solvers.hpp:382
Reordering
The reordering method used by the BlockILU factorization.
Definition: solvers.hpp:497
virtual double CalcObjective(const Vector &x) const =0
Objective F(x). In parallel, the result should be reduced over tasks.
void add(const Vector &v1, const Vector &v2, Vector &v)
Definition: vector.cpp:238
virtual ~OptimizationSolver()
Definition: solvers.hpp:404
double Norm(const Vector &x) const
Definition: solvers.hpp:54
const Vector * d_hi
Definition: solvers.hpp:364
double GetFinalNorm() const
Definition: solvers.hpp:70
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:1107
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:180
void SetKDim(int dim)
Definition: solvers.hpp:215
const Vector * x_lo
Definition: solvers.hpp:364
FGMRESSolver(MPI_Comm _comm)
Definition: solvers.hpp:212
Direct sparse solver using UMFPACK.
Definition: solvers.hpp:580
const Vector * x_hi
Definition: solvers.hpp:364
UMFPackSolver(bool _use_long_ints=false)
For larger matrices, if the solver fails, set the parameter _use_long_ints = true.
Definition: solvers.hpp:596
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:594
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:1144
void SetLinearConstraint(const Vector &_w, double _a)
Definition: solvers.cpp:1632
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.hpp:275
void SetPrintLevel(int print_lvl)
Definition: solvers.cpp:70
const Vector * GetBoundsVec_Lo() const
Definition: solvers.hpp:387
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.cpp:97
double solve(double l, const Vector &xt, Vector &x, int &nclip) const
Solve QP at fixed lambda.
Definition: solvers.hpp:433
int GetNumConstraints() const
Definition: solvers.cpp:1605
Data type sparse matrix.
Definition: sparsemat.hpp:40
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:1302
virtual ~UMFPackSolver()
Definition: solvers.cpp:2448
Jacobi smoothing for a given bilinear form (no matrix necessary).
Definition: solvers.hpp:84
double b
Definition: lissajous.cpp:42
void SetMaxIter(int max_it)
Definition: solvers.hpp:65
void median(const Vector &lo, const Vector &hi)
v = median(v,lo,hi) entrywise. Implementation assumes lo &lt;= hi.
Definition: vector.cpp:448
SuiteSparse_long * AJ
Definition: solvers.hpp:586
void SetKDim(int dim)
Set the number of iteration to perform between restarts, default is 50.
Definition: solvers.hpp:197
Newton&#39;s method for solving F(x)=b for a given operator F.
Definition: solvers.hpp:300
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:516
void print_iteration(int it, double r, double l) const
Definition: solvers.cpp:1638
void Setup(const Vector &diag)
Definition: solvers.cpp:140
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:529
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: solvers.cpp:2470
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.cpp:1130
double Dot(const Vector &x, const Vector &y) const
Definition: solvers.cpp:54
int GetConverged() const
Definition: solvers.hpp:69
BiCGSTABSolver(MPI_Comm _comm)
Definition: solvers.hpp:242
klu_symbolic * Symbolic
Definition: solvers.hpp:625
const Vector * GetBoundsVec_Hi() const
Definition: solvers.hpp:388
SLBQPOptimizer(MPI_Comm _comm)
Definition: solvers.hpp:460
BlockILU(int block_size_, Reordering reordering_=Reordering::MINIMUM_DISCARDED_FILL, int k_fill_=0)
Definition: solvers.cpp:1969
virtual void Mult(const Vector &b, Vector &x) const
Solve the nonlinear system with right-hand side b.
Definition: solvers.cpp:1339
OperatorJacobiSmoother(const BilinearForm &a, const Array< int > &ess_tdof_list, const double damping=1.0)
Definition: solvers.cpp:109
CGSolver(MPI_Comm _comm)
Definition: solvers.hpp:163
void SetPrintLevel(int print_lvl)
Set the print level field in the Control data member.
Definition: solvers.hpp:612
Stationary linear iteration: x &lt;- x + B (b - A x)
Definition: solvers.hpp:120
void SetInequalityConstraint(const Vector &dl, const Vector &dh)
Definition: solvers.cpp:1587
void SetEqualityConstraint(const Vector &c)
Definition: solvers.cpp:1579
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.hpp:413
void SetAbsTol(double atol)
Definition: solvers.hpp:64
void SetRelTol(double rtol)
Definition: solvers.hpp:63
double Control[UMFPACK_CONTROL]
Definition: solvers.hpp:591
Abstract base class for iterative solver.
Definition: solvers.hpp:34
OptimizationSolver(MPI_Comm _comm)
Definition: solvers.hpp:402
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:978
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:2501
double * GetBlockData()
Definition: solvers.hpp:541
SLISolver(MPI_Comm _comm)
Definition: solvers.hpp:131
GMRES method.
Definition: solvers.hpp:184
GMRESSolver(MPI_Comm _comm)
Definition: solvers.hpp:193
const Vector * GetEqualityVec() const
Definition: solvers.hpp:384
virtual void SetOptimizationProblem(const OptimizationProblem &prob)
Definition: solvers.cpp:1613
void SetBounds(const Vector &_lo, const Vector &_hi)
Definition: solvers.cpp:1626
const Operator * GetD() const
Definition: solvers.hpp:383
double a
Definition: lissajous.cpp:41
void UpdateVectors()
Definition: solvers.cpp:174
Abstract solver for OptimizationProblems.
Definition: solvers.hpp:394
klu_common Common
Definition: solvers.hpp:646
virtual double ComputeScalingFactor(const Vector &x, const Vector &b) const
This method can be overloaded in derived classes to implement line search algorithms.
Definition: solvers.hpp:325
OptimizationProblem(int insize, const Operator *C_, const Operator *D_)
In parallel, insize is the number of the local true dofs.
Definition: solvers.cpp:1569
void SetSolutionBounds(const Vector &xl, const Vector &xh)
Definition: solvers.cpp:1597
const Vector * GetInequalityVec_Lo() const
Definition: solvers.hpp:385
virtual void Mult(const Vector &xt, Vector &x) const =0
Operator application: y=A(x).
void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: solvers.hpp:105
int dim
Definition: ex24.cpp:43
int * GetBlockJ()
Definition: solvers.hpp:536
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:245
virtual void Mult(const Vector &xt, Vector &x) const
Definition: solvers.cpp:1645
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:166
virtual void CalcObjectiveGrad(const Vector &x, Vector &grad) const
The result grad is expected to enter with the correct size.
Definition: solvers.hpp:375
Vector data type.
Definition: vector.hpp:48
MINRESSolver(MPI_Comm _comm)
Definition: solvers.hpp:272
const Vector * d_lo
Definition: solvers.hpp:364
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:91
KLUSolver(SparseMatrix &A)
Definition: solvers.hpp:634
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:1418
virtual ~KLUSolver()
Definition: solvers.cpp:2529
double * Data()
Definition: densemat.hpp:819
void UpdateVectors()
Definition: solvers.cpp:351
BiCGSTAB method.
Definition: solvers.hpp:231
SparseMatrix * mat
Definition: solvers.hpp:624
Base class for solvers.
Definition: operator.hpp:500
const OptimizationProblem * problem
Definition: solvers.hpp:397
Abstract operator.
Definition: operator.hpp:24
const Operator * D
Definition: solvers.hpp:363
Rank 3 tensor (array of matrices)
Definition: densemat.hpp:725
virtual void SetOptimizationProblem(const OptimizationProblem &prob)
Definition: solvers.hpp:408
const Vector * c_e
Definition: solvers.hpp:364
klu_numeric * Numeric
Definition: solvers.hpp:626
int width
Dimension of the input / number of columns in the matrix.
Definition: operator.hpp:28
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:2382
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:134
const Operator * C
Not owned, some can remain unused (NULL).
Definition: solvers.hpp:363
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:323
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:942
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.cpp:1328
virtual void SetOperator(const Operator &op)
Factorize the given Operator op which must be a SparseMatrix.
Definition: solvers.cpp:2285