Generate and return the roots of a polynomial of arbitrary degree specified by its coefficients coef
.
More...
Generate and return the roots of a polynomial of arbitrary degree specified by its coefficients coef
.
See the documentation of pm_polynomial for details of the root-finding methods.
- Parameters
-
[in] | coef | : The input contiguous vector of,
-
type
complex of kind any supported by the processor (e.g., CK, CK32, CK64, or CK128), or
-
type
real of kind any supported by the processor (e.g., RK, RK32, RK64, or RK128),
containing the coefficients of the polynomial in the order of increasing power.
By definition, the degree of the polynomial is size(coef) - 1 .
|
[in] | method | : The input scalar constant that can be any of the following:
-
the scalar constant eigen or a constant object of type eigen_type implying the use of the Eigenvalue root-finding method.
-
the scalar constant jenkins or a constant object of type jenkins_type implying the use of the Jenkins-Traub root-finding method.
-
the scalar constant laguerre or a constant object of type laguerre_type implying the use of the Laguerre root-finding method.
-
the scalar constant sgl or a constant object of type sgl_type implying the use of the Skowron-Gould root-finding method.
This method is yet to be fully implemented.
Which polynomial root-finding method should I use?
-
If you have a polynomial of highly varying coefficients, then the Eigenvalue method as specified by eigen_type is likely more reliable.
-
The Jenkins-Traub is also considered a relatively reliable fast Sure-Fire technique for finding the roots of polynomials.
(optional, default = eigen_type) |
- Returns
root
: The output allocatable
vector of type complex
of the same kind as the input coef
, containing the roots of the polynomial specified by its coefficients coef
.
By definition, the size of the output root is the same as the degree of the polynomial (i.e., size(coef) - 1
).
However, if the algorithm fails to find any of the roots at any stage, it will return only the roots found.
An output root
of zero size indicates the failure of the algorithm to converge or find any of the polynomial roots.
Possible calling interfaces ⛓
complex(kind(coef)), allocatable :: root(:)
root
= getPolyRoot(coef(:), method) ! allocatable output.
Generate and return the roots of a polynomial of arbitrary degree specified by its coefficients coef.
This module contains procedures and generic interfaces for performing various mathematical operations...
type(eigen_type), parameter eigen
This is a scalar parameter object of type eigen_type that is exclusively used to signify the use of E...
type(jenkins_type), parameter jenkins
This is a scalar parameter object of type jenkins_type that is exclusively used to signify the use of...
type(laguerre_type), parameter laguerre
This is a scalar parameter object of type laguerre_type that is exclusively used to signify the use o...
type(sgl_type), parameter sgl
This is a scalar parameter object of type sgl_type that is exclusively used to signify the use of Sko...
- Warning
- All warnings and conditions associated with setPolyRoot also apply to this generic interface.
These conditions are verified only if the library is built with the preprocessor macro CHECK_ENABLED=1
.
-
The
pure
procedure(s) documented herein become impure
when the ParaMonte library is compiled with preprocessor macro CHECK_ENABLED=1
.
By default, these procedures are pure
in release
build and impure
in debug
and testing
builds. The procedures of this generic interface are always impure
when the input argument method
is set to an object of type eigen_type.
- See also
- getRoot
setRoot
getPolyRoot
setPolyRoot
Example usage ⛓
11 complex(CKG),
allocatable :: coef(:), root(:)
13 type(display_type) :: disp
17 call disp%show(
"coef = [complex(CKG) :: -4, 1, -4, 1]")
18 coef
= [
complex(CKG) ::
-4,
1,
-4,
1]
21 call disp%show(
"root = getPolyRoot(coef%re) ! polynomial with real coefficients.")
25 call disp%show(
"getPolyVal(coef, root)")
28 call disp%show(
"root = getPolyRoot(coef%re, method = eigen) ! polynomial with real coefficients.")
32 call disp%show(
"getPolyVal(coef, root)")
35 call disp%show(
"root = getPolyRoot(coef%re, method = jenkins) ! polynomial with real coefficients.")
39 call disp%show(
"getPolyVal(coef, root)")
42 call disp%show(
"root = getPolyRoot(coef%re, method = laguerre) ! polynomial with real coefficients.")
46 call disp%show(
"getPolyVal(coef, root)")
51 call disp%show(
"coef = [complex(CKG) :: -120, 274, -225, 85, -15, 1]")
52 coef
= [
complex(CKG) ::
-120,
274,
-225,
85,
-15,
1]
55 call disp%show(
"root = getPolyRoot(coef%re) ! polynomial with real coefficients.")
59 call disp%show(
"getPolyVal(coef, root)")
62 call disp%show(
"root = getPolyRoot(coef%re, method = eigen) ! polynomial with real coefficients.")
66 call disp%show(
"getPolyVal(coef, root)")
69 call disp%show(
"root = getPolyRoot(coef%re, method = jenkins) ! polynomial with real coefficients.")
73 call disp%show(
"getPolyVal(coef, root)")
76 call disp%show(
"root = getPolyRoot(coef%re, method = laguerre) ! polynomial with real coefficients.")
80 call disp%show(
"getPolyVal(coef, root)")
85 call disp%show(
"coef = [complex(CKG) :: 8, -8, 16, -16, 8, -8]")
86 coef
= [
complex(CKG) ::
8,
-8,
16,
-16,
8,
-8]
89 call disp%show(
"root = getPolyRoot(coef%re) ! polynomial with real coefficients.")
93 call disp%show(
"getPolyVal(coef, root)")
96 call disp%show(
"root = getPolyRoot(coef%re, method = eigen) ! polynomial with real coefficients.")
100 call disp%show(
"getPolyVal(coef, root)")
103 call disp%show(
"root = getPolyRoot(coef%re, method = jenkins) ! polynomial with real coefficients.")
107 call disp%show(
"getPolyVal(coef, root)")
110 call disp%show(
"root = getPolyRoot(coef%re, method = laguerre) ! polynomial with real coefficients.")
114 call disp%show(
"getPolyVal(coef, root)")
119 call disp%show(
"coef = [complex(CKG) :: -120, 274, -225, 85, -15, 1]")
120 coef
= [
complex(CKG) ::
-120,
274,
-225,
85,
-15,
1]
123 call disp%show(
"root = getPolyRoot(coef)")
127 call disp%show(
"getPolyVal(coef, root)")
130 call disp%show(
"root = getPolyRoot(coef, method = eigen)")
134 call disp%show(
"getPolyVal(coef, root)")
137 call disp%show(
"root = getPolyRoot(coef, method = jenkins)")
141 call disp%show(
"getPolyVal(coef, root)")
144 call disp%show(
"root = getPolyRoot(coef, method = laguerre)")
148 call disp%show(
"getPolyVal(coef, root)")
This is a generic method of the derived type display_type with pass attribute.
This is a generic method of the derived type display_type with pass attribute.
Generate and return the value of the polynomial of arbitrary degree whose coefficients are specified ...
This module contains classes and procedures for input/output (IO) or generic display operations on st...
type(display_type) disp
This is a scalar module variable an object of type display_type for general display.
This module defines the relevant Fortran kind type-parameters frequently used in the ParaMonte librar...
integer, parameter LK
The default logical kind in the ParaMonte library: kind(.true.) in Fortran, kind(....
integer, parameter CKS
The single-precision complex kind in Fortran mode. On most platforms, this is a 32-bit real kind.
integer, parameter IK
The default integer kind in the ParaMonte library: int32 in Fortran, c_int32_t in C-Fortran Interoper...
integer, parameter SK
The default character kind in the ParaMonte library: kind("a") in Fortran, c_char in C-Fortran Intero...
Generate and return an object of type display_type.
Example Unix compile command via Intel ifort
compiler ⛓
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
Example Windows Batch compile command via Intel ifort
compiler ⛓
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
Example Unix / MinGW compile command via GNU gfortran
compiler ⛓
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
Example output ⛓
2coef
= [
complex(CKG) ::
-4,
1,
-4,
1]
4(
-4.00000000,
+0.00000000), (
+1.00000000,
+0.00000000), (
-4.00000000,
+0.00000000), (
+1.00000000,
+0.00000000)
7(
+3.99999952,
+0.00000000), (
+0.00000000,
+0.999999821), (
+0.00000000,
-0.999999821)
9(
-0.810623169E-5,
+0.00000000), (
-0.143051147E-5,
+0.357627812E-6), (
-0.143051147E-5,
-0.357627812E-6)
13(
+3.99999952,
+0.00000000), (
+0.00000000,
+0.999999821), (
+0.00000000,
-0.999999821)
15(
-0.810623169E-5,
+0.00000000), (
-0.143051147E-5,
+0.357627812E-6), (
-0.143051147E-5,
-0.357627812E-6)
19(
-0.894069672E-6,
+1.00000000), (
-0.894069672E-6,
-1.00000000), (
+4.00000191,
+0.00000000)
21(
+0.190734863E-5,
+0.715255919E-5), (
+0.190734863E-5,
-0.715255919E-5), (
+0.324249268E-4,
+0.00000000)
25(
+0.977888703E-8,
+1.00000000), (
-0.467116479E-8,
-1.00000000), (
+0.160798663E-8,
-1.00000000)
27(
+0.00000000,
-0.391155481E-7), (
+0.00000000,
-0.186846592E-7), (
+0.00000000,
+0.643194653E-8)
30coef
= [
complex(CKG) ::
-120,
274,
-225,
85,
-15,
1]
32(
-120.000000,
+0.00000000), (
+274.000000,
+0.00000000), (
-225.000000,
+0.00000000), (
+85.0000000,
+0.00000000), (
-15.0000000,
+0.00000000), (
+1.00000000,
+0.00000000)
35(
+4.99990797,
+0.00000000), (
+4.00018120,
+0.00000000), (
+2.99988103,
+0.00000000), (
+2.00003195,
+0.00000000), (
+0.999997020,
+0.00000000)
37(
-0.205993652E-2,
+0.00000000), (
-0.915527344E-3,
+0.00000000), (
-0.411987305E-3,
+0.00000000), (
-0.190734863E-3,
+0.00000000), (
-0.839233398E-4,
+0.00000000)
41(
+4.99990797,
+0.00000000), (
+4.00018120,
+0.00000000), (
+2.99988103,
+0.00000000), (
+2.00003195,
+0.00000000), (
+0.999997020,
+0.00000000)
43(
-0.205993652E-2,
+0.00000000), (
-0.915527344E-3,
+0.00000000), (
-0.411987305E-3,
+0.00000000), (
-0.190734863E-3,
+0.00000000), (
-0.839233398E-4,
+0.00000000)
47(
+0.999999523,
+0.00000000), (
+2.00003052,
+0.00000000), (
+2.99990487,
+0.00000000), (
+4.00010586,
+0.00000000), (
+4.99995899,
+0.00000000)
49(
-0.152587891E-4,
+0.00000000), (
-0.183105469E-3,
+0.00000000), (
-0.373840332E-3,
+0.00000000), (
-0.671386719E-3,
+0.00000000), (
-0.106048584E-2,
+0.00000000)
53(
+0.999993324,
+0.00000000), (
+2.00000691,
+0.931322575E-9), (
+1.99999917,
+0.291038305E-10), (
+0.999999583,
+0.349245965E-9), (
+1.00000000,
+0.593718141E-8)
55(
-0.175476074E-3,
+0.00000000), (
-0.457763672E-4,
-0.558787860E-8), (
+0.762939453E-5,
-0.174622983E-9), (
-0.228881836E-4,
+0.838191383E-8), (
+0.00000000,
+0.142492354E-6)
58coef
= [
complex(CKG) ::
8,
-8,
16,
-16,
8,
-8]
60(
+8.00000000,
+0.00000000), (
-8.00000000,
+0.00000000), (
+16.0000000,
+0.00000000), (
-16.0000000,
+0.00000000), (
+8.00000000,
+0.00000000), (
-8.00000000,
+0.00000000)
63(
+0.999999642,
+0.00000000), (
-0.429153442E-4,
+1.00021696), (
-0.429153442E-4,
-1.00021696), (
+0.427067280E-4,
+0.999782741), (
+0.427067280E-4,
-0.999782741)
65(
+0.114440918E-4,
+0.00000000), (
+0.238418579E-5,
-0.804837327E-6), (
+0.238418579E-5,
+0.804837327E-6), (
+0.190734863E-5,
-0.566709787E-6), (
+0.190734863E-5,
+0.566709787E-6)
69(
+0.999999642,
+0.00000000), (
-0.429153442E-4,
+1.00021696), (
-0.429153442E-4,
-1.00021696), (
+0.427067280E-4,
+0.999782741), (
+0.427067280E-4,
-0.999782741)
71(
+0.114440918E-4,
+0.00000000), (
+0.238418579E-5,
-0.804837327E-6), (
+0.238418579E-5,
+0.804837327E-6), (
+0.190734863E-5,
-0.566709787E-6), (
+0.190734863E-5,
+0.566709787E-6)
75(
-0.339653343E-5,
+1.00000632), (
-0.339653343E-5,
-1.00000632), (
+0.339158305E-5,
+0.999993563), (
+0.339158305E-5,
-0.999993563), (
+1.00000000,
+0.00000000)
77(
+0.00000000,
-0.469048246E-6), (
+0.00000000,
+0.469048246E-6), (
+0.00000000,
+0.430129148E-6), (
+0.00000000,
-0.430129148E-6), (
+0.00000000,
+0.00000000)
81(
+0.111623813E-3,
+0.999268234), (
-0.395235606E-3,
-0.999610364), (
+0.317461498E-3,
+0.999678612), (
+1.00000000,
-0.395812094E-8), (
+0.999999940,
+0.162981451E-8)
83(
+0.219345093E-4,
-0.112091075E-4), (
-0.953674316E-5,
+0.959518366E-5), (
+0.667572021E-5,
+0.645453110E-5), (
+0.00000000,
+0.126659870E-6), (
+0.143051147E-5,
-0.521540500E-7)
86coef
= [
complex(CKG) ::
-120,
274,
-225,
85,
-15,
1]
88(
-120.000000,
+0.00000000), (
+274.000000,
+0.00000000), (
-225.000000,
+0.00000000), (
+85.0000000,
+0.00000000), (
-15.0000000,
+0.00000000), (
+1.00000000,
+0.00000000)
91(
+4.99997854,
-0.120472396E-5), (
+4.00004864,
+0.389384786E-5), (
+2.99997115,
-0.376791718E-5), (
+2.00000644,
+0.113765248E-5), (
+0.999999404,
+0.159270712E-6)
93(
-0.137329102E-3,
-0.289111576E-4), (
-0.129699707E-3,
-0.233645478E-4), (
-0.190734863E-3,
-0.150715496E-4), (
-0.762939453E-5,
-0.682582322E-5), (
-0.762939453E-5,
+0.382250619E-5)
97(
+4.99997854,
-0.120472396E-5), (
+4.00004864,
+0.389384786E-5), (
+2.99997115,
-0.376791718E-5), (
+2.00000644,
+0.113765248E-5), (
+0.999999404,
+0.159270712E-6)
99(
-0.137329102E-3,
-0.289111576E-4), (
-0.129699707E-3,
-0.233645478E-4), (
-0.190734863E-3,
-0.150715496E-4), (
-0.762939453E-5,
-0.682582322E-5), (
-0.762939453E-5,
+0.382250619E-5)
103(
+1.00000048,
-0.684522092E-7), (
+2.00026107,
-0.101812184E-3), (
+2.99921584,
+0.304903486E-3), (
+4.00078106,
-0.304773916E-3), (
+4.99974155,
+0.101751066E-3)
105(
+0.762939453E-5,
-0.164284984E-5), (
-0.154876709E-2,
+0.610602554E-3), (
-0.317382812E-2,
+0.121960416E-2), (
-0.476837158E-2,
+0.183103979E-2), (
-0.605010986E-2,
+0.243940251E-2)
109(
+0.999993324,
+0.00000000), (
+2.00000691,
+0.931322575E-9), (
+1.99999917,
+0.291038305E-10), (
+0.999999583,
+0.349245965E-9), (
+1.00000000,
+0.593718141E-8)
111(
-0.175476074E-3,
+0.00000000), (
-0.457763672E-4,
-0.558787860E-8), (
+0.762939453E-5,
-0.174622983E-9), (
-0.228881836E-4,
+0.838191383E-8), (
+0.00000000,
+0.142492354E-6)
- Benchmarks:
Benchmark :: The effects of method
on runtime efficiency ⛓
- The following program compares the runtime performance of setPolyRoot using different polynomial root finding algorithms.
8 use iso_fortran_env,
only:
error_unit
12 integer(IK) ,
parameter :: degmax
= 9_IK
13 integer(IK) ,
allocatable :: rootCount(:)
17 integer(IK) :: fileUnit
18 integer(IK) :: rootUnit
19 real(RKG) :: dumsum
= 0._RKG
20 complex(RKG) :: coef(
0:
2**degmax)
21 complex(RKG) :: root(
2**degmax)
22 type(bench_type) ,
allocatable :: bench(:)
24 bench
= [
bench_type(name
= SK_
"Eigen", exec
= setPolyRootEigen, overhead
= setOverhead)
&
25 ,
bench_type(name
= SK_
"Jenkins", exec
= setPolyRootJenkins, overhead
= setOverhead)
&
26 ,
bench_type(name
= SK_
"Laguerre", exec
= setPolyRootLaguerre, overhead
= setOverhead)
&
27 ,
bench_type(name
= SK_
"SGL", exec
= setPolyRootSGL, overhead
= setOverhead)
&
29 allocate(rootCount(
size(bench)))
31 write(
*,
"(*(g0,:,' '))")
32 write(
*,
"(*(g0,:,' '))")
"Benchmarking setPolyRoot()"
33 write(
*,
"(*(g0,:,' '))")
35 open(newunit
= fileUnit, file
= "main.out", status
= "replace")
36 open(newunit
= rootUnit, file
= "root.count", status
= "replace")
38 write(fileUnit,
"(*(g0,:,','))")
"Polynomial Degree", (bench(ibench)
%name, ibench
= 1,
size(bench))
39 write(rootUnit,
"(*(g0,:,','))")
"Polynomial Degree", (bench(ibench)
%name, ibench
= 1,
size(bench))
40 loopOverArraySize:
do ideg
= 0, degmax
43 call setUnifRand(coef(
0 : degree), (
1._RKG,
1._RKG), (
2._RKG,
2._RKG))
44 write(
*,
"(*(g0,:,' '))")
"Benchmarking with coef size", degree
45 do ibench
= 1,
size(bench)
46 bench(ibench)
%timing
= bench(ibench)
%getTiming(minsec
= 0.07_RK)
48 write(rootUnit,
"(*(g0,:,','))") degree, rootCount
49 write(fileUnit,
"(*(g0,:,','))") degree, (bench(ibench)
%timing
%mean, ibench
= 1,
size(bench))
51 end do loopOverArraySize
52 write(
*,
"(*(g0,:,' '))") dumsum
53 write(
*,
"(*(g0,:,' '))")
63 subroutine setOverhead()
68 dumsum
= dumsum
+ rootCount(ibench)
71 subroutine setPolyRootEigen()
73 call setPolyRoot(root(
1 : degree), rootCount(ibench), coef(
0 : degree), method)
77 subroutine setPolyRootJenkins()
79 call setPolyRoot(root(
1 : degree), rootCount(ibench), coef(
0 : degree), method)
83 subroutine setPolyRootLaguerre()
85 call setPolyRoot(root(
1 : degree), rootCount(ibench), coef(
0 : degree), method)
89 subroutine setPolyRootSGL()
91 call setPolyRoot(root(
1 : degree), rootCount(ibench), coef(
0 : degree), method)
92 rootCount(ibench)
= root(
1)
%re
Allocate or resize (shrink or expand) an input allocatable scalar string or array of rank 1....
Generate and return an object of type timing_type containing the benchmark timing information and sta...
Return a uniform random scalar or contiguous array of arbitrary rank of randomly uniformly distribute...
Return the roots of a polynomial of arbitrary degree specified by its coefficients coef.
This module contains procedures and generic interfaces for resizing allocatable arrays of various typ...
This module contains abstract interfaces and types that facilitate benchmarking of different procedur...
This module contains classes and procedures for computing various statistical quantities related to t...
integer, parameter RK
The default real kind in the ParaMonte library: real64 in Fortran, c_double in C-Fortran Interoperati...
integer, parameter RKD
The double precision real kind in Fortran mode. On most platforms, this is an 64-bit real kind.
integer, parameter RKH
The scalar integer constant of intrinsic default kind, representing the highest-precision real kind t...
This is the class for creating benchmark and performance-profiling objects.
Example Unix compile command via Intel ifort
compiler ⛓
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
Example Windows Batch compile command via Intel ifort
compiler ⛓
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
Example Unix / MinGW compile command via GNU gfortran
compiler ⛓
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
Postprocessing of the benchmark output ⛓
3import matplotlib.pyplot
as plt
8dirname = os.path.basename(os.getcwd())
12df = pd.read_csv(
"main.out", delimiter =
",")
13colnames = list(df.columns.values)
19ax = plt.figure(figsize = 1.25 * np.array([6.4,4.6]), dpi = 200)
22for colname
in colnames[1:]:
23 plt.plot( df[colnames[0]].values
28plt.xticks(fontsize = fontsize)
29plt.yticks(fontsize = fontsize)
30ax.set_xlabel(colnames[0], fontsize = fontsize)
31ax.set_ylabel(
"Runtime [ seconds ]", fontsize = fontsize)
32ax.set_title(
" vs. ".join(colnames[1:])+
"\nLower is better.", fontsize = fontsize)
36plt.grid(visible =
True, which =
"both", axis =
"both", color =
"0.85", linestyle =
"-")
37ax.tick_params(axis =
"y", which =
"minor")
38ax.tick_params(axis =
"x", which =
"minor")
39ax.legend ( colnames[1:]
46plt.savefig(
"benchmark." + dirname +
".runtime.png")
52ax = plt.figure(figsize = 1.25 * np.array([6.4,4.6]), dpi = 200)
55plt.plot( df[colnames[0]].values
56 , np.ones(len(df[colnames[0]].values))
61for colname
in colnames[2:]:
62 plt.plot( df[colnames[0]].values
63 , df[colname].values / df[colnames[1]].values
67plt.xticks(fontsize = fontsize)
68plt.yticks(fontsize = fontsize)
69ax.set_xlabel(colnames[0], fontsize = fontsize)
70ax.set_ylabel(
"Runtime compared to {}".format(colnames[1]), fontsize = fontsize)
71ax.set_title(
"Runtime Ratio Comparison. Lower means faster.\nLower than 1 means faster than {}.".format(colnames[1]), fontsize = fontsize)
75plt.grid(visible =
True, which =
"both", axis =
"both", color =
"0.85", linestyle =
"-")
76ax.tick_params(axis =
"y", which =
"minor")
77ax.tick_params(axis =
"x", which =
"minor")
78ax.legend ( colnames[1:]
85plt.savefig(
"benchmark." + dirname +
".runtime.ratio.png")
91df = pd.read_csv(
"root.count", delimiter =
",")
92colnames = list(df.columns.values)
94ax = plt.figure(figsize = 1.25 * np.array([6.4,4.6]), dpi = 200)
97plt.plot( range(1, df[colnames[0]].values[-1] * 2)
98 , range(1, df[colnames[0]].values[-1] * 2)
103for colname
in colnames[1:-1]:
104 plt.plot( df[colnames[0]].values
109plt.xticks(fontsize = fontsize)
110plt.yticks(fontsize = fontsize)
111ax.set_xlabel(colnames[0], fontsize = fontsize)
112ax.set_ylabel(
"Number of Roots Found", fontsize = fontsize)
113ax.set_title(
" vs. ".join(colnames[1:])+
"\nHigher is better.", fontsize = fontsize)
117plt.grid(visible =
True, which =
"both", axis =
"both", color =
"0.85", linestyle =
"-")
118ax.tick_params(axis =
"y", which =
"minor")
119ax.tick_params(axis =
"x", which =
"minor")
120ax.legend ( [
"Equality"] + list(colnames[1:])
123 , fontsize = fontsize
127plt.savefig(
"benchmark." + dirname +
".root.count.png")
Visualization of the benchmark output ⛓
Benchmark moral ⛓
- Among all summation algorithms, jenkins_type appears to offer the fastest root finding algorithms.
- The eigen_type method also tends to offer excellent performance.
- Unlike the above two, laguerre_type algorithm tends to significantly trail behind both in performance and reliability in finding all roots of the polynomial.
- Test:
- test_pm_polynomial
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.
-
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.
-
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.
- Copyright
- Computational Data Science Lab
- Author:
- Amir Shahmoradi, Oct 16, 2009, 11:14 AM, Michigan
Definition at line 3750 of file pm_polynomial.F90.