ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
pm_mathRoot::getRoot Interface Reference

Generate and return a root of a specified continuous real-valued one-dimensional mathematical function such that \(f(\mathrm{root}) = 0\) with the user-specified or the default root-finding method.
More...

Detailed Description

Generate and return a root of a specified continuous real-valued one-dimensional mathematical function such that \(f(\mathrm{root}) = 0\) with the user-specified or the default root-finding method.

See the documentation of pm_mathRoot for details of the various root-finding methods.

Note
Which root-finding method should I use among all available?
The brent method is the recommended algorithm of choice for general one-dimensional root-finding problems without knowledge of function derivative (gradient).
The Newton-Raphson method is the recommended algorithm of choice for general one-dimensional root-finding problems with knowledge of function derivative (gradient).
Parameters
[in]method: The input scalar constant that can be one of the following:
  1. The constant false or an object of type false_type, signifying the use of the False-Position root-finding method within the algorithm.
  2. The constant bisection or an object of type bisection_type, signifying the use of the Bisection root-finding method within the algorithm.
  3. The constant secant or an object of type secant_type, signifying the use of the Secant root-finding method within the algorithm.
  4. The constant brent or an object of type brent_type, signifying the use of the Brent root-finding method within the algorithm.
  5. The constant ridders or an object of type ridders_type, signifying the use of the Ridders root-finding method within the algorithm.
  6. The constant toms748 or an object of type toms748_type, signifying the use of the TOMS748 root-finding method within the algorithm.
  7. The constant newton or an object of type newton_type, signifying the use of the Newton root-finding method within the algorithm.
  8. The constant halley or an object of type halley_type, signifying the use of the Halley root-finding method within the algorithm.
  9. The constant schroder or an object of type schroder_type, signifying the use of the Schroder root-finding method within the algorithm.
(optional, default = brent. Note the implicit constraint this default option sets on the input getFunc() interface below.)
[in]getFunc: The external user-specified function whose interface depends on the specified value for the input argument method.
  1. If the specified method is any of the following,
    1. the constant brent or an object of type brent_type,
    2. the constant false or an object of type false_type,
    3. the constant secant or an object of type secant_type,
    4. the constant ridders or an object of type ridders_type,
    5. the constant toms748 or an object of type toms748_type,
    6. the constant bisection or an object of type bisection_type,
    none of which require the derivative of the target function to find the function root, then getFunc() must take a single scalar input argument x of the same type and kind as the output argument root (below).
    On output, getFunc() must return a scalar of the same type and kind as the function input argument x, containing the value of the target function evaluated at the specified input x.
    The following illustrates the generic interface of getFunc() for the above values of method,
    function getFunc(x) result(func)
    real(RKG) , intent(in) :: x
    real(RKG) :: func
    end function
    where RKG refers to the kind of the output argument root.
  2. If the specified method is any of the following,
    1. the constant newton or an object of type newton_type,
    2. the constant halley or an object of type halley_type,
    3. the constant schroder or an object of type schroder_type,
    all of which either require the first (Newton) or higher derivatives (Halley/Schroder) of the target function to find the function root, then getFunc() must take two arguments:
    1. A scalar input argument x of the same type and kind as the output argument root (below).
    2. A scalar input argument order of type integer of default kind IK.
    On output, getFunc() must return a scalar of the same type and kind as the function input argument x, containing the derivative of the specified input order of the target function evaluated at the specified input x:
    1. An input order value of 0 must yield the target function value at the input x.
    2. An input order value of 1 must yield the target function first derivative at the input x.
    3. An input order value of 2 must yield the target function second derivative at the input x.
    The following illustrates the generic interface of getFunc() for the above values of method,
    function getFunc(x, order) result(func)
    real(RKG) , intent(in) :: x
    integer(IK) , intent(in) :: order
    real(RKG) :: func
    end function
    where RKG refers to the kind of the output argument root.
[in,out]root: The input/output scalar of type real of kind any supported by the processor (e.g., RK, RK32, RK64, or RK128).
On input,
  1. If the specified method is of type newton_type or halley_type, then root must contain the best guess initial value for the function root within the bracket specified by the input arguments lb and ub.
  2. Otherwise, any input value for root is ignored for all other root-finding methods.
On output, root contains the best approximation to the root of the user-specified function getFunc() within the specified search interval [lb, ub] and tolerance threshold abstol.
[in]lb: The input scalar of same type and kind as the output root, representing the lower bound of the bracket (search) interval.
[in]ub: The input scalar of same type and kind as the output root, representing the upper bound of the bracket (search) interval.
[in]abstol: The input scalar of same type and kind as the output root, representing the absolute tolerance used as the stopping criterion of the search.
The iterations of the specified input method continue until the search interval becomes smaller than abstol in absolute units.
Care must be taken for specifying a reasonable value for abstol (see the warnings below).
If no suitable value for abstol is known a priori, try abstol = epsilon(0._RKG)**.8 * (abs(lb) + abs(ub)) where RKG refers to the kind of the output argument root.
(optional, default = epsilon(0._RKG)**.8 * (abs(lb) + abs(ub)))
[out]neval: The output scalar argument of type integer of default kind IK, containing the total number of getFunc() function calls made by the algorithm.
  1. A positive output neval implies the successful convergence of the algorithm to the function root after +neval function evaluations.
  2. A negative output neval implies the failure of convergence of the algorithm to the function root after -neval function evaluations.
  3. A zero output neval only occurs if the root lies at the search bracket boundaries specified by the input argument lb or ub.
Convergence failures rarely occur. If they do, setting the input arguments abstol and/or niter to larger values may resolve the failure.
[in]niter: The input scalar of type integer of default kind IK, representing the maximum number of iterations allowed for the specified method to find the root.
The default number of steps niter within the algorithm is a compile-time constant that depends on the kind of the real input arguments.
(optional, default = ceiling(2 * precision(lb) / log10(2.)))


Possible calling interfaces

use pm_mathRoot, only: getRoot
root = getRoot(getFunc, lb, ub, abstol = abstol, neval = neval, niter = niter)
root = getRoot(method, getFunc, lb, ub, abstol = abstol, neval = neval, niter = niter)
Generate and return a root of a specified continuous real-valued one-dimensional mathematical functio...
This module contains classes and procedures for computing the roots of one-dimensional continuous mat...
Warning
The condition lb < ub must hold for the corresponding input arguments.
The condition 0 < abstol must hold for the corresponding input arguments.
The condition abstol < (ub - lb) must hold for the corresponding input arguments.
The condition abs(getFunc(lb), lf) < abstol must hold for the corresponding input arguments.
The condition abs(getFunc(ub), uf) < abstol must hold for the corresponding input arguments.
The condition sign(lf, 1.) /= sign(uf, 1.) must hold for the corresponding input arguments.
The condition 0 < niter must hold for the corresponding input arguments.
These conditions are verified only if the library is built with the preprocessor macro CHECK_ENABLED=1.
It is crucial to keep in mind that computers use a fixed number of binary digits to represent floating-point numbers.
While the user-specified function might analytically pass through zero, it is possible that its computed value is never zero, for any floating-point argument.
One must also decide on what accuracy on the root is attainable.
For example, convergence to within \(10^{−6}\) in absolute value is reasonable when the root lies near \(1\), but unachievable if the root lies near an extremely large number near floating-point overflow.
Remarks
The procedures under discussion are impure.
The procedures under discussion are recursive.
See also
getRoot
setRoot
getPolyRoot
setPolyRoot


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK
4 use pm_kind, only: RKG => RKH ! all processor kinds are supported.
6 use pm_io, only: display_type
7
8 implicit none
9
10 real(RKG) :: root
11 integer(IK) :: neval
12 type(display_type) :: disp
13 disp = display_type(file = "main.out.F90")
14
15 call disp%skip()
16 call disp%show("root = getRoot(getSin, lb = 2._RKG, ub = 4._RKG, neval = neval)")
17 root = getRoot(getSin, lb = 2._RKG, ub = 4._RKG, neval = neval)
18 call disp%show("[root, getSin(root)]")
19 call disp%show( [root, getSin(root)] )
20 call disp%show("neval")
21 call disp%show( neval )
22 call disp%skip()
23 call disp%show("root = getRoot(false, getSin, lb = 2._RKG, ub = 4._RKG, neval = neval)")
24 root = getRoot(false, getSin, lb = 2._RKG, ub = 4._RKG, neval = neval)
25 call disp%show("[root, getSin(root)]")
26 call disp%show( [root, getSin(root)] )
27 call disp%show("neval")
28 call disp%show( neval )
29 call disp%skip()
30 call disp%show("root = getRoot(bisection, getSin, lb = 2._RKG, ub = 4._RKG, neval = neval)")
31 root = getRoot(bisection, getSin, lb = 2._RKG, ub = 4._RKG, neval = neval)
32 call disp%show("[root, getSin(root)]")
33 call disp%show( [root, getSin(root)] )
34 call disp%show("neval")
35 call disp%show( neval )
36 call disp%skip()
37 call disp%show("root = getRoot(secant, getSin, lb = 2._RKG, ub = 4._RKG, neval = neval)")
38 root = getRoot(secant, getSin, lb = 2._RKG, ub = 4._RKG, neval = neval)
39 call disp%show("[root, getSin(root)]")
40 call disp%show( [root, getSin(root)] )
41 call disp%show("neval")
42 call disp%show( neval )
43 call disp%skip()
44 call disp%show("root = getRoot(brent, getSin, lb = 2._RKG, ub = 4._RKG, neval = neval)")
45 root = getRoot(brent, getSin, lb = 2._RKG, ub = 4._RKG, neval = neval)
46 call disp%show("[root, getSin(root)]")
47 call disp%show( [root, getSin(root)] )
48 call disp%show("neval")
49 call disp%show( neval )
50 call disp%skip()
51 call disp%show("root = getRoot(ridders, getSin, lb = 2._RKG, ub = 4._RKG, neval = neval)")
52 root = getRoot(ridders, getSin, lb = 2._RKG, ub = 4._RKG, neval = neval)
53 call disp%show("[root, getSin(root)]")
54 call disp%show( [root, getSin(root)] )
55 call disp%show("neval")
56 call disp%show( neval )
57 call disp%skip()
58 call disp%show("root = getRoot(toms748, getSin, lb = 2._RKG, ub = 4._RKG, neval = neval)")
59 root = getRoot(toms748, getSin, lb = 2._RKG, ub = 4._RKG, neval = neval)
60 call disp%show("[root, getSin(root)]")
61 call disp%show( [root, getSin(root)] )
62 call disp%show("neval")
63 call disp%show( neval )
64 call disp%skip()
65 call disp%show("root = getRoot(newton, getSinDiff, lb = 2._RKG, ub = 4._RKG, neval = neval)")
66 root = getRoot(newton, getSinDiff, lb = 2._RKG, ub = 4._RKG, neval = neval)
67 call disp%show("[root, getSin(root)]")
68 call disp%show( [root, getSin(root)] )
69 call disp%show("neval")
70 call disp%show( neval )
71 call disp%skip()
72 call disp%show("root = getRoot(halley, getSinDiff, lb = 2._RKG, ub = 4._RKG, neval = neval)")
73 root = getRoot(halley, getSinDiff, lb = 2._RKG, ub = 4._RKG, neval = neval)
74 call disp%show("[root, getSin(root)]")
75 call disp%show( [root, getSin(root)] )
76 call disp%show("neval")
77 call disp%show( neval )
78 call disp%skip()
79 call disp%show("root = getRoot(schroder, getSinDiff, lb = 2._RKG, ub = 4._RKG, neval = neval)")
80 root = getRoot(schroder, getSinDiff, lb = 2._RKG, ub = 4._RKG, neval = neval)
81 call disp%show("[root, getSin(root)]")
82 call disp%show( [root, getSin(root)] )
83 call disp%show("neval")
84 call disp%show( neval )
85 call disp%skip()
86
87
88
89 call disp%skip()
90 call disp%show("root = getRoot(getCos, lb = 1._RKG, ub = 3._RKG, neval = neval)")
91 root = getRoot(getCos, lb = 1._RKG, ub = 3._RKG, neval = neval)
92 call disp%show("[root, getCos(root)]")
93 call disp%show( [root, getCos(root)] )
94 call disp%show("neval")
95 call disp%show( neval )
96 call disp%skip()
97 call disp%show("root = getRoot(false, getCos, lb = 1._RKG, ub = 3._RKG, neval = neval)")
98 root = getRoot(false, getCos, lb = 1._RKG, ub = 3._RKG, neval = neval)
99 call disp%show("[root, getCos(root)]")
100 call disp%show( [root, getCos(root)] )
101 call disp%show("neval")
102 call disp%show( neval )
103 call disp%skip()
104 call disp%show("root = getRoot(bisection, getCos, lb = 1._RKG, ub = 3._RKG, neval = neval)")
105 root = getRoot(bisection, getCos, lb = 1._RKG, ub = 3._RKG, neval = neval)
106 call disp%show("[root, getCos(root)]")
107 call disp%show( [root, getCos(root)] )
108 call disp%show("neval")
109 call disp%show( neval )
110 call disp%skip()
111 call disp%show("root = getRoot(secant, getCos, lb = 1._RKG, ub = 3._RKG, neval = neval)")
112 root = getRoot(secant, getCos, lb = 1._RKG, ub = 3._RKG, neval = neval)
113 call disp%show("[root, getCos(root)]")
114 call disp%show( [root, getCos(root)] )
115 call disp%show("neval")
116 call disp%show( neval )
117 call disp%skip()
118 call disp%show("root = getRoot(brent, getCos, lb = 1._RKG, ub = 3._RKG, neval = neval)")
119 root = getRoot(brent, getCos, lb = 1._RKG, ub = 3._RKG, neval = neval)
120 call disp%show("[root, getCos(root)]")
121 call disp%show( [root, getCos(root)] )
122 call disp%show("neval")
123 call disp%show( neval )
124 call disp%skip()
125 call disp%show("root = getRoot(ridders, getCos, lb = 1._RKG, ub = 3._RKG, neval = neval)")
126 root = getRoot(ridders, getCos, lb = 1._RKG, ub = 3._RKG, neval = neval)
127 call disp%show("[root, getCos(root)]")
128 call disp%show( [root, getCos(root)] )
129 call disp%show("neval")
130 call disp%show( neval )
131 call disp%skip()
132 call disp%show("root = getRoot(toms748, getCos, lb = 1._RKG, ub = 3._RKG, neval = neval)")
133 root = getRoot(toms748, getCos, lb = 1._RKG, ub = 3._RKG, neval = neval)
134 call disp%show("[root, getCos(root)]")
135 call disp%show( [root, getCos(root)] )
136 call disp%show("neval")
137 call disp%show( neval )
138 call disp%skip()
139 call disp%show("root = getRoot(newton, getCosDiff, lb = 1._RKG, ub = 3._RKG, neval = neval)")
140 root = getRoot(newton, getCosDiff, lb = 1._RKG, ub = 3._RKG, neval = neval)
141 call disp%show("[root, getCos(root)]")
142 call disp%show( [root, getCos(root)] )
143 call disp%show("neval")
144 call disp%show( neval )
145 call disp%skip()
146 call disp%show("root = getRoot(halley, getCosDiff, lb = 1._RKG, ub = 3._RKG, neval = neval)")
147 root = getRoot(halley, getCosDiff, lb = 1._RKG, ub = 3._RKG, neval = neval)
148 call disp%show("[root, getCos(root)]")
149 call disp%show( [root, getCos(root)] )
150 call disp%show("neval")
151 call disp%show( neval )
152 call disp%skip()
153 call disp%show("root = getRoot(schroder, getCosDiff, lb = 1._RKG, ub = 3._RKG, neval = neval)")
154 root = getRoot(schroder, getCosDiff, lb = 1._RKG, ub = 3._RKG, neval = neval)
155 call disp%show("[root, getCos(root)]")
156 call disp%show( [root, getCos(root)] )
157 call disp%show("neval")
158 call disp%show( neval )
159 call disp%skip()
160
161
162
163 call disp%skip()
164 call disp%show("root = getRoot(getQuad, lb = -1._RKG, ub = 4._RKG, neval = neval)")
165 root = getRoot(getQuad, lb = -1._RKG, ub = 4._RKG, neval = neval)
166 call disp%show("[root, getQuad(root)]")
167 call disp%show( [root, getQuad(root)] )
168 call disp%show("neval")
169 call disp%show( neval )
170 call disp%skip()
171 call disp%show("root = getRoot(false, getQuad, lb = -1._RKG, ub = 4._RKG, neval = neval)")
172 root = getRoot(false, getQuad, lb = -1._RKG, ub = 4._RKG, neval = neval)
173 call disp%show("[root, getQuad(root)]")
174 call disp%show( [root, getQuad(root)] )
175 call disp%show("neval")
176 call disp%show( neval )
177 call disp%skip()
178 call disp%show("root = getRoot(bisection, getQuad, lb = -1._RKG, ub = 4._RKG, neval = neval)")
179 root = getRoot(bisection, getQuad, lb = -1._RKG, ub = 4._RKG, neval = neval)
180 call disp%show("[root, getQuad(root)]")
181 call disp%show( [root, getQuad(root)] )
182 call disp%show("neval")
183 call disp%show( neval )
184 call disp%skip()
185 call disp%show("root = getRoot(secant, getQuad, lb = -1._RKG, ub = 4._RKG, neval = neval)")
186 root = getRoot(secant, getQuad, lb = -1._RKG, ub = 4._RKG, neval = neval)
187 call disp%show("[root, getQuad(root)]")
188 call disp%show( [root, getQuad(root)] )
189 call disp%show("neval")
190 call disp%show( neval )
191 call disp%skip()
192 call disp%show("root = getRoot(brent, getQuad, lb = -1._RKG, ub = 4._RKG, neval = neval)")
193 root = getRoot(brent, getQuad, lb = -1._RKG, ub = 4._RKG, neval = neval)
194 call disp%show("[root, getQuad(root)]")
195 call disp%show( [root, getQuad(root)] )
196 call disp%show("neval")
197 call disp%show( neval )
198 call disp%skip()
199 call disp%show("root = getRoot(ridders, getQuad, lb = -1._RKG, ub = 4._RKG, neval = neval)")
200 root = getRoot(ridders, getQuad, lb = -1._RKG, ub = 4._RKG, neval = neval)
201 call disp%show("[root, getQuad(root)]")
202 call disp%show( [root, getQuad(root)] )
203 call disp%show("neval")
204 call disp%show( neval )
205 call disp%skip()
206 call disp%show("root = getRoot(toms748, getQuad, lb = -1._RKG, ub = 4._RKG, neval = neval)")
207 root = getRoot(toms748, getQuad, lb = -1._RKG, ub = 4._RKG, neval = neval)
208 call disp%show("[root, getQuad(root)]")
209 call disp%show( [root, getQuad(root)] )
210 call disp%show("neval")
211 call disp%show( neval )
212 call disp%skip()
213 call disp%show("root = getRoot(newton, getQuadDiff, lb = -1._RKG, ub = 4._RKG, neval = neval)")
214 root = getRoot(newton, getQuadDiff, lb = -1._RKG, ub = 4._RKG, neval = neval)
215 call disp%show("[root, getQuad(root)]")
216 call disp%show( [root, getQuad(root)] )
217 call disp%show("neval")
218 call disp%show( neval )
219 call disp%skip()
220 call disp%show("root = getRoot(halley, getQuadDiff, lb = -1._RKG, ub = 4._RKG, neval = neval)")
221 root = getRoot(halley, getQuadDiff, lb = -1._RKG, ub = 4._RKG, neval = neval)
222 call disp%show("[root, getQuad(root)]")
223 call disp%show( [root, getQuad(root)] )
224 call disp%show("neval")
225 call disp%show( neval )
226 call disp%skip()
227 call disp%show("root = getRoot(schroder, getQuadDiff, lb = -1._RKG, ub = 4._RKG, neval = neval)")
228 root = getRoot(schroder, getQuadDiff, lb = -1._RKG, ub = 4._RKG, neval = neval)
229 call disp%show("[root, getQuad(root)]")
230 call disp%show( [root, getQuad(root)] )
231 call disp%show("neval")
232 call disp%show( neval )
233 call disp%skip()
234
235contains
236
237 pure function getSin(x) result(func)
238 real(RKG), intent(in) :: x
239 real(RKG) :: func
240 func = sin(x)
241 end function
242
243 pure function getCos(x) result(func)
244 real(RKG), intent(in) :: x
245 real(RKG) :: func
246 func = cos(x)
247 end function
248
249 pure function getQuad(x) result(func)
250 real(RKG), intent(in) :: x
251 real(RKG) :: func
252 func = x * (x - 1._RKG) * (x - 2._RKG)
253 end function
254
255 pure function getSinDiff(x, order) result(func)
256 integer(IK), intent(in) :: order
257 real(RKG), intent(in) :: x
258 real(RKG) :: func
259 if (order == 0) func = +getSin(x)
260 if (order == 1) func = +cos(x)
261 if (order == 2) func = -sin(x)
262 end function
263
264 pure function getCosDiff(x, order) result(func)
265 integer(IK), intent(in) :: order
266 real(RKG), intent(in) :: x
267 real(RKG) :: func
268 if (order == 0) func = +getCos(x)
269 if (order == 1) func = -sin(x)
270 if (order == 2) func = -cos(x)
271 end function
272
273 pure function getQuadDiff(x, order) result(func)
274 integer(IK), intent(in) :: order
275 real(RKG), intent(in) :: x
276 real(RKG) :: func
277 if (order == 0) func = getQuad(x)
278 if (order == 1) func = 3._RKG * x**2 - 6._RKG * x + 2._RKG
279 if (order == 2) func = 6._RKG * x - 6._RKG
280 end function
281
282end program example
This is a generic method of the derived type display_type with pass attribute.
Definition: pm_io.F90:11726
This is a generic method of the derived type display_type with pass attribute.
Definition: pm_io.F90:11508
This module contains classes and procedures for input/output (IO) or generic display operations on st...
Definition: pm_io.F90:252
type(display_type) disp
This is a scalar module variable an object of type display_type for general display.
Definition: pm_io.F90:11393
This module defines the relevant Fortran kind type-parameters frequently used in the ParaMonte librar...
Definition: pm_kind.F90:268
integer, parameter LK
The default logical kind in the ParaMonte library: kind(.true.) in Fortran, kind(....
Definition: pm_kind.F90:541
integer, parameter IK
The default integer kind in the ParaMonte library: int32 in Fortran, c_int32_t in C-Fortran Interoper...
Definition: pm_kind.F90:540
integer, parameter SK
The default character kind in the ParaMonte library: kind("a") in Fortran, c_char in C-Fortran Intero...
Definition: pm_kind.F90:539
integer, parameter RKH
The scalar integer constant of intrinsic default kind, representing the highest-precision real kind t...
Definition: pm_kind.F90:858
type(brent_type), parameter brent
This is a scalar parameter object of type brent_type that is exclusively used to signify the use of B...
type(ridders_type), parameter ridders
This is a scalar parameter object of type ridders_type that is exclusively used to signify the use of...
type(toms748_type), parameter toms748
This is a scalar parameter object of type toms748_type that is exclusively used to signify the use of...
type(newton_type), parameter newton
This is a scalar parameter object of type newton_type that is exclusively used to signify the use of ...
type(bisection_type), parameter bisection
This is a scalar parameter object of type bisection_type that is exclusively used to signify the use ...
type(false_type), parameter false
This is a scalar parameter object of type false_type that is exclusively used to signify the use of F...
type(halley_type), parameter halley
This is a scalar parameter object of type halley_type that is exclusively used to signify the use of ...
type(schroder_type), parameter schroder
This is a scalar parameter object of type schroder_type that is exclusively used to signify the use o...
type(secant_type), parameter secant
This is a scalar parameter object of type secant_type that is exclusively used to signify the use of ...
Generate and return an object of type display_type.
Definition: pm_io.F90:10282

Example Unix compile command via Intel ifort compiler
1#!/usr/bin/env sh
2rm main.exe
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example Windows Batch compile command via Intel ifort compiler
1del main.exe
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
4main.exe

Example Unix / MinGW compile command via GNU gfortran compiler
1#!/usr/bin/env sh
2rm main.exe
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example output
1
2root = getRoot(getSin, lb = 2._RKG, ub = 4._RKG, neval = neval)
3[root, getSin(root)]
4+3.14159265358979323846264338327950280, +0.867181013012378102479704402604335225E-34
5neval
6+9
7
8root = getRoot(false, getSin, lb = 2._RKG, ub = 4._RKG, neval = neval)
9[root, getSin(root)]
10+3.14159265358979323846264338327950280, +0.867181013012378102479704402604335225E-34
11neval
12+9
13
14root = getRoot(bisection, getSin, lb = 2._RKG, ub = 4._RKG, neval = neval)
15[root, getSin(root)]
16+3.14159265358979323846264338592990596, -0.265040307298939228819052880532077411E-26
17neval
18+93
19
20root = getRoot(secant, getSin, lb = 2._RKG, ub = 4._RKG, neval = neval)
21[root, getSin(root)]
22+3.14159265358979323846264338327950280, +0.867181013012378102479704402604335225E-34
23neval
24+9
25
26root = getRoot(brent, getSin, lb = 2._RKG, ub = 4._RKG, neval = neval)
27[root, getSin(root)]
28+3.14159265358979323846264338327950280, +0.867181013012378102479704402604335225E-34
29neval
30+9
31
32root = getRoot(ridders, getSin, lb = 2._RKG, ub = 4._RKG, neval = neval)
33[root, getSin(root)]
34+3.14159265358979323846264338327950280, +0.867181013012378102479704402604335225E-34
35neval
36+22
37
38root = getRoot(toms748, getSin, lb = 2._RKG, ub = 4._RKG, neval = neval)
39[root, getSin(root)]
40+3.14159265358979323846264338327950280, +0.867181013012378102479704402604335225E-34
41neval
42+9
43
44root = getRoot(newton, getSinDiff, lb = 2._RKG, ub = 4._RKG, neval = neval)
45[root, getSin(root)]
46+3.14159265358979323846264338327950280, +0.867181013012378102479704402604335225E-34
47neval
48+12
49
50root = getRoot(halley, getSinDiff, lb = 2._RKG, ub = 4._RKG, neval = neval)
51[root, getSin(root)]
52+3.14159265358979323846264338327950280, +0.867181013012378102479704402604335225E-34
53neval
54+14
55
56root = getRoot(schroder, getSinDiff, lb = 2._RKG, ub = 4._RKG, neval = neval)
57[root, getSin(root)]
58+3.14159265358979323846264338327950280, +0.867181013012378102479704402604335225E-34
59neval
60+14
61
62
63root = getRoot(getCos, lb = 1._RKG, ub = 3._RKG, neval = neval)
64[root, getCos(root)]
65+1.57079632679489661923132169163975140, +0.433590506506189051239852201302167613E-34
66neval
67+9
68
69root = getRoot(false, getCos, lb = 1._RKG, ub = 3._RKG, neval = neval)
70[root, getCos(root)]
71+1.57079632679489661923132169163975140, +0.433590506506189051239852201302167613E-34
72neval
73+9
74
75root = getRoot(bisection, getCos, lb = 1._RKG, ub = 3._RKG, neval = neval)
76[root, getCos(root)]
77+1.57079632679489661923132169458054011, -0.294078867038732832157848450435953325E-26
78neval
79+93
80
81root = getRoot(secant, getCos, lb = 1._RKG, ub = 3._RKG, neval = neval)
82[root, getCos(root)]
83+1.57079632679489661923132169163975140, +0.433590506506189051239852201302167613E-34
84neval
85+9
86
87root = getRoot(brent, getCos, lb = 1._RKG, ub = 3._RKG, neval = neval)
88[root, getCos(root)]
89+1.57079632679489661923132169163975140, +0.433590506506189051239852201302167613E-34
90neval
91+9
92
93root = getRoot(ridders, getCos, lb = 1._RKG, ub = 3._RKG, neval = neval)
94[root, getCos(root)]
95+1.57079632679489661923132169164015931, -0.407868603170565934772132143019357364E-30
96neval
97+22
98
99root = getRoot(toms748, getCos, lb = 1._RKG, ub = 3._RKG, neval = neval)
100[root, getCos(root)]
101+1.57079632679489661923132169163975140, +0.433590506506189051239852201302167613E-34
102neval
103+10
104
105root = getRoot(newton, getCosDiff, lb = 1._RKG, ub = 3._RKG, neval = neval)
106[root, getCos(root)]
107+1.57079632679489661923132169163975140, +0.433590506506189051239852201302167613E-34
108neval
109+12
110
111root = getRoot(halley, getCosDiff, lb = 1._RKG, ub = 3._RKG, neval = neval)
112[root, getCos(root)]
113+1.57079632679489661923132169163975140, +0.433590506506189051239852201302167613E-34
114neval
115+17
116
117root = getRoot(schroder, getCosDiff, lb = 1._RKG, ub = 3._RKG, neval = neval)
118[root, getCos(root)]
119+1.57079632679489661923132169163975140, +0.433590506506189051239852201302167613E-34
120neval
121+17
122
123
124root = getRoot(getQuad, lb = -1._RKG, ub = 4._RKG, neval = neval)
125[root, getQuad(root)]
126+0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000
127neval
128+3
129
130root = getRoot(false, getQuad, lb = -1._RKG, ub = 4._RKG, neval = neval)
131[root, getQuad(root)]
132+0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000
133neval
134+3
135
136root = getRoot(bisection, getQuad, lb = -1._RKG, ub = 4._RKG, neval = neval)
137[root, getQuad(root)]
138+1.99999999999999999999999999838441287, -0.323117426778526435496644019556792704E-26
139neval
140+94
141
142root = getRoot(secant, getQuad, lb = -1._RKG, ub = 4._RKG, neval = neval)
143[root, getQuad(root)]
144+0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000
145neval
146+3
147
148root = getRoot(brent, getQuad, lb = -1._RKG, ub = 4._RKG, neval = neval)
149[root, getQuad(root)]
150+0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000
151neval
152+3
153
154root = getRoot(ridders, getQuad, lb = -1._RKG, ub = 4._RKG, neval = neval)
155[root, getQuad(root)]
156+2.00000000000000000000000000000000000, +0.00000000000000000000000000000000000
157neval
158+24
159
160root = getRoot(toms748, getQuad, lb = -1._RKG, ub = 4._RKG, neval = neval)
161[root, getQuad(root)]
162+2.00000000000000000000000000000000000, +0.00000000000000000000000000000000000
163neval
164+4
165
166root = getRoot(newton, getQuadDiff, lb = -1._RKG, ub = 4._RKG, neval = neval)
167[root, getQuad(root)]
168+0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000
169neval
170+6
171
172root = getRoot(halley, getQuadDiff, lb = -1._RKG, ub = 4._RKG, neval = neval)
173[root, getQuad(root)]
174+1.00000000000000000000000000000000000, -0.00000000000000000000000000000000000
175neval
176+20
177
178root = getRoot(schroder, getQuadDiff, lb = -1._RKG, ub = 4._RKG, neval = neval)
179[root, getQuad(root)]
180+1.00000000000000000000000000000000000, -0.00000000000000000000000000000000000
181neval
182+20
183
184
Test:
test_pm_mathRoot


Final Remarks


If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.

  1. If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
  2. If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.

This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.

Author:
Amir Shahmoradi, Oct 16, 2009, 11:14 AM, Michigan

Definition at line 1589 of file pm_mathRoot.F90.


The documentation for this interface was generated from the following file: