GetFEM  5.4.2
gmm_iter_solvers.h
Go to the documentation of this file.
1 /* -*- c++ -*- (enables emacs c++ mode) */
2 /*===========================================================================
3 
4  Copyright (C) 2002-2020 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, write to the Free Software Foundation,
19  Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301, USA.
20 
21  As a special exception, you may use this file as it is a part of a free
22  software library without restriction. Specifically, if other files
23  instantiate templates or use macros or inline functions from this file,
24  or you compile this file and link it with other files to produce an
25  executable, this file does not by itself cause the resulting executable
26  to be covered by the GNU Lesser General Public License. This exception
27  does not however invalidate any other reasons why the executable file
28  might be covered by the GNU Lesser General Public License.
29 
30 ===========================================================================*/
31 
32 /**@file gmm_iter_solvers.h
33  @author Yves Renard <Yves.Renard@insa-lyon.fr>
34  @date October 13, 2002.
35  @brief Include standard gmm iterative solvers (cg, gmres, ...)
36 */
37 #ifndef GMM_ITER_SOLVERS_H__
38 #define GMM_ITER_SOLVERS_H__
39 
40 #include "gmm_iter.h"
41 
42 
43 namespace gmm {
44 
45  /** mixed method to find a zero of a real function G, a priori
46  * between a and b. If the zero is not between a and b, iterations
47  * of secant are applied. When a convenient interval is found,
48  * iterations of dichotomie and regula falsi are applied.
49  */
50  template <typename FUNC, typename T>
51  T find_root(const FUNC &G, T a = T(0), T b = T(1),
52  T tol = gmm::default_tol(T())) {
53  T c, Ga = G(a), Gb = G(b), Gc, d;
54  d = gmm::abs(b - a);
55 #if 0
56  for (int i = 0; i < 4; i++) { /* secant iterations. */
57  if (d < tol) return (b + a) / 2.0;
58  c = b - Gb * (b - a) / (Gb - Ga); Gc = G(c);
59  a = b; b = c; Ga = Gb; Gb = Gc;
60  d = gmm::abs(b - a);
61  }
62 #endif
63  while (Ga * Gb > 0.0) { /* secant iterations. */
64  if (d < tol) return (b + a) / 2.0;
65  c = b - Gb * (b - a) / (Gb - Ga); Gc = G(c);
66  a = b; b = c; Ga = Gb; Gb = Gc;
67  d = gmm::abs(b - a);
68  }
69 
70  c = std::max(a, b); a = std::min(a, b); b = c;
71  while (d > tol) {
72  c = b - (b - a) * (Gb / (Gb - Ga)); /* regula falsi. */
73  if (c > b) c = b;
74  if (c < a) c = a;
75  Gc = G(c);
76  if (Gc*Gb > 0) { b = c; Gb = Gc; } else { a = c; Ga = Gc; }
77  c = (b + a) / 2.0 ; Gc = G(c); /* Dichotomie. */
78  if (Gc*Gb > 0) { b = c; Gb = Gc; } else { a = c; Ga = Gc; }
79  d = gmm::abs(b - a); c = (b + a) / 2.0; if ((c == a) || (c == b)) d = 0.0;
80  }
81  return (b + a) / 2.0;
82  }
83 
84 }
85 
86 #include "gmm_precond_diagonal.h"
87 #include "gmm_precond_ildlt.h"
88 #include "gmm_precond_ildltt.h"
90 #include "gmm_precond_ilu.h"
91 #include "gmm_precond_ilut.h"
92 #include "gmm_precond_ilutp.h"
93 
94 
95 
96 #include "gmm_solver_cg.h"
97 #include "gmm_solver_bicgstab.h"
98 #include "gmm_solver_qmr.h"
102 #include "gmm_tri_solve.h"
103 #include "gmm_solver_gmres.h"
104 #include "gmm_solver_bfgs.h"
105 #include "gmm_least_squares_cg.h"
106 
107 // #include "gmm_solver_idgmres.h"
108 
109 
110 
111 #endif // GMM_ITER_SOLVERS_H__
gmm_tri_solve.h
Solve triangular linear system for dense matrices.
gmm_solver_constrained_cg.h
Constrained conjugate gradient.
gmm_solver_bicgstab.h
BiCGStab iterative solver.
gmm_solver_bfgs.h
Implements BFGS (Broyden, Fletcher, Goldfarb, Shanno) algorithm.
gmm_modified_gram_schmidt.h
Modified Gram-Schmidt orthogonalization.
gmm_precond_diagonal.h
Diagonal matrix preconditoner.
gmm::find_root
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.
Definition: gmm_iter_solvers.h:51
gmm_solver_cg.h
Conjugate gradient iterative solver.
gmm_precond_mr_approx_inverse.h
Approximate inverse via MR iteration.
gmm_precond_ildltt.h
incomplete LDL^t (cholesky) preconditioner with fill-in and threshold.
gmm_precond_ilut.h
ILUT: Incomplete LU with threshold and K fill-in Preconditioner.
gmm_solver_gmres.h
GMRES (Generalized Minimum Residual) iterative solver.
gmm_solver_qmr.h
Quasi-Minimal Residual iterative solver.
gmm_precond_ildlt.h
Incomplete Level 0 ILDLT Preconditioner.
gmm_precond_ilutp.h
ILUTP: Incomplete LU with threshold and K fill-in Preconditioner and column pivoting.
gmm_precond_ilu.h
Incomplete LU without fill-in Preconditioner.
gmm_iter.h
Iteration object.
gmm_solver_Schwarz_additive.h

Rabisu Mirror Service We provide mirrors to support Open source communities. Our mirror server is located in Istanbul/Turkey region.

Please do not hesitate to contact mirror@rabisu.com for new open source mirror submissions.