Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RootFinderAlgorithms.h
Go to the documentation of this file.
1// @(#)root/mathmore:$Id$
2// Author: L. Moneta, A. Zsenei 08/2005
3
4 /**********************************************************************
5 * *
6 * Copyright (c) 2004 ROOT Foundation, CERN/PH-SFT *
7 * *
8 * This library is free software; you can redistribute it and/or *
9 * modify it under the terms of the GNU General Public License *
10 * as published by the Free Software Foundation; either version 2 *
11 * of the License, or (at your option) any later version. *
12 * *
13 * This library is distributed in the hope that it will be useful, *
14 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
15 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU *
16 * General Public License for more details. *
17 * *
18 * You should have received a copy of the GNU General Public License *
19 * along with this library (see file COPYING); if not, write *
20 * to the Free Software Foundation, Inc., 59 Temple Place, Suite *
21 * 330, Boston, MA 02111-1307 USA, or contact the author. *
22 * *
23 **********************************************************************/
24
25// Header file for GSL ROOT Finder Algorithms
26//
27// Created by: moneta at Sun Nov 14 14:07:50 2004
28//
29// Last update: Sun Nov 14 14:07:50 2004
30//
31#ifndef ROOT_Math_GSLRootFinderAlgorithms
32#define ROOT_Math_GSLRootFinderAlgorithms
33
34
35#include "Math/GSLRootFinder.h"
36
38
39namespace ROOT {
40namespace Math {
41
42 /**
43 Root-Finding Algorithms
44
45 */
46
47namespace Roots {
48
49//________________________________________________________________________________________________________
50 /**
51 Roots::Bisection
52 Bisection algorithm, simplest algorithm for bracketing the roots of a function, but slowest one.
53 See the <A HREF="http://www.gnu.org/software/gsl/manual/html_node/Root-Bracketing-Algorithms.html">GSL manual</A> for more information
54 @ingroup RootFinders
55 */
56
57 class Bisection : public GSLRootFinder {
58
59 public:
60
61 Bisection();
62 ~Bisection() override;
63
64 private:
65 // usually copying is non trivial, so we make this unaccessible
66
67 Bisection(const Bisection &);
69
70 };
71
72//________________________________________________________________________________________________________
73 /**
74 False Position algorithm based on linear interpolation.
75 See the <A HREF="http://www.gnu.org/software/gsl/manual/html_node/Root-Bracketing-Algorithms.html">GSL manual</A> for more information
76 @ingroup RootFinders
77 */
78
79 class FalsePos : public GSLRootFinder {
80
81 public:
82
83 FalsePos();
84 ~FalsePos() override;
85
86 private:
87 // usually copying is non trivial, so we make this unaccessible
88 FalsePos(const FalsePos &);
89 FalsePos & operator = (const FalsePos &);
90
91 };
92
93
94
95//________________________________________________________________________________________________________
96 /**
97 Brent-Dekker algorithm which combines an interpolation strategy with the bisection algorithm
98 See the <A HREF="http://www.gnu.org/software/gsl/manual/html_node/Root-Bracketing-Algorithms.html">
99 GSL manual</A> for more information
100
101 @ingroup RootFinders
102 */
103
104 class Brent : public GSLRootFinder {
105
106 public:
107
108 Brent();
109 ~Brent() override;
110
111 private:
112 // usually copying is non trivial, so we make this unaccessible
113 Brent(const Brent &);
114 Brent & operator = (const Brent &);
115
116 };
117
118
119 //----------------------------------------------------------------------
120 // algorithm with derivatives
121 //----------------------------------------------------------------------
122
123//________________________________________________________________________________________________________
124 /**
125 a Newton algorithm, which computes the derivative at each iteration
126 See the <A HREF="http://www.gnu.org/software/gsl/manual/html_node/Root-Finding-Algorithms-using-Derivatives.html">
127 GSL manual</A> for more information
128
129 @ingroup RootFinders
130 */
131
132 class Newton : public GSLRootFinderDeriv {
133
134 public:
135
136 Newton();
137 ~Newton() override;
138
139 private:
140 // usually copying is non trivial, so we make this unaccessible
141 Newton(const Newton &);
142 Newton & operator = (const Newton &);
143
144 };
145
146
147//________________________________________________________________________________________________________
148 /**
149 \a Secant algorithm, simplified version of Newton method, which does not require the derivative at every step.
150 See the <A HREF="http://www.gnu.org/software/gsl/manual/html_node/Root-Finding-Algorithms-using-Derivatives.html">
151 GSL manual</A> for more information
152 @ingroup RootFinders
153 */
154
155 class Secant : public GSLRootFinderDeriv {
156
157 public:
158
159 Secant();
160 ~Secant() override;
161
162 private:
163 // usually copying is non trivial, so we make this unaccessible
164 Secant(const Secant &);
165 Secant & operator = (const Secant &);
166
167 };
168
169//________________________________________________________________________________________________________
170 /**
171 \a Steffenson method, providing the fastes convergence.
172 See the <A HREF="http://www.gnu.org/software/gsl/manual/html_node/Root-Finding-Algorithms-using-Derivatives.html">
173 GSL manual</A> for more information
174
175 @ingroup RootFinders
176 */
177
179
180 public:
181
182 Steffenson();
183 ~Steffenson() override;
184
185 private:
186 // usually copying is non trivial, so we make this unaccessible
187 Steffenson(const Steffenson &);
189
190 };
191
192
193}
194
195} // namespace Math
196} // namespace ROOT
197
198
199#endif /* ROOT_Math_GSLRootFinderAlgorithms */
Base class for GSL Root-Finding algorithms for one dimensional functions which use function derivativ...
Base class for GSL Root-Finding algorithms for one dimensional functions which do not use function de...
Roots::Bisection Bisection algorithm, simplest algorithm for bracketing the roots of a function,...
Bisection & operator=(const Bisection &)
Brent-Dekker algorithm which combines an interpolation strategy with the bisection algorithm See the ...
Brent & operator=(const Brent &)
False Position algorithm based on linear interpolation.
FalsePos & operator=(const FalsePos &)
a Newton algorithm, which computes the derivative at each iteration See the GSL manual for more infor...
Newton & operator=(const Newton &)
Secant algorithm, simplified version of Newton method, which does not require the derivative at every...
Secant & operator=(const Secant &)
Steffenson method, providing the fastes convergence.
Steffenson & operator=(const Steffenson &)
Namespace for new Math classes and functions.
tbb::task_arena is an alias of tbb::interface7::task_arena, which doesn't allow to forward declare tb...