GetFEM  5.5
gmm_iter_solvers.h
Go to the documentation of this file.
1 /* -*- c++ -*- (enables emacs c++ mode) */
2 /*===========================================================================
3 
4  Copyright (C) 2002-2026 Yves Renard
5 
6  This file is a part of GetFEM
7 
8  GetFEM is free software; you can redistribute it and/or modify it
9  under the terms of the GNU Lesser General Public License as published
10  by the Free Software Foundation; either version 3 of the License, or
11  (at your option) any later version along with the GCC Runtime Library
12  Exception either version 3.1 or (at your option) any later version.
13  This program is distributed in the hope that it will be useful, but
14  WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
15  or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
16  License and GCC Runtime Library Exception for more details.
17  You should have received a copy of the GNU Lesser General Public License
18  along with this program. If not, see https://www.gnu.org/licenses/.
19 
20  As a special exception, you may use this file as it is a part of a free
21  software library without restriction. Specifically, if other files
22  instantiate templates or use macros or inline functions from this file,
23  or you compile this file and link it with other files to produce an
24  executable, this file does not by itself cause the resulting executable
25  to be covered by the GNU Lesser General Public License. This exception
26  does not however invalidate any other reasons why the executable file
27  might be covered by the GNU Lesser General Public License.
28 
29 ===========================================================================*/
30 
31 /**@file gmm_iter_solvers.h
32  @author Yves Renard <Yves.Renard@insa-lyon.fr>
33  @date October 13, 2002.
34  @brief Include standard gmm iterative solvers (cg, gmres, ...)
35 */
36 #ifndef GMM_ITER_SOLVERS_H__
37 #define GMM_ITER_SOLVERS_H__
38 
39 #include "gmm_iter.h"
40 
41 
42 namespace gmm {
43 
44  /** mixed method to find a zero of a real function G, a priori
45  * between a and b. If the zero is not between a and b, iterations
46  * of secant are applied. When a convenient interval is found,
47  * iterations of dichotomie and regula falsi are applied.
48  */
49  template <typename FUNC, typename T>
50  T find_root(const FUNC &G, T a = T(0), T b = T(1),
51  T tol = gmm::default_tol(T())) {
52  T c, Ga = G(a), Gb = G(b), Gc, d;
53  d = gmm::abs(b - a);
54 #if 0
55  for (int i = 0; i < 4; i++) { /* secant iterations. */
56  if (d < tol) return (b + a) / 2.0;
57  c = b - Gb * (b - a) / (Gb - Ga); Gc = G(c);
58  a = b; b = c; Ga = Gb; Gb = Gc;
59  d = gmm::abs(b - a);
60  }
61 #endif
62  while (Ga * Gb > 0.0) { /* secant iterations. */
63  if (d < tol) return (b + a) / 2.0;
64  c = b - Gb * (b - a) / (Gb - Ga); Gc = G(c);
65  a = b; b = c; Ga = Gb; Gb = Gc;
66  d = gmm::abs(b - a);
67  }
68 
69  c = std::max(a, b); a = std::min(a, b); b = c;
70  while (d > tol) {
71  c = b - (b - a) * (Gb / (Gb - Ga)); /* regula falsi. */
72  if (c > b) c = b;
73  if (c < a) c = a;
74  Gc = G(c);
75  if (Gc*Gb > 0) { b = c; Gb = Gc; } else { a = c; Ga = Gc; }
76  c = (b + a) / 2.0 ; Gc = G(c); /* Dichotomie. */
77  if (Gc*Gb > 0) { b = c; Gb = Gc; } else { a = c; Ga = Gc; }
78  d = gmm::abs(b - a); c = (b + a) / 2.0; if ((c == a) || (c == b)) d = 0.0;
79  }
80  return (b + a) / 2.0;
81  }
82 
83 }
84 
85 #include "gmm_precond_diagonal.h"
86 #include "gmm_precond_ildlt.h"
87 #include "gmm_precond_ildltt.h"
89 #include "gmm_precond_ilu.h"
90 #include "gmm_precond_ilut.h"
91 #include "gmm_precond_ilutp.h"
92 
93 
94 
95 #include "gmm_solver_cg.h"
96 #include "gmm_solver_bicgstab.h"
97 #include "gmm_solver_qmr.h"
101 #include "gmm_tri_solve.h"
102 #include "gmm_solver_gmres.h"
103 #include "gmm_solver_bfgs.h"
104 #include "gmm_least_squares_cg.h"
105 
106 // #include "gmm_solver_idgmres.h"
107 
108 
109 
110 #endif // GMM_ITER_SOLVERS_H__
Iteration object.
T find_root(const FUNC &G, T a=T(0), T b=T(1), T tol=gmm::default_tol(T()))
mixed method to find a zero of a real function G, a priori between a and b.
Modified Gram-Schmidt orthogonalization.
Diagonal matrix preconditoner.
Incomplete Level 0 ILDLT Preconditioner.
incomplete LDL^t (cholesky) preconditioner with fill-in and threshold.
Incomplete LU without fill-in Preconditioner.
ILUT: Incomplete LU with threshold and K fill-in Preconditioner.
ILUTP: Incomplete LU with threshold and K fill-in Preconditioner and column pivoting.
Approximate inverse via MR iteration.
Implements BFGS (Broyden, Fletcher, Goldfarb, Shanno) algorithm.
BiCGStab iterative solver.
Conjugate gradient iterative solver.
Constrained conjugate gradient.
GMRES (Generalized Minimum Residual) iterative solver.
Quasi-Minimal Residual iterative solver.
Solve triangular linear system for dense matrices.