MFEM  v4.5.2
Finite element discretization library
solvers.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2023, 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 #include "handle.hpp"
18 #include <memory>
19 
20 #ifdef MFEM_USE_MPI
21 #include <mpi.h>
22 #endif
23 
24 #ifdef MFEM_USE_SUITESPARSE
25 #include "sparsemat.hpp"
26 #include <umfpack.h>
27 #include <klu.h>
28 #endif
29 
30 namespace mfem
31 {
32 
33 class BilinearForm;
34 
35 /// Abstract base class for an iterative solver monitor
37 {
38 protected:
39  /// The last IterativeSolver to which this monitor was attached.
41 
42 public:
44 
46 
47  /// Monitor the residual vector r
48  virtual void MonitorResidual(int it, double norm, const Vector &r,
49  bool final)
50  {
51  }
52 
53  /// Monitor the solution vector x
54  virtual void MonitorSolution(int it, double norm, const Vector &x,
55  bool final)
56  {
57  }
58 
59  /** @brief This method is invoked by IterativeSolver::SetMonitor, informing
60  the monitor which IterativeSolver is using it. */
61  void SetIterativeSolver(const IterativeSolver &solver)
62  { iter_solver = &solver; }
63 };
64 
65 /// Abstract base class for iterative solver
66 class IterativeSolver : public Solver
67 {
68 public:
69  /** @brief Settings for the output behavior of the IterativeSolver.
70 
71  By default, all output is suppressed. The construction of the desired
72  print level can be achieved through a builder pattern, for example
73 
74  PrintLevel().Errors().Warnings()
75 
76  constructs the print level with only errors and warnings enabled.
77  */
78  struct PrintLevel
79  {
80  /** @brief If a fatal problem has been detected the failure will be
81  reported to @ref mfem::err. */
82  bool errors = false;
83  /** @brief If a non-fatal problem has been detected some context-specific
84  information will be reported to @ref mfem::out */
85  bool warnings = false;
86  /** @brief Detailed information about each iteration will be reported to
87  @ref mfem::out */
88  bool iterations = false;
89  /** @brief A summary of the solver process will be reported after the last
90  iteration to @ref mfem::out */
91  bool summary = false;
92  /** @brief Information about the first and last iteration will be printed
93  to @ref mfem::out */
94  bool first_and_last = false;
95 
96  /// Initializes the print level to suppress
97  PrintLevel() = default;
98 
99  /** @name Builder
100  These methods are utilized to construct PrintLevel objects through a
101  builder approach by chaining the function calls in this group. */
102  ///@{
103  PrintLevel &None() { *this = PrintLevel(); return *this; }
104  PrintLevel &Warnings() { warnings=true; return *this; }
105  PrintLevel &Errors() { errors=true; return *this; }
106  PrintLevel &Iterations() { iterations=true; return *this; }
107  PrintLevel &FirstAndLast() { first_and_last=true; return *this; }
108  PrintLevel &Summary() { summary=true; return *this; }
110  ///@}
111  };
112 
113 #ifdef MFEM_USE_MPI
114 private:
115  int dot_prod_type; // 0 - local, 1 - global over 'comm'
116  MPI_Comm comm = MPI_COMM_NULL;
117 #endif
118 
119 protected:
120  const Operator *oper;
123 
124  /// @name Reporting (protected attributes and member functions)
125  ///@{
126 
127  /** @brief (DEPRECATED) Legacy print level definition, which is left for
128  compatibility with custom iterative solvers.
129  @deprecated #print_options should be used instead. */
130  int print_level = -1;
131 
132  /** @brief Output behavior for the iterative solver.
133 
134  This primarily controls the output behavior of the iterative solvers
135  provided by this library. This member must be synchronized with
136  #print_level to ensure compatibility with custom iterative solvers. */
138 
139  /// Convert a legacy print level integer to a PrintLevel object
141 
142  /// @brief Use some heuristics to guess a legacy print level corresponding to
143  /// the given PrintLevel.
144  static int GuessLegacyPrintLevel(PrintLevel);
145  ///@}
146 
147  /// @name Convergence (protected attributes)
148  ///@{
149 
150  /// Limit for the number of iterations the solver is allowed to do
151  int max_iter;
152 
153  /// Relative tolerance.
154  double rel_tol;
155 
156  /// Absolute tolerance.
157  double abs_tol;
158 
159  ///@}
160 
161  /// @name Solver statistics (protected attributes)
162  /// Every IterativeSolver is expected to define these in its Mult() call.
163  ///@{
164 
165  mutable int final_iter = -1;
166  mutable bool converged = false;
167  mutable double initial_norm = -1.0, final_norm = -1.0;
168 
169  ///@}
170 
171  double Dot(const Vector &x, const Vector &y) const;
172  double Norm(const Vector &x) const { return sqrt(Dot(x, x)); }
173  void Monitor(int it, double norm, const Vector& r, const Vector& x,
174  bool final=false) const;
175 
176 public:
177  IterativeSolver();
178 
179 #ifdef MFEM_USE_MPI
180  IterativeSolver(MPI_Comm comm_);
181 #endif
182 
183  /** @name Convergence
184  @brief Termination criteria for the iterative solvers.
185 
186  @details While the convergence criterion is solver specific, most of the
187  provided iterative solvers use one of the following criteria
188 
189  \f$ ||r||_X \leq tol_{rel}||r_0||_X \f$,
190 
191  \f$ ||r||_X \leq tol_{abs} \f$,
192 
193  \f$ ||r||_X \leq \max\{ tol_{abs}, tol_{rel} ||r_0||_X \} \f$,
194 
195  where X denotes the space in which the norm is measured. The choice of
196  X depends on the specific iterative solver.
197  */
198  ///@{
199  void SetRelTol(double rtol) { rel_tol = rtol; }
200  void SetAbsTol(double atol) { abs_tol = atol; }
201  void SetMaxIter(int max_it) { max_iter = max_it; }
202  ///@}
203 
204  /** @name Reporting
205  These options control the internal reporting behavior into ::mfem::out
206  and ::mfem::err of the iterative solvers.
207  */
208  ///@{
209 
210  /// @brief Legacy method to set the level of verbosity of the solver output.
211  /** This is the old way to control what information will be printed to
212  ::mfem::out and ::mfem::err. The behavior for the print level for all
213  iterative solvers is:
214 
215  - -1: Suppress all outputs.
216  - 0: Print information about all detected issues (e.g. no convergence).
217  - 1: Same as level 0, but with detailed information about each
218  iteration.
219  - 2: Print detected issues and a summary when the solver terminates.
220  - 3: Same as 2, but print also the first and last iterations.
221  - >3: Custom print options which are dependent on the specific solver.
222 
223  In parallel, only rank 0 produces output.
224 
225  @note It is recommended to use @ref SetPrintLevel(PrintLevel) instead.
226 
227  @note Some derived classes, like KINSolver, redefine this method and use
228  their own set of print level constants. */
229  virtual void SetPrintLevel(int print_lvl);
230 
231  /// @brief Set the level of verbosity of the solver output.
232  /** In parallel, only rank 0 produces outputs. Errors are output to
233  ::mfem::err and all other information to ::mfem::out.
234 
235  @note Not all subclasses of IterativeSolver support all possible options.
236 
237  @note Some derived classes, like KINSolver, disable this method in favor
238  of SetPrintLevel(int).
239 
240  @sa PrintLevel for possible options.
241  */
242  virtual void SetPrintLevel(PrintLevel);
243  ///@}
244 
245  /// @name Solver statistics.
246  /// These are valid after the call to Mult().
247  ///@{
248 
249  /// Returns the number of iterations taken during the last call to Mult()
250  int GetNumIterations() const { return final_iter; }
251  /// Returns true if the last call to Mult() converged successfully.
252  bool GetConverged() const { return converged; }
253  /// @brief Returns the initial residual norm from the last call to Mult().
254  ///
255  /// This function returns the norm of the residual (or preconditioned
256  /// residual, depending on the solver), computed before the start of the
257  /// iteration.
258  double GetInitialNorm() const { return initial_norm; }
259  /// @brief Returns the final residual norm after termination of the solver
260  /// during the last call to Mult().
261  ///
262  /// This function returns the norm of the residual (or preconditioned
263  /// residual, depending on the solver), corresponding to the returned
264  /// solution.
265  double GetFinalNorm() const { return final_norm; }
266  /// @brief Returns the final residual norm after termination of the solver
267  /// during the last call to Mult(), divided by the initial residual norm.
268  /// Returns -1 if one of these norms is left undefined by the solver.
269  ///
270  /// @sa GetFinalNorm(), GetInitialNorm()
271  double GetFinalRelNorm() const
272  {
273  if (final_norm < 0.0 || initial_norm < 0.0) { return -1.0; }
274  return final_norm / initial_norm;
275  }
276 
277  ///@}
278 
279  /// This should be called before SetOperator
280  virtual void SetPreconditioner(Solver &pr);
281 
282  /// Also calls SetOperator for the preconditioner if there is one
283  virtual void SetOperator(const Operator &op) override;
284 
285  /// Set the iterative solver monitor
287  { monitor = &m; m.SetIterativeSolver(*this); }
288 
289 #ifdef MFEM_USE_MPI
290  /** @brief Return the associated MPI communicator, or MPI_COMM_NULL if no
291  communicator is set. */
292  MPI_Comm GetComm() const
293  { return dot_prod_type == 0 ? MPI_COMM_NULL : comm; }
294 #endif
295 };
296 
297 
298 /// Jacobi smoothing for a given bilinear form (no matrix necessary).
299 /** Useful with tensorized, partially assembled operators. Can also be defined
300  by given diagonal vector. This is basic Jacobi iteration; for tolerances,
301  iteration control, etc. wrap with SLISolver. */
303 {
304 public:
305  /** @brief Default constructor: the diagonal will be computed by subsequent
306  calls to SetOperator() using the Operator method AssembleDiagonal. */
307  /** In this case the array of essential tdofs will be empty. */
308  OperatorJacobiSmoother(const double damping=1.0);
309 
310  /** Setup a Jacobi smoother with the diagonal of @a a obtained by calling
311  a.AssembleDiagonal(). It is assumed that the underlying operator acts as
312  the identity on entries in ess_tdof_list, corresponding to (assembled)
313  DIAG_ONE policy or ConstrainedOperator in the matrix-free setting.
314 
315  @note For objects created with this constructor, calling SetOperator()
316  will only set the internal Operator pointer to the given new Operator
317  without any other changes to the object. This is done to preserve the
318  original behavior of this class. */
320  const Array<int> &ess_tdof_list,
321  const double damping=1.0);
322 
323  /** Application is by the *inverse* of the given vector. It is assumed that
324  the underlying operator acts as the identity on entries in ess_tdof_list,
325  corresponding to (assembled) DIAG_ONE policy or ConstrainedOperator in
326  the matrix-free setting.
327 
328  @note For objects created with this constructor, calling SetOperator()
329  will only set the internal Operator pointer to the given new Operator
330  without any other changes to the object. This is done to preserve the
331  original behavior of this class. */
333  const Array<int> &ess_tdof_list,
334  const double damping=1.0);
335 
337 
338  /// Replace diagonal entries with their absolute values.
339  void SetPositiveDiagonal(bool pos_diag = true) { use_abs_diag = pos_diag; }
340 
341  void Mult(const Vector &x, Vector &y) const;
342  void MultTranspose(const Vector &x, Vector &y) const { Mult(x, y); }
343 
344  /** @brief Recompute the diagonal using the method AssembleDiagonal of the
345  given new Operator, @a op. */
346  /** Note that (Par)BilinearForm operators are treated similar to the way they
347  are treated in the constructor that takes a BilinearForm parameter.
348  Specifically, this means that the OperatorJacobiSmoother will work with
349  true-dof vectors even though the size of the BilinearForm may be
350  different.
351 
352  When the new Operator, @a op, is not a (Par)BilinearForm, any previously
353  set array of essential true-dofs will be thrown away because in this case
354  any essential b.c. will be handled by the AssembleDiagonal method. */
355  void SetOperator(const Operator &op);
356 
357 private:
358  Vector dinv;
359  const double damping;
360  const Array<int> *ess_tdof_list; // not owned; may be NULL
361  mutable Vector residual;
362  /// Uses absolute values of the diagonal entries.
363  bool use_abs_diag = false;
364 
365  const Operator *oper; // not owned
366 
367  // To preserve the original behavior, some constructors set this flag to
368  // false to disallow updating the OperatorJacobiSmoother with SetOperator.
369  const bool allow_updates;
370 
371 public:
372  void Setup(const Vector &diag);
373 };
374 
375 /// Chebyshev accelerated smoothing with given vector, no matrix necessary
376 /** Potentially useful with tensorized operators, for example. This is just a
377  very basic Chebyshev iteration, if you want tolerances, iteration control,
378  etc. wrap this with SLISolver. */
380 {
381 public:
382  /** Application is by *inverse* of the given vector. It is assumed the
383  underlying operator acts as the identity on entries in ess_tdof_list,
384  corresponding to (assembled) DIAG_ONE policy or ConstrainedOperator in
385  the matrix-free setting. The estimated largest eigenvalue of the
386  diagonally preconditoned operator must be provided via
387  max_eig_estimate. */
388  OperatorChebyshevSmoother(const Operator &oper_, const Vector &d,
389  const Array<int>& ess_tdof_list,
390  int order, double max_eig_estimate);
391 
392  /// Deprecated: see pass-by-reference version above
393  MFEM_DEPRECATED
394  OperatorChebyshevSmoother(const Operator* oper_, const Vector &d,
395  const Array<int>& ess_tdof_list,
396  int order, double max_eig_estimate);
397 
398  /** Application is by *inverse* of the given vector. It is assumed the
399  underlying operator acts as the identity on entries in ess_tdof_list,
400  corresponding to (assembled) DIAG_ONE policy or ConstrainedOperator in
401  the matrix-free setting. The largest eigenvalue of the diagonally
402  preconditoned operator is estimated internally via a power method. The
403  accuracy of the estimated eigenvalue may be controlled via
404  power_iterations and power_tolerance. */
405 #ifdef MFEM_USE_MPI
406  OperatorChebyshevSmoother(const Operator &oper_, const Vector &d,
407  const Array<int>& ess_tdof_list,
408  int order, MPI_Comm comm = MPI_COMM_NULL,
409  int power_iterations = 10,
410  double power_tolerance = 1e-8);
411 
412  /// Deprecated: see pass-by-reference version above
413  MFEM_DEPRECATED
414  OperatorChebyshevSmoother(const Operator* oper_, const Vector &d,
415  const Array<int>& ess_tdof_list,
416  int order, MPI_Comm comm = MPI_COMM_NULL,
417  int power_iterations = 10,
418  double power_tolerance = 1e-8);
419 #else
420  OperatorChebyshevSmoother(const Operator &oper_, const Vector &d,
421  const Array<int>& ess_tdof_list,
422  int order, int power_iterations = 10,
423  double power_tolerance = 1e-8);
424 
425  /// Deprecated: see pass-by-reference version above
426  MFEM_DEPRECATED
427  OperatorChebyshevSmoother(const Operator* oper_, const Vector &d,
428  const Array<int>& ess_tdof_list,
429  int order, int power_iterations = 10,
430  double power_tolerance = 1e-8);
431 #endif
432 
434 
435  void Mult(const Vector&x, Vector &y) const;
436 
437  void MultTranspose(const Vector &x, Vector &y) const { Mult(x, y); }
438 
439  void SetOperator(const Operator &op_)
440  {
441  oper = &op_;
442  }
443 
444  void Setup();
445 
446 private:
447  const int order;
448  double max_eig_estimate;
449  const int N;
450  Vector dinv;
451  const Vector &diag;
452  Array<double> coeffs;
453  const Array<int>& ess_tdof_list;
454  mutable Vector residual;
455  mutable Vector helperVector;
456  const Operator* oper;
457 };
458 
459 
460 /// Stationary linear iteration: x <- x + B (b - A x)
462 {
463 protected:
464  mutable Vector r, z;
465 
466  void UpdateVectors();
467 
468 public:
469  SLISolver() { }
470 
471 #ifdef MFEM_USE_MPI
472  SLISolver(MPI_Comm comm_) : IterativeSolver(comm_) { }
473 #endif
474 
475  virtual void SetOperator(const Operator &op)
477 
478  virtual void Mult(const Vector &b, Vector &x) const;
479 };
480 
481 /// Stationary linear iteration. (tolerances are squared)
482 void SLI(const Operator &A, const Vector &b, Vector &x,
483  int print_iter = 0, int max_num_iter = 1000,
484  double RTOLERANCE = 1e-12, double ATOLERANCE = 1e-24);
485 
486 /// Preconditioned stationary linear iteration. (tolerances are squared)
487 void SLI(const Operator &A, Solver &B, const Vector &b, Vector &x,
488  int print_iter = 0, int max_num_iter = 1000,
489  double RTOLERANCE = 1e-12, double ATOLERANCE = 1e-24);
490 
491 
492 /// Conjugate gradient method
493 class CGSolver : public IterativeSolver
494 {
495 protected:
496  mutable Vector r, d, z;
497 
498  void UpdateVectors();
499 
500 public:
501  CGSolver() { }
502 
503 #ifdef MFEM_USE_MPI
504  CGSolver(MPI_Comm comm_) : IterativeSolver(comm_) { }
505 #endif
506 
507  virtual void SetOperator(const Operator &op)
509 
510  virtual void Mult(const Vector &b, Vector &x) const;
511 };
512 
513 /// Conjugate gradient method. (tolerances are squared)
514 void CG(const Operator &A, const Vector &b, Vector &x,
515  int print_iter = 0, int max_num_iter = 1000,
516  double RTOLERANCE = 1e-12, double ATOLERANCE = 1e-24);
517 
518 /// Preconditioned conjugate gradient method. (tolerances are squared)
519 void PCG(const Operator &A, Solver &B, const Vector &b, Vector &x,
520  int print_iter = 0, int max_num_iter = 1000,
521  double RTOLERANCE = 1e-12, double ATOLERANCE = 1e-24);
522 
523 
524 /// GMRES method
526 {
527 protected:
528  int m; // see SetKDim()
529 
530 public:
531  GMRESSolver() { m = 50; }
532 
533 #ifdef MFEM_USE_MPI
534  GMRESSolver(MPI_Comm comm_) : IterativeSolver(comm_) { m = 50; }
535 #endif
536 
537  /// Set the number of iteration to perform between restarts, default is 50.
538  void SetKDim(int dim) { m = dim; }
539 
540  virtual void Mult(const Vector &b, Vector &x) const;
541 };
542 
543 /// FGMRES method
545 {
546 protected:
547  int m;
548 
549 public:
550  FGMRESSolver() { m = 50; }
551 
552 #ifdef MFEM_USE_MPI
553  FGMRESSolver(MPI_Comm comm_) : IterativeSolver(comm_) { m = 50; }
554 #endif
555 
556  void SetKDim(int dim) { m = dim; }
557 
558  virtual void Mult(const Vector &b, Vector &x) const;
559 };
560 
561 /// GMRES method. (tolerances are squared)
562 int GMRES(const Operator &A, Vector &x, const Vector &b, Solver &M,
563  int &max_iter, int m, double &tol, double atol, int printit);
564 
565 /// GMRES method. (tolerances are squared)
566 void GMRES(const Operator &A, Solver &B, const Vector &b, Vector &x,
567  int print_iter = 0, int max_num_iter = 1000, int m = 50,
568  double rtol = 1e-12, double atol = 1e-24);
569 
570 
571 /// BiCGSTAB method
573 {
574 protected:
575  mutable Vector p, phat, s, shat, t, v, r, rtilde;
576 
577  void UpdateVectors();
578 
579 public:
581 
582 #ifdef MFEM_USE_MPI
583  BiCGSTABSolver(MPI_Comm comm_) : IterativeSolver(comm_) { }
584 #endif
585 
586  virtual void SetOperator(const Operator &op)
588 
589  virtual void Mult(const Vector &b, Vector &x) const;
590 };
591 
592 /// BiCGSTAB method. (tolerances are squared)
593 int BiCGSTAB(const Operator &A, Vector &x, const Vector &b, Solver &M,
594  int &max_iter, double &tol, double atol, int printit);
595 
596 /// BiCGSTAB method. (tolerances are squared)
597 void BiCGSTAB(const Operator &A, Solver &B, const Vector &b, Vector &x,
598  int print_iter = 0, int max_num_iter = 1000,
599  double rtol = 1e-12, double atol = 1e-24);
600 
601 
602 /// MINRES method
604 {
605 protected:
606  mutable Vector v0, v1, w0, w1, q;
607  mutable Vector u1; // used in the preconditioned version
608 
609 public:
611 
612 #ifdef MFEM_USE_MPI
613  MINRESSolver(MPI_Comm comm_) : IterativeSolver(comm_) { }
614 #endif
615 
616  virtual void SetPreconditioner(Solver &pr)
617  {
619  if (oper) { u1.SetSize(width); }
620  }
621 
622  virtual void SetOperator(const Operator &op);
623 
624  virtual void Mult(const Vector &b, Vector &x) const;
625 };
626 
627 /// MINRES method without preconditioner. (tolerances are squared)
628 void MINRES(const Operator &A, const Vector &b, Vector &x, int print_it = 0,
629  int max_it = 1000, double rtol = 1e-12, double atol = 1e-24);
630 
631 /// MINRES method with preconditioner. (tolerances are squared)
632 void MINRES(const Operator &A, Solver &B, const Vector &b, Vector &x,
633  int print_it = 0, int max_it = 1000,
634  double rtol = 1e-12, double atol = 1e-24);
635 
636 
637 /// Newton's method for solving F(x)=b for a given operator F.
638 /** The method GetGradient() must be implemented for the operator F.
639  The preconditioner is used (in non-iterative mode) to evaluate
640  the action of the inverse gradient of the operator. */
642 {
643 protected:
644  mutable Vector r, c;
645  mutable Operator *grad;
646 
647  // Adaptive linear solver rtol variables
648 
649  // Method to determine rtol, 0 means the adaptive algorithm is deactivated.
650  int lin_rtol_type = 0;
651  // rtol to use in first iteration
652  double lin_rtol0;
653  // Maximum rtol
654  double lin_rtol_max;
655  // Function norm ||F(x)|| of the previous iterate
656  mutable double fnorm_last = 0.0;
657  // Linear residual norm of the previous iterate
658  mutable double lnorm_last = 0.0;
659  // Forcing term (linear residual rtol) from the previous iterate
660  mutable double eta_last = 0.0;
661  // Eisenstat-Walker factor gamma
662  double gamma;
663  // Eisenstat-Walker factor alpha
664  double alpha;
665 
666  /** @brief Method for the adaptive linear solver rtol invoked before the
667  linear solve. */
668  void AdaptiveLinRtolPreSolve(const Vector &x,
669  const int it,
670  const double fnorm) const;
671 
672  /** @brief Method for the adaptive linear solver rtol invoked after the
673  linear solve. */
674  void AdaptiveLinRtolPostSolve(const Vector &x,
675  const Vector &b,
676  const int it,
677  const double fnorm) const;
678 
679 public:
681 
682 #ifdef MFEM_USE_MPI
683  NewtonSolver(MPI_Comm comm_) : IterativeSolver(comm_) { }
684 #endif
685  virtual void SetOperator(const Operator &op);
686 
687  /// Set the linear solver for inverting the Jacobian.
688  /** This method is equivalent to calling SetPreconditioner(). */
689  virtual void SetSolver(Solver &solver) { prec = &solver; }
690 
691  /// Solve the nonlinear system with right-hand side @a b.
692  /** If `b.Size() != Height()`, then @a b is assumed to be zero. */
693  virtual void Mult(const Vector &b, Vector &x) const;
694 
695  /** @brief This method can be overloaded in derived classes to implement line
696  search algorithms. */
697  /** The base class implementation (NewtonSolver) simply returns 1. A return
698  value of 0 indicates a failure, interrupting the Newton iteration. */
699  virtual double ComputeScalingFactor(const Vector &x, const Vector &b) const
700  { return 1.0; }
701 
702  /** @brief This method can be overloaded in derived classes to perform
703  computations that need knowledge of the newest Newton state. */
704  virtual void ProcessNewState(const Vector &x) const { }
705 
706  /// Enable adaptive linear solver relative tolerance algorithm.
707  /** Compute a relative tolerance for the Krylov method after each nonlinear
708  iteration, based on the algorithm presented in [1].
709 
710  The maximum linear solver relative tolerance @a rtol_max should be < 1. For
711  @a type 1 the parameters @a alpha and @a gamma are ignored. For @a type 2
712  @a alpha has to be between 0 and 1 and @a gamma between 1 and 2.
713 
714  [1] Eisenstat, Stanley C., and Homer F. Walker. "Choosing the forcing terms
715  in an inexact Newton method."
716  */
717  void SetAdaptiveLinRtol(const int type = 2,
718  const double rtol0 = 0.5,
719  const double rtol_max = 0.9,
720  const double alpha = 0.5 * (1.0 + sqrt(5.0)),
721  const double gamma = 1.0);
722 };
723 
724 /** L-BFGS method for solving F(x)=b for a given operator F, by minimizing
725  the norm of F(x) - b. Requires only the action of the operator F. */
726 class LBFGSSolver : public NewtonSolver
727 {
728 protected:
729  int m = 10;
731 
733  {
734  for (int i = 0; i < skArray.Size(); i++)
735  {
736  delete skArray[i];
737  delete ykArray[i];
738  }
739  }
740 
742  {
744  skArray.SetSize(m);
745  ykArray.SetSize(m);
746  for (int i = 0; i < m; i++)
747  {
748  skArray[i] = new Vector(width);
749  ykArray[i] = new Vector(width);
750  skArray[i]->UseDevice(true);
751  ykArray[i]->UseDevice(true);
752  }
753  }
754 
755 public:
757 
758 #ifdef MFEM_USE_MPI
759  LBFGSSolver(MPI_Comm comm_) : NewtonSolver(comm_) { }
760 #endif
761 
762  virtual void SetOperator(const Operator &op)
763  {
766  }
767 
768  void SetHistorySize(int dim)
769  {
770  m = dim;
772  }
773 
774  /// Solve the nonlinear system with right-hand side @a b.
775  /** If `b.Size() != Height()`, then @a b is assumed to be zero. */
776  virtual void Mult(const Vector &b, Vector &x) const;
777 
778  virtual void SetPreconditioner(Solver &pr)
779  { MFEM_WARNING("L-BFGS won't use the given preconditioner."); }
780  virtual void SetSolver(Solver &solver)
781  { MFEM_WARNING("L-BFGS won't use the given solver."); }
782 
784 };
785 
786 
787 /** Adaptive restarted GMRES.
788  m_max and m_min(=1) are the maximal and minimal restart parameters.
789  m_step(=1) is the step to use for going from m_max and m_min.
790  cf(=0.4) is a desired convergence factor. */
791 int aGMRES(const Operator &A, Vector &x, const Vector &b,
792  const Operator &M, int &max_iter,
793  int m_max, int m_min, int m_step, double cf,
794  double &tol, double &atol, int printit);
795 
796 
797 /** Defines operators and constraints for the following optimization problem:
798  *
799  * Find x that minimizes the objective function F(x), subject to
800  * C(x) = c_e,
801  * d_lo <= D(x) <= d_hi,
802  * x_lo <= x <= x_hi.
803  *
804  * The operators F, C, D must take input of the same size (same width).
805  * Gradients of F, C, D might be needed, depending on the OptimizationSolver.
806  * When used with Hiop, gradients of C and D must be DenseMatrices.
807  * F always returns a scalar value, see CalcObjective(), CalcObjectiveGrad().
808  * C and D can have arbitrary heights.
809  * C and D can be NULL, meaning that their constraints are not used.
810  *
811  * When used in parallel, all Vectors are assumed to be true dof vectors, and
812  * the operators are expected to be defined for tdof vectors. */
814 {
815 protected:
816  /// Not owned, some can remain unused (NULL).
817  const Operator *C, *D;
818  const Vector *c_e, *d_lo, *d_hi, *x_lo, *x_hi;
819 
820 public:
821  const int input_size;
822 
823  /// In parallel, insize is the number of the local true dofs.
824  OptimizationProblem(int insize, const Operator *C_, const Operator *D_);
825 
826  /// Objective F(x). In parallel, the result should be reduced over tasks.
827  virtual double CalcObjective(const Vector &x) const = 0;
828  /// The result grad is expected to enter with the correct size.
829  virtual void CalcObjectiveGrad(const Vector &x, Vector &grad) const
830  { MFEM_ABORT("The objective gradient is not implemented."); }
831 
832  void SetEqualityConstraint(const Vector &c);
833  void SetInequalityConstraint(const Vector &dl, const Vector &dh);
834  void SetSolutionBounds(const Vector &xl, const Vector &xh);
835 
836  const Operator *GetC() const { return C; }
837  const Operator *GetD() const { return D; }
838  const Vector *GetEqualityVec() const { return c_e; }
839  const Vector *GetInequalityVec_Lo() const { return d_lo; }
840  const Vector *GetInequalityVec_Hi() const { return d_hi; }
841  const Vector *GetBoundsVec_Lo() const { return x_lo; }
842  const Vector *GetBoundsVec_Hi() const { return x_hi; }
843 
844  int GetNumConstraints() const;
845 };
846 
847 /// Abstract solver for OptimizationProblems.
849 {
850 protected:
852 
853 public:
855 #ifdef MFEM_USE_MPI
856  OptimizationSolver(MPI_Comm comm_): IterativeSolver(comm_), problem(NULL) { }
857 #endif
858  virtual ~OptimizationSolver() { }
859 
860  /** This function is virtual as solvers might need to perform some initial
861  * actions (e.g. validation) with the OptimizationProblem. */
863  { problem = &prob; }
864 
865  virtual void Mult(const Vector &xt, Vector &x) const = 0;
866 
867  virtual void SetPreconditioner(Solver &pr)
868  { MFEM_ABORT("Not meaningful for this solver."); }
869  virtual void SetOperator(const Operator &op)
870  { MFEM_ABORT("Not meaningful for this solver."); }
871 };
872 
873 /** SLBQP optimizer:
874  * (S)ingle (L)inearly Constrained with (B)ounds (Q)uadratic (P)rogram
875  *
876  * Minimize || x-x_t ||, subject to
877  * sum w_i x_i = a,
878  * x_lo <= x <= x_hi.
879  */
881 {
882 protected:
884  double a;
885 
886  /// Solve QP at fixed lambda
887  inline double solve(double l, const Vector &xt, Vector &x, int &nclip) const
888  {
889  add(xt, l, w, x);
890  if (problem == NULL) { x.median(lo,hi); }
891  else
892  {
895  }
896  nclip++;
897  if (problem == NULL) { return Dot(w, x) - a; }
898  else
899  {
900  Vector c(1);
901  // Includes parallel communication.
902  problem->GetC()->Mult(x, c);
903 
904  return c(0) - (*problem->GetEqualityVec())(0);
905  }
906  }
907 
908  inline void print_iteration(int it, double r, double l) const;
909 
910 public:
912 
913 #ifdef MFEM_USE_MPI
914  SLBQPOptimizer(MPI_Comm comm_) : OptimizationSolver(comm_) { }
915 #endif
916 
917  /** Setting an OptimizationProblem will overwrite the Vectors given by
918  * SetBounds and SetLinearConstraint. The objective function remains
919  * unchanged. */
920  virtual void SetOptimizationProblem(const OptimizationProblem &prob);
921 
922  void SetBounds(const Vector &lo_, const Vector &hi_);
923  void SetLinearConstraint(const Vector &w_, double a_);
924 
925  /** We let the target values play the role of the initial vector xt, from
926  * which the operator generates the optimal vector x. */
927  virtual void Mult(const Vector &xt, Vector &x) const;
928 };
929 
930 /** Block ILU solver:
931  * Performs a block ILU(k) approximate factorization with specified block
932  * size. Currently only k=0 is supported. This is useful as a preconditioner
933  * for DG-type discretizations, where the system matrix has a natural
934  * (elemental) block structure.
935  *
936  * In the case of DG discretizations, the block size should usually be set to
937  * either ndofs_per_element or vdim*ndofs_per_element (if the finite element
938  * space has Ordering::byVDIM). The block size must evenly divide the size of
939  * the matrix.
940  *
941  * Renumbering the blocks is also supported by specifying a reordering method.
942  * Currently greedy minimum discarded fill ordering and no reordering are
943  * supported. Renumbering the blocks can lead to a much better approximate
944  * factorization.
945  */
946 class BlockILU : public Solver
947 {
948 public:
949 
950  /// The reordering method used by the BlockILU factorization.
951  enum class Reordering
952  {
954  NONE
955  };
956 
957  /** Create an "empty" BlockILU solver. SetOperator must be called later to
958  * actually form the factorization
959  */
960  BlockILU(int block_size_,
962  int k_fill_ = 0);
963 
964  /** Create a block ILU approximate factorization for the matrix @a op.
965  * @a op should be of type either SparseMatrix or HypreParMatrix. In the
966  * case that @a op is a HypreParMatrix, the ILU factorization is performed
967  * on the diagonal blocks of the parallel decomposition.
968  */
969  BlockILU(const Operator &op, int block_size_ = 1,
971  int k_fill_ = 0);
972 
973  /** Perform the block ILU factorization for the matrix @a op.
974  * As in the constructor, @a op must either be a SparseMatrix or
975  * HypreParMatrix
976  */
977  void SetOperator(const Operator &op);
978 
979  /// Solve the system `LUx = b`, where `L` and `U` are the block ILU factors.
980  void Mult(const Vector &b, Vector &x) const;
981 
982  /** Get the I array for the block CSR representation of the factorization.
983  * Similar to SparseMatrix::GetI(). Mostly used for testing.
984  */
985  int *GetBlockI() { return IB.GetData(); }
986 
987  /** Get the J array for the block CSR representation of the factorization.
988  * Similar to SparseMatrix::GetJ(). Mostly used for testing.
989  */
990  int *GetBlockJ() { return JB.GetData(); }
991 
992  /** Get the data array for the block CSR representation of the factorization.
993  * Similar to SparseMatrix::GetData(). Mostly used for testing.
994  */
995  double *GetBlockData() { return AB.Data(); }
996 
997 private:
998  /// Set up the block CSR structure corresponding to a sparse matrix @a A
999  void CreateBlockPattern(const class SparseMatrix &A);
1000 
1001  /// Perform the block ILU factorization
1002  void Factorize();
1003 
1004  int block_size;
1005 
1006  /// Fill level for block ILU(k) factorizations. Only k=0 is supported.
1007  int k_fill;
1008 
1009  Reordering reordering;
1010 
1011  /// Temporary vector used in the Mult() function.
1012  mutable Vector y;
1013 
1014  /// Permutation and inverse permutation vectors for the block reordering.
1015  Array<int> P, Pinv;
1016 
1017  /** Block CSR storage of the factorization. The block upper triangular part
1018  * stores the U factor. The L factor implicitly has identity on the diagonal
1019  * blocks, and the rest of L is given by the strictly block lower triangular
1020  * part.
1021  */
1022  Array<int> IB, ID, JB;
1023  DenseTensor AB;
1024 
1025  /// DB(i) stores the LU factorization of the i'th diagonal block
1026  mutable DenseTensor DB;
1027  /// Pivot arrays for the LU factorizations given by #DB
1028  mutable Array<int> ipiv;
1029 };
1030 
1031 
1032 /// Monitor that checks whether the residual is zero at a given set of dofs.
1033 /** This monitor is useful for checking if the initial guess, rhs, operator, and
1034  preconditioner are properly setup for solving in the subspace with imposed
1035  essential boundary conditions. */
1037 {
1038 protected:
1039  const Array<int> *ess_dofs_list; ///< Not owned
1040 
1041 public:
1042  ResidualBCMonitor(const Array<int> &ess_dofs_list_)
1043  : ess_dofs_list(&ess_dofs_list_) { }
1044 
1045  void MonitorResidual(int it, double norm, const Vector &r,
1046  bool final) override;
1047 };
1048 
1049 
1050 #ifdef MFEM_USE_SUITESPARSE
1051 
1052 /// Direct sparse solver using UMFPACK
1053 class UMFPackSolver : public Solver
1054 {
1055 protected:
1058  void *Numeric;
1059  SuiteSparse_long *AI, *AJ;
1060 
1061  void Init();
1062 
1063 public:
1064  double Control[UMFPACK_CONTROL];
1065  mutable double Info[UMFPACK_INFO];
1066 
1067  /** @brief For larger matrices, if the solver fails, set the parameter @a
1068  use_long_ints_ = true. */
1069  UMFPackSolver(bool use_long_ints_ = false)
1070  : use_long_ints(use_long_ints_) { Init(); }
1071  /** @brief Factorize the given SparseMatrix using the defaults. For larger
1072  matrices, if the solver fails, set the parameter @a use_long_ints_ =
1073  true. */
1074  UMFPackSolver(SparseMatrix &A, bool use_long_ints_ = false)
1075  : use_long_ints(use_long_ints_) { Init(); SetOperator(A); }
1076 
1077  /** @brief Factorize the given Operator @a op which must be a SparseMatrix.
1078 
1079  The factorization uses the parameters set in the #Control data member.
1080  @note This method calls SparseMatrix::SortColumnIndices() with @a op,
1081  modifying the matrix if the column indices are not already sorted. */
1082  virtual void SetOperator(const Operator &op);
1083 
1084  /// Set the print level field in the #Control data member.
1085  void SetPrintLevel(int print_lvl) { Control[UMFPACK_PRL] = print_lvl; }
1086 
1087  virtual void Mult(const Vector &b, Vector &x) const;
1088  virtual void MultTranspose(const Vector &b, Vector &x) const;
1089 
1090  virtual ~UMFPackSolver();
1091 };
1092 
1093 /// Direct sparse solver using KLU
1094 class KLUSolver : public Solver
1095 {
1096 protected:
1098  klu_symbolic *Symbolic;
1099  klu_numeric *Numeric;
1100 
1101  void Init();
1102 
1103 public:
1105  : mat(0),Symbolic(0),Numeric(0)
1106  { Init(); }
1108  : mat(0),Symbolic(0),Numeric(0)
1109  { Init(); SetOperator(A); }
1110 
1111  // Works on sparse matrices only; calls SparseMatrix::SortColumnIndices().
1112  virtual void SetOperator(const Operator &op);
1113 
1114  virtual void Mult(const Vector &b, Vector &x) const;
1115  virtual void MultTranspose(const Vector &b, Vector &x) const;
1116 
1117  virtual ~KLUSolver();
1118 
1119  mutable klu_common Common;
1120 };
1121 
1122 #endif // MFEM_USE_SUITESPARSE
1123 
1124 /// Block diagonal solver for A, each block is inverted by direct solver
1126 {
1127  SparseMatrix& block_dof;
1128  mutable Array<int> local_dofs;
1129  mutable Vector sub_rhs;
1130  mutable Vector sub_sol;
1131  std::unique_ptr<DenseMatrixInverse[]> block_solvers;
1132 public:
1133  /// block_dof is a boolean matrix, block_dof(i, j) = 1 if j-th dof belongs to
1134  /// i-th block, block_dof(i, j) = 0 otherwise.
1135  DirectSubBlockSolver(const SparseMatrix& A, const SparseMatrix& block_dof);
1136  virtual void Mult(const Vector &x, Vector &y) const;
1137  virtual void SetOperator(const Operator &op) { }
1138 };
1139 
1140 /// Solver S such that I - A * S = (I - A * S1) * (I - A * S0).
1141 /// That is, S = S0 + S1 - S1 * A * S0.
1142 class ProductSolver : public Solver
1143 {
1144  OperatorPtr A;
1145  OperatorPtr S0;
1146  OperatorPtr S1;
1147 public:
1149  bool ownA, bool ownS0, bool ownS1)
1150  : Solver(A_->NumRows()), A(A_, ownA), S0(S0_, ownS0), S1(S1_, ownS1) { }
1151  virtual void Mult(const Vector &x, Vector &y) const;
1152  virtual void MultTranspose(const Vector &x, Vector &y) const;
1153  virtual void SetOperator(const Operator &op) { }
1154 };
1155 
1156 /// Solver wrapper which orthogonalizes the input and output vector
1157 /**
1158  * OrthoSolver wraps an existing Solver and orthogonalizes the input vector
1159  * before passing it to the Mult() method of the Solver. This is a convenience
1160  * implementation to handle e.g. a Poisson problem with pure Neumann boundary
1161  * conditions, where this procedure removes the Nullspace.
1162  */
1163 class OrthoSolver : public Solver
1164 {
1165 private:
1166 #ifdef MFEM_USE_MPI
1167  MPI_Comm mycomm;
1168  mutable HYPRE_BigInt global_size;
1169  const bool parallel;
1170 #else
1171  mutable int global_size;
1172 #endif
1173 
1174 public:
1175  OrthoSolver();
1176 #ifdef MFEM_USE_MPI
1177  OrthoSolver(MPI_Comm mycomm_);
1178 #endif
1179 
1180  /// Set the solver used by the OrthoSolver.
1181  /** The action of the OrthoSolver is given by P * s * P where P is the
1182  projection to the subspace of vectors with zero sum. Calling this method
1183  is required before calling SetOperator() or Mult(). */
1184  void SetSolver(Solver &s);
1185 
1186  /// Set the Operator that is the OrthoSolver is to invert (approximately).
1187  /** The Operator @a op is simply forwarded to the solver object given by
1188  SetSolver() which needs to be called before this method. Calling this
1189  method is optional when the solver already has an associated Operator. */
1190  virtual void SetOperator(const Operator &op);
1191 
1192  /** @brief Perform the action of the OrthoSolver: P * solver * P where P is
1193  the projection to the subspace of vectors with zero sum. */
1194  /** @note The projection P can be written as P = I - 1 1^T / (1^T 1) where
1195  I is the identity matrix and 1 is the column-vector with all components
1196  equal to 1. */
1197  void Mult(const Vector &b, Vector &x) const;
1198 
1199 private:
1200  Solver *solver = nullptr;
1201 
1202  mutable Vector b_ortho;
1203 
1204  void Orthogonalize(const Vector &v, Vector &v_ortho) const;
1205 };
1206 
1207 #ifdef MFEM_USE_MPI
1208 /** This smoother does relaxations on an auxiliary space (determined by a map
1209  from the original space to the auxiliary space provided by the user).
1210  The smoother on the auxiliary space is a HypreSmoother. Its options can be
1211  modified through GetSmoother.
1212  For example, the space can be the nullspace of div/curl, in which case the
1213  smoother can be used to construct a Hiptmair smoother. */
1214 class AuxSpaceSmoother : public Solver
1215 {
1216  OperatorPtr aux_map_;
1217  OperatorPtr aux_system_;
1218  OperatorPtr aux_smoother_;
1219  void Mult(const Vector &x, Vector &y, bool transpose) const;
1220 public:
1221  AuxSpaceSmoother(const HypreParMatrix &op, HypreParMatrix *aux_map,
1222  bool op_is_symmetric = true, bool own_aux_map = false);
1223  virtual void Mult(const Vector &x, Vector &y) const { Mult(x, y, false); }
1224  virtual void MultTranspose(const Vector &x, Vector &y) const { Mult(x, y, true); }
1225  virtual void SetOperator(const Operator &op) { }
1226  HypreSmoother& GetSmoother() { return *aux_smoother_.As<HypreSmoother>(); }
1227  using Operator::Mult;
1228 };
1229 #endif // MFEM_USE_MPI
1230 
1231 }
1232 
1233 #endif // MFEM_SOLVERS
const Array< int > * ess_dofs_list
Not owned.
Definition: solvers.hpp:1039
MINRESSolver(MPI_Comm comm_)
Definition: solvers.hpp:613
GMRESSolver(MPI_Comm comm_)
Definition: solvers.hpp:534
Conjugate gradient method.
Definition: solvers.hpp:493
CGSolver(MPI_Comm comm_)
Definition: solvers.hpp:504
HypreSmoother & GetSmoother()
Definition: solvers.hpp:1226
void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
Definition: solvers.hpp:437
double Info[UMFPACK_INFO]
Definition: solvers.hpp:1065
Chebyshev accelerated smoothing with given vector, no matrix necessary.
Definition: solvers.hpp:379
SparseMatrix * mat
Definition: solvers.hpp:1057
const Operator * oper
Definition: solvers.hpp:120
Array< Vector * > ykArray
Definition: solvers.hpp:730
void SetPositiveDiagonal(bool pos_diag=true)
Replace diagonal entries with their absolute values.
Definition: solvers.hpp:339
const Operator * GetC() const
Definition: solvers.hpp:836
SuiteSparse_long * AI
Definition: solvers.hpp:1059
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:762
double GetInitialNorm() const
Returns the initial residual norm from the last call to Mult().
Definition: solvers.hpp:258
void SetHistorySize(int dim)
Definition: solvers.hpp:768
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:869
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:512
bool summary
A summary of the solver process will be reported after the last iteration to mfem::out.
Definition: solvers.hpp:91
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:475
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: solvers.cpp:300
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:704
virtual void Mult(const Vector &b, Vector &x) const
Solve the nonlinear system with right-hand side b.
Definition: solvers.cpp:2005
bool first_and_last
Information about the first and last iteration will be printed to mfem::out.
Definition: solvers.hpp:94
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:1162
const Vector * GetBoundsVec_Hi() const
Definition: solvers.hpp:842
OperatorJacobiSmoother(const double damping=1.0)
Default constructor: the diagonal will be computed by subsequent calls to SetOperator() using the Ope...
Definition: solvers.cpp:200
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:1612
T * GetData()
Returns the data.
Definition: array.hpp:115
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:712
FGMRES method.
Definition: solvers.hpp:544
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:974
int * GetBlockI()
Definition: solvers.hpp:985
void InitializeStorageVectors()
Definition: solvers.hpp:741
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:525
Direct sparse solver using KLU.
Definition: solvers.hpp:1094
bool errors
If a fatal problem has been detected the failure will be reported to mfem::err.
Definition: solvers.hpp:82
int NumRows() const
Get the number of rows (size of output) of the Operator. Synonym with Height().
Definition: operator.hpp:69
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:689
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:3317
void SetOperator(const Operator &op)
Definition: solvers.cpp:2763
MINRES method.
Definition: solvers.hpp:603
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:1074
Reordering
The reordering method used by the BlockILU factorization.
Definition: solvers.hpp:951
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:2994
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: solvers.hpp:1223
virtual ~IterativeSolverMonitor()
Definition: solvers.hpp:45
ResidualBCMonitor(const Array< int > &ess_dofs_list_)
Definition: solvers.hpp:1042
int print_level
(DEPRECATED) Legacy print level definition, which is left for compatibility with custom iterative sol...
Definition: solvers.hpp:130
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: solvers.hpp:1153
PrintLevel()=default
Initializes the print level to suppress.
virtual double CalcObjective(const Vector &x) const =0
Objective F(x). In parallel, the result should be reduced over tasks.
SLISolver(MPI_Comm comm_)
Definition: solvers.hpp:472
void add(const Vector &v1, const Vector &v2, Vector &v)
Definition: vector.cpp:314
virtual ~OptimizationSolver()
Definition: solvers.hpp:858
void Monitor(int it, double norm, const Vector &r, const Vector &x, bool final=false) const
Definition: solvers.cpp:190
const Vector * d_hi
Definition: solvers.hpp:818
double GetFinalNorm() const
Returns the final residual norm after termination of the solver during the last call to Mult()...
Definition: solvers.hpp:265
virtual void SetOperator(const Operator &op)
Set the Operator that is the OrthoSolver is to invert (approximately).
Definition: solvers.cpp:3440
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: solvers.cpp:3370
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:1568
void SetKDim(int dim)
Definition: solvers.hpp:556
const Vector * GetInequalityVec_Lo() const
Definition: solvers.hpp:839
virtual void SetOperator(const Operator &op) override
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.cpp:179
const Vector * x_lo
Definition: solvers.hpp:818
double Norm(const Vector &x) const
Definition: solvers.hpp:172
Direct sparse solver using UMFPACK.
Definition: solvers.hpp:1053
Operator * grad
Definition: solvers.hpp:645
const Vector * x_hi
Definition: solvers.hpp:818
const class IterativeSolver * iter_solver
The last IterativeSolver to which this monitor was attached.
Definition: solvers.hpp:40
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.hpp:616
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
Definition: solvers.cpp:71
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:3331
int GetNumConstraints() const
Definition: solvers.cpp:2338
Data type sparse matrix.
Definition: sparsemat.hpp:50
double solve(double l, const Vector &xt, Vector &x, int &nclip) const
Solve QP at fixed lambda.
Definition: solvers.hpp:887
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:1778
virtual ~UMFPackSolver()
Definition: solvers.cpp:3264
Jacobi smoothing for a given bilinear form (no matrix necessary).
Definition: solvers.hpp:302
int GetNumIterations() const
Returns the number of iterations taken during the last call to Mult()
Definition: solvers.hpp:250
double b
Definition: lissajous.cpp:42
virtual void CalcObjectiveGrad(const Vector &x, Vector &grad) const
The result grad is expected to enter with the correct size.
Definition: solvers.hpp:829
void SetMaxIter(int max_it)
Definition: solvers.hpp:201
OperatorChebyshevSmoother(const Operator &oper_, const Vector &d, const Array< int > &ess_tdof_list, int order, double max_eig_estimate)
Definition: solvers.cpp:328
void SetLinearConstraint(const Vector &w_, double a_)
Definition: solvers.cpp:2365
void median(const Vector &lo, const Vector &hi)
v = median(v,lo,hi) entrywise. Implementation assumes lo <= hi.
Definition: vector.cpp:524
DirectSubBlockSolver(const SparseMatrix &A, const SparseMatrix &block_dof)
Definition: solvers.cpp:3355
SuiteSparse_long * AJ
Definition: solvers.hpp:1059
void SetKDim(int dim)
Set the number of iteration to perform between restarts, default is 50.
Definition: solvers.hpp:538
Newton&#39;s method for solving F(x)=b for a given operator F.
Definition: solvers.hpp:641
void print_iteration(int it, double r, double l) const
Definition: solvers.cpp:2371
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:892
void Setup(const Vector &diag)
Definition: solvers.cpp:277
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
Definition: solvers.hpp:1224
Parallel smoothers in hypre.
Definition: hypre.hpp:971
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:907
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: solvers.cpp:3286
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.cpp:1591
const Vector * GetInequalityVec_Hi() const
Definition: solvers.hpp:840
Array< Vector * > skArray
Definition: solvers.hpp:730
prob_type prob
Definition: ex25.cpp:153
klu_symbolic * Symbolic
Definition: solvers.hpp:1098
void SetIterativeSolver(const IterativeSolver &solver)
This method is invoked by IterativeSolver::SetMonitor, informing the monitor which IterativeSolver is...
Definition: solvers.hpp:61
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: solvers.hpp:1225
BlockILU(int block_size_, Reordering reordering_=Reordering::MINIMUM_DISCARDED_FILL, int k_fill_=0)
Definition: solvers.cpp:2745
void SetPrintLevel(int print_lvl)
Set the print level field in the Control data member.
Definition: solvers.hpp:1085
void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
Definition: solvers.hpp:342
Stationary linear iteration: x <- x + B (b - A x)
Definition: solvers.hpp:461
void SetInequalityConstraint(const Vector &dl, const Vector &dh)
Definition: solvers.cpp:2320
void SetEqualityConstraint(const Vector &c)
Definition: solvers.cpp:2312
Abstract base class for an iterative solver monitor.
Definition: solvers.hpp:36
const Vector * GetBoundsVec_Lo() const
Definition: solvers.hpp:841
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.hpp:867
void SetOperator(const Operator &op_)
Set/update the solver for the given operator.
Definition: solvers.hpp:439
ProductSolver(Operator *A_, Solver *S0_, Solver *S1_, bool ownA, bool ownS0, bool ownS1)
Definition: solvers.hpp:1148
virtual void Mult(const Vector &b, Vector &x) const
Solve the nonlinear system with right-hand side b.
Definition: solvers.cpp:1817
void SetAbsTol(double atol)
Definition: solvers.hpp:200
PrintLevel print_options
Output behavior for the iterative solver.
Definition: solvers.hpp:137
void SetRelTol(double rtol)
Definition: solvers.hpp:199
double Control[UMFPACK_CONTROL]
Definition: solvers.hpp:1064
Monitor that checks whether the residual is zero at a given set of dofs.
Definition: solvers.hpp:1036
Abstract base class for iterative solver.
Definition: solvers.hpp:66
UMFPackSolver(bool use_long_ints_=false)
For larger matrices, if the solver fails, set the parameter use_long_ints_ = true.
Definition: solvers.hpp:1069
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:3227
void DeleteStorageVectors()
Definition: solvers.hpp:732
void SetAdaptiveLinRtol(const int type=2, const double rtol0=0.5, const double rtol_max=0.9, const double alpha=0.5 *(1.0+sqrt(5.0)), const double gamma=1.0)
Enable adaptive linear solver relative tolerance algorithm.
Definition: solvers.cpp:1925
bool GetConverged() const
Returns true if the last call to Mult() converged successfully.
Definition: solvers.hpp:252
int max_iter
Limit for the number of iterations the solver is allowed to do.
Definition: solvers.hpp:151
LBFGSSolver(MPI_Comm comm_)
Definition: solvers.hpp:759
HYPRE_Int HYPRE_BigInt
double * GetBlockData()
Definition: solvers.hpp:995
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:3192
GMRES method.
Definition: solvers.hpp:525
virtual void SetSolver(Solver &solver)
Set the linear solver for inverting the Jacobian.
Definition: solvers.hpp:780
void AdaptiveLinRtolPostSolve(const Vector &x, const Vector &b, const int it, const double fnorm) const
Method for the adaptive linear solver rtol invoked after the linear solve.
Definition: solvers.cpp:1985
static int GuessLegacyPrintLevel(PrintLevel)
Use some heuristics to guess a legacy print level corresponding to the given PrintLevel.
Definition: solvers.cpp:149
FGMRESSolver(MPI_Comm comm_)
Definition: solvers.hpp:553
virtual void SetOptimizationProblem(const OptimizationProblem &prob)
Definition: solvers.cpp:2346
void SetMonitor(IterativeSolverMonitor &m)
Set the iterative solver monitor.
Definition: solvers.hpp:286
double a
Definition: lissajous.cpp:41
void UpdateVectors()
Definition: solvers.cpp:519
Abstract solver for OptimizationProblems.
Definition: solvers.hpp:848
A "square matrix" operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
Solver wrapper which orthogonalizes the input and output vector.
Definition: solvers.hpp:1163
virtual ~LBFGSSolver()
Definition: solvers.hpp:783
klu_common Common
Definition: solvers.hpp:1119
double rel_tol
Relative tolerance.
Definition: solvers.hpp:154
OptimizationProblem(int insize, const Operator *C_, const Operator *D_)
In parallel, insize is the number of the local true dofs.
Definition: solvers.cpp:2302
void SetSolutionBounds(const Vector &xl, const Vector &xh)
Definition: solvers.cpp:2330
virtual void Mult(const Vector &xt, Vector &x) const
Definition: solvers.cpp:2380
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get()...
Definition: handle.hpp:104
AuxSpaceSmoother(const HypreParMatrix &op, HypreParMatrix *aux_map, bool op_is_symmetric=true, bool own_aux_map=false)
Definition: solvers.cpp:3507
virtual void Mult(const Vector &xt, Vector &x) const =0
Operator application: y=A(x).
BiCGSTABSolver(MPI_Comm comm_)
Definition: solvers.hpp:583
void SetOperator(const Operator &op)
Recompute the diagonal using the method AssembleDiagonal of the given new Operator, op.
Definition: solvers.cpp:241
int dim
Definition: ex24.cpp:53
SLBQPOptimizer(MPI_Comm comm_)
Definition: solvers.hpp:914
const Vector * GetEqualityVec() const
Definition: solvers.hpp:838
PrintLevel FromLegacyPrintLevel(int)
Convert a legacy print level integer to a PrintLevel object.
Definition: solvers.cpp:114
int * GetBlockJ()
Definition: solvers.hpp:990
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:586
double Dot(const Vector &x, const Vector &y) const
Definition: solvers.cpp:55
MPI_Comm GetComm() const
Return the associated MPI communicator, or MPI_COMM_NULL if no communicator is set.
Definition: solvers.hpp:292
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:507
Vector data type.
Definition: vector.hpp:60
void SetSolver(Solver &s)
Set the solver used by the OrthoSolver.
Definition: solvers.cpp:3431
void Mult(const Vector &b, Vector &x) const
Perform the action of the OrthoSolver: P * solver * P where P is the projection to the subspace of ve...
Definition: solvers.cpp:3450
void MonitorResidual(int it, double norm, const Vector &r, bool final) override
Monitor the residual vector r.
Definition: solvers.cpp:3045
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.hpp:778
const Vector * d_lo
Definition: solvers.hpp:818
virtual void MonitorSolution(int it, double norm, const Vector &x, bool final)
Monitor the solution vector x.
Definition: solvers.hpp:54
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:173
KLUSolver(SparseMatrix &A)
Definition: solvers.hpp:1107
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:2147
virtual ~KLUSolver()
Definition: solvers.cpp:3345
RefCoord s[3]
double * Data()
Definition: densemat.hpp:1095
void UpdateVectors()
Definition: solvers.cpp:703
BiCGSTAB method.
Definition: solvers.hpp:572
bool warnings
If a non-fatal problem has been detected some context-specific information will be reported to mfem::...
Definition: solvers.hpp:85
bool iterations
Detailed information about each iteration will be reported to mfem::out.
Definition: solvers.hpp:88
OptimizationSolver(MPI_Comm comm_)
Definition: solvers.hpp:856
SparseMatrix * mat
Definition: solvers.hpp:1097
Base class for solvers.
Definition: operator.hpp:682
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: solvers.hpp:1137
const OptimizationProblem * problem
Definition: solvers.hpp:851
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:1376
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:699
Abstract operator.
Definition: operator.hpp:24
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:343
Block diagonal solver for A, each block is inverted by direct solver.
Definition: solvers.hpp:1125
void SetBounds(const Vector &lo_, const Vector &hi_)
Definition: solvers.cpp:2359
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: solvers.cpp:3385
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
Definition: solvers.cpp:3402
void AdaptiveLinRtolPreSolve(const Vector &x, const int it, const double fnorm) const
Method for the adaptive linear solver rtol invoked before the linear solve.
Definition: solvers.cpp:1938
const Operator * D
Definition: solvers.hpp:817
Rank 3 tensor (array of matrices)
Definition: densemat.hpp:978
virtual void SetOptimizationProblem(const OptimizationProblem &prob)
Definition: solvers.hpp:862
NewtonSolver(MPI_Comm comm_)
Definition: solvers.hpp:683
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: solvers.cpp:479
const Operator * GetD() const
Definition: solvers.hpp:837
const Vector * c_e
Definition: solvers.hpp:818
IterativeSolverMonitor * monitor
Definition: solvers.hpp:122
klu_numeric * Numeric
Definition: solvers.hpp:1099
virtual void MonitorResidual(int it, double norm, const Vector &r, bool final)
Monitor the residual vector r.
Definition: solvers.hpp:48
int width
Dimension of the input / number of columns in the matrix.
Definition: operator.hpp:28
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:475
Settings for the output behavior of the IterativeSolver.
Definition: solvers.hpp:78
const Operator * C
Not owned, some can remain unused (NULL).
Definition: solvers.hpp:817
double abs_tol
Absolute tolerance.
Definition: solvers.hpp:157
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:671
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:1338
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.cpp:1806
virtual void SetOperator(const Operator &op)
Factorize the given Operator op which must be a SparseMatrix.
Definition: solvers.cpp:3097
double GetFinalRelNorm() const
Returns the final residual norm after termination of the solver during the last call to Mult()...
Definition: solvers.hpp:271