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

Generate and return the roots of a polynomial of arbitrary degree specified by its coefficients coef.
More...

Detailed Description

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,
  1. type complex of kind any supported by the processor (e.g., CK, CK32, CK64, or CK128), or
  2. 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:
  1. the scalar constant eigen or a constant object of type eigen_type implying the use of the Eigenvalue root-finding method.
  2. the scalar constant jenkins or a constant object of type jenkins_type implying the use of the Jenkins-Traub root-finding method.
  3. the scalar constant laguerre or a constant object of type laguerre_type implying the use of the Laguerre root-finding method.
  4. 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?
  1. If you have a polynomial of highly varying coefficients, then the Eigenvalue method as specified by eigen_type is likely more reliable.
  2. 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(:)) ! allocatable output.
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

1program example
2
3 use pm_kind, only: SK, IK, LK
4 use pm_kind, only: CKG => CKS ! all processor real and complex kinds are supported.
5 use pm_io, only: display_type
8
9 implicit none
10
11 complex(CKG), allocatable :: coef(:), root(:)
12
13 type(display_type) :: disp
14 disp = display_type(file = "main.out.F90")
15
16 call disp%skip()
17 call disp%show("coef = [complex(CKG) :: -4, 1, -4, 1]")
18 coef = [complex(CKG) :: -4, 1, -4, 1]
19 call disp%show("coef")
20 call disp%show( coef )
21 call disp%show("root = getPolyRoot(coef%re) ! polynomial with real coefficients.")
22 root = getPolyRoot(coef%re)
23 call disp%show("root")
24 call disp%show( root )
25 call disp%show("getPolyVal(coef, root)")
26 call disp%show( getPolyVal(coef, root) )
27 call disp%skip()
28 call disp%show("root = getPolyRoot(coef%re, method = eigen) ! polynomial with real coefficients.")
29 root = getPolyRoot(coef%re, method = eigen)
30 call disp%show("root")
31 call disp%show( root )
32 call disp%show("getPolyVal(coef, root)")
33 call disp%show( getPolyVal(coef, root) )
34 call disp%skip()
35 call disp%show("root = getPolyRoot(coef%re, method = jenkins) ! polynomial with real coefficients.")
36 root = getPolyRoot(coef%re, method = jenkins)
37 call disp%show("root")
38 call disp%show( root )
39 call disp%show("getPolyVal(coef, root)")
40 call disp%show( getPolyVal(coef, root) )
41 call disp%skip()
42 call disp%show("root = getPolyRoot(coef%re, method = laguerre) ! polynomial with real coefficients.")
43 root = getPolyRoot(coef%re, method = laguerre)
44 call disp%show("root")
45 call disp%show( root )
46 call disp%show("getPolyVal(coef, root)")
47 call disp%show( getPolyVal(coef, root) )
48 call disp%skip()
49
50 call disp%skip()
51 call disp%show("coef = [complex(CKG) :: -120, 274, -225, 85, -15, 1]")
52 coef = [complex(CKG) :: -120, 274, -225, 85, -15, 1]
53 call disp%show("coef")
54 call disp%show( coef )
55 call disp%show("root = getPolyRoot(coef%re) ! polynomial with real coefficients.")
56 root = getPolyRoot(coef%re)
57 call disp%show("root")
58 call disp%show( root )
59 call disp%show("getPolyVal(coef, root)")
60 call disp%show( getPolyVal(coef, root) )
61 call disp%skip()
62 call disp%show("root = getPolyRoot(coef%re, method = eigen) ! polynomial with real coefficients.")
63 root = getPolyRoot(coef%re, method = eigen)
64 call disp%show("root")
65 call disp%show( root )
66 call disp%show("getPolyVal(coef, root)")
67 call disp%show( getPolyVal(coef, root) )
68 call disp%skip()
69 call disp%show("root = getPolyRoot(coef%re, method = jenkins) ! polynomial with real coefficients.")
70 root = getPolyRoot(coef%re, method = jenkins)
71 call disp%show("root")
72 call disp%show( root )
73 call disp%show("getPolyVal(coef, root)")
74 call disp%show( getPolyVal(coef, root) )
75 call disp%skip()
76 call disp%show("root = getPolyRoot(coef%re, method = laguerre) ! polynomial with real coefficients.")
77 root = getPolyRoot(coef%re, method = laguerre)
78 call disp%show("root")
79 call disp%show( root )
80 call disp%show("getPolyVal(coef, root)")
81 call disp%show( getPolyVal(coef, root) )
82 call disp%skip()
83
84 call disp%skip()
85 call disp%show("coef = [complex(CKG) :: 8, -8, 16, -16, 8, -8]")
86 coef = [complex(CKG) :: 8, -8, 16, -16, 8, -8]
87 call disp%show("coef")
88 call disp%show( coef )
89 call disp%show("root = getPolyRoot(coef%re) ! polynomial with real coefficients.")
90 root = getPolyRoot(coef%re)
91 call disp%show("root")
92 call disp%show( root )
93 call disp%show("getPolyVal(coef, root)")
94 call disp%show( getPolyVal(coef, root) )
95 call disp%skip()
96 call disp%show("root = getPolyRoot(coef%re, method = eigen) ! polynomial with real coefficients.")
97 root = getPolyRoot(coef%re, method = eigen)
98 call disp%show("root")
99 call disp%show( root )
100 call disp%show("getPolyVal(coef, root)")
101 call disp%show( getPolyVal(coef, root) )
102 call disp%skip()
103 call disp%show("root = getPolyRoot(coef%re, method = jenkins) ! polynomial with real coefficients.")
104 root = getPolyRoot(coef%re, method = jenkins)
105 call disp%show("root")
106 call disp%show( root )
107 call disp%show("getPolyVal(coef, root)")
108 call disp%show( getPolyVal(coef, root) )
109 call disp%skip()
110 call disp%show("root = getPolyRoot(coef%re, method = laguerre) ! polynomial with real coefficients.")
111 root = getPolyRoot(coef%re, method = laguerre)
112 call disp%show("root")
113 call disp%show( root )
114 call disp%show("getPolyVal(coef, root)")
115 call disp%show( getPolyVal(coef, root) )
116 call disp%skip()
117
118 call disp%skip()
119 call disp%show("coef = [complex(CKG) :: -120, 274, -225, 85, -15, 1]")
120 coef = [complex(CKG) :: -120, 274, -225, 85, -15, 1]
121 call disp%show("coef")
122 call disp%show( coef )
123 call disp%show("root = getPolyRoot(coef)")
124 root = getPolyRoot(coef)
125 call disp%show("root")
126 call disp%show( root )
127 call disp%show("getPolyVal(coef, root)")
128 call disp%show( getPolyVal(coef, root) )
129 call disp%skip()
130 call disp%show("root = getPolyRoot(coef, method = eigen)")
131 root = getPolyRoot(coef, method = eigen)
132 call disp%show("root")
133 call disp%show( root )
134 call disp%show("getPolyVal(coef, root)")
135 call disp%show( getPolyVal(coef, root) )
136 call disp%skip()
137 call disp%show("root = getPolyRoot(coef, method = jenkins)")
138 root = getPolyRoot(coef, method = jenkins)
139 call disp%show("root")
140 call disp%show( root )
141 call disp%show("getPolyVal(coef, root)")
142 call disp%show( getPolyVal(coef, root) )
143 call disp%skip()
144 call disp%show("root = getPolyRoot(coef, method = laguerre)")
145 root = getPolyRoot(coef, method = laguerre)
146 call disp%show("root")
147 call disp%show( root )
148 call disp%show("getPolyVal(coef, root)")
149 call disp%show( getPolyVal(coef, root) )
150 call disp%skip()
151
152end 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
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...
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 CKS
The single-precision complex kind in Fortran mode. On most platforms, this is a 32-bit real kind.
Definition: pm_kind.F90:570
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
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
2coef = [complex(CKG) :: -4, 1, -4, 1]
3coef
4(-4.00000000, +0.00000000), (+1.00000000, +0.00000000), (-4.00000000, +0.00000000), (+1.00000000, +0.00000000)
5root = getPolyRoot(coef%re) ! polynomial with real coefficients.
6root
7(+3.99999952, +0.00000000), (+0.00000000, +0.999999821), (+0.00000000, -0.999999821)
8getPolyVal(coef, root)
9(-0.810623169E-5, +0.00000000), (-0.143051147E-5, +0.357627812E-6), (-0.143051147E-5, -0.357627812E-6)
10
11root = getPolyRoot(coef%re, method = eigen) ! polynomial with real coefficients.
12root
13(+3.99999952, +0.00000000), (+0.00000000, +0.999999821), (+0.00000000, -0.999999821)
14getPolyVal(coef, root)
15(-0.810623169E-5, +0.00000000), (-0.143051147E-5, +0.357627812E-6), (-0.143051147E-5, -0.357627812E-6)
16
17root = getPolyRoot(coef%re, method = jenkins) ! polynomial with real coefficients.
18root
19(-0.894069672E-6, +1.00000000), (-0.894069672E-6, -1.00000000), (+4.00000191, +0.00000000)
20getPolyVal(coef, root)
21(+0.190734863E-5, +0.715255919E-5), (+0.190734863E-5, -0.715255919E-5), (+0.324249268E-4, +0.00000000)
22
23root = getPolyRoot(coef%re, method = laguerre) ! polynomial with real coefficients.
24root
25(+0.977888703E-8, +1.00000000), (-0.467116479E-8, -1.00000000), (+0.160798663E-8, -1.00000000)
26getPolyVal(coef, root)
27(+0.00000000, -0.391155481E-7), (+0.00000000, -0.186846592E-7), (+0.00000000, +0.643194653E-8)
28
29
30coef = [complex(CKG) :: -120, 274, -225, 85, -15, 1]
31coef
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)
33root = getPolyRoot(coef%re) ! polynomial with real coefficients.
34root
35(+4.99990797, +0.00000000), (+4.00018120, +0.00000000), (+2.99988103, +0.00000000), (+2.00003195, +0.00000000), (+0.999997020, +0.00000000)
36getPolyVal(coef, root)
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)
38
39root = getPolyRoot(coef%re, method = eigen) ! polynomial with real coefficients.
40root
41(+4.99990797, +0.00000000), (+4.00018120, +0.00000000), (+2.99988103, +0.00000000), (+2.00003195, +0.00000000), (+0.999997020, +0.00000000)
42getPolyVal(coef, root)
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)
44
45root = getPolyRoot(coef%re, method = jenkins) ! polynomial with real coefficients.
46root
47(+0.999999523, +0.00000000), (+2.00003052, +0.00000000), (+2.99990487, +0.00000000), (+4.00010586, +0.00000000), (+4.99995899, +0.00000000)
48getPolyVal(coef, root)
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)
50
51root = getPolyRoot(coef%re, method = laguerre) ! polynomial with real coefficients.
52root
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)
54getPolyVal(coef, root)
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)
56
57
58coef = [complex(CKG) :: 8, -8, 16, -16, 8, -8]
59coef
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)
61root = getPolyRoot(coef%re) ! polynomial with real coefficients.
62root
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)
64getPolyVal(coef, root)
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)
66
67root = getPolyRoot(coef%re, method = eigen) ! polynomial with real coefficients.
68root
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)
70getPolyVal(coef, root)
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)
72
73root = getPolyRoot(coef%re, method = jenkins) ! polynomial with real coefficients.
74root
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)
76getPolyVal(coef, root)
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)
78
79root = getPolyRoot(coef%re, method = laguerre) ! polynomial with real coefficients.
80root
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)
82getPolyVal(coef, root)
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)
84
85
86coef = [complex(CKG) :: -120, 274, -225, 85, -15, 1]
87coef
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)
89root = getPolyRoot(coef)
90root
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)
92getPolyVal(coef, root)
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)
94
95root = getPolyRoot(coef, method = eigen)
96root
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)
98getPolyVal(coef, root)
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)
100
101root = getPolyRoot(coef, method = jenkins)
102root
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)
104getPolyVal(coef, root)
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)
106
107root = getPolyRoot(coef, method = laguerre)
108root
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)
110getPolyVal(coef, root)
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)
112
113
Benchmarks:


Benchmark :: The effects of method on runtime efficiency

The following program compares the runtime performance of setPolyRoot using different polynomial root finding algorithms.
1! Test the performance of `setPolyRoot()` with and without the selection `control` argument.
2program benchmark
3
4 use pm_bench, only: bench_type
5 use pm_distUnif, only: setUnifRand
7 use pm_kind, only: SK, IK, LK, RKH, RK, RKG => RKD
8 use iso_fortran_env, only: error_unit
9
10 implicit none
11
12 integer(IK) , parameter :: degmax = 9_IK
13 integer(IK) , allocatable :: rootCount(:)
14 integer(IK) :: ibench
15 integer(IK) :: ideg
16 integer(IK) :: degree
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(:)
23
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) &
28 ]
29 allocate(rootCount(size(bench)))
30
31 write(*,"(*(g0,:,' '))")
32 write(*,"(*(g0,:,' '))") "Benchmarking setPolyRoot()"
33 write(*,"(*(g0,:,' '))")
34
35 open(newunit = fileUnit, file = "main.out", status = "replace")
36 open(newunit = rootUnit, file = "root.count", status = "replace")
37
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
41
42 degree = 2**ideg
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)
47 end do
48 write(rootUnit,"(*(g0,:,','))") degree, rootCount
49 write(fileUnit,"(*(g0,:,','))") degree, (bench(ibench)%timing%mean, ibench = 1, size(bench))
50
51 end do loopOverArraySize
52 write(*,"(*(g0,:,' '))") dumsum
53 write(*,"(*(g0,:,' '))")
54
55 close(fileUnit)
56
57contains
58
59 !%%%%%%%%%%%%%%%%%%%%
60 ! procedure wrappers.
61 !%%%%%%%%%%%%%%%%%%%%
62
63 subroutine setOverhead()
64 call getDummy()
65 end subroutine
66
67 subroutine getDummy()
68 dumsum = dumsum + rootCount(ibench)
69 end subroutine
70
71 subroutine setPolyRootEigen()
72 use pm_polynomial, only: setPolyRoot, method => eigen
73 call setPolyRoot(root(1 : degree), rootCount(ibench), coef(0 : degree), method)
74 call getDummy()
75 end subroutine
76
77 subroutine setPolyRootJenkins()
78 use pm_polynomial, only: setPolyRoot, method => jenkins
79 call setPolyRoot(root(1 : degree), rootCount(ibench), coef(0 : degree), method)
80 call getDummy()
81 end subroutine
82
83 subroutine setPolyRootLaguerre()
84 use pm_polynomial, only: setPolyRoot, method => laguerre
85 call setPolyRoot(root(1 : degree), rootCount(ibench), coef(0 : degree), method)
86 call getDummy()
87 end subroutine
88
89 subroutine setPolyRootSGL()
90 use pm_polynomial, only: setPolyRoot, method => sgl
91 call setPolyRoot(root(1 : degree), rootCount(ibench), coef(0 : degree), method)
92 rootCount(ibench) = root(1)%re
93 call getDummy()
94 end subroutine
95
96end program benchmark
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...
Definition: pm_bench.F90:574
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...
Definition: pm_bench.F90:41
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...
Definition: pm_kind.F90:543
integer, parameter RKD
The double precision real kind in Fortran mode. On most platforms, this is an 64-bit real kind.
Definition: pm_kind.F90:568
integer, parameter RKH
The scalar integer constant of intrinsic default kind, representing the highest-precision real kind t...
Definition: pm_kind.F90:858
This is the class for creating benchmark and performance-profiling objects.
Definition: pm_bench.F90:386


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


Postprocessing of the benchmark output

1#!/usr/bin/env python
2
3import matplotlib.pyplot as plt
4import pandas as pd
5import numpy as np
6
7import os
8dirname = os.path.basename(os.getcwd())
9
10fontsize = 14
11
12df = pd.read_csv("main.out", delimiter = ",")
13colnames = list(df.columns.values)
14
15
18
19ax = plt.figure(figsize = 1.25 * np.array([6.4,4.6]), dpi = 200)
20ax = plt.subplot()
21
22for colname in colnames[1:]:
23 plt.plot( df[colnames[0]].values
24 , df[colname].values
25 , linewidth = 2
26 )
27
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)
33ax.set_xscale("log")
34ax.set_yscale("log")
35plt.minorticks_on()
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:]
40 #, loc='center left'
41 #, bbox_to_anchor=(1, 0.5)
42 , fontsize = fontsize
43 )
44
45plt.tight_layout()
46plt.savefig("benchmark." + dirname + ".runtime.png")
47
48
51
52ax = plt.figure(figsize = 1.25 * np.array([6.4,4.6]), dpi = 200)
53ax = plt.subplot()
54
55plt.plot( df[colnames[0]].values
56 , np.ones(len(df[colnames[0]].values))
57 , linestyle = "--"
58 #, color = "black"
59 , linewidth = 2
60 )
61for colname in colnames[2:]:
62 plt.plot( df[colnames[0]].values
63 , df[colname].values / df[colnames[1]].values
64 , linewidth = 2
65 )
66
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)
72ax.set_xscale("log")
73ax.set_yscale("log")
74plt.minorticks_on()
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:]
79 #, bbox_to_anchor = (1, 0.5)
80 #, loc = "center left"
81 , fontsize = fontsize
82 )
83
84plt.tight_layout()
85plt.savefig("benchmark." + dirname + ".runtime.ratio.png")
86
87
90
91df = pd.read_csv("root.count", delimiter = ",")
92colnames = list(df.columns.values)
93
94ax = plt.figure(figsize = 1.25 * np.array([6.4,4.6]), dpi = 200)
95ax = plt.subplot()
96
97plt.plot( range(1, df[colnames[0]].values[-1] * 2)
98 , range(1, df[colnames[0]].values[-1] * 2)
99 , linestyle = "--"
100 , color = "black"
101 , linewidth = 2
102 )
103for colname in colnames[1:-1]:
104 plt.plot( df[colnames[0]].values
105 , df[colname].values
106 , linewidth = 2
107 )
108
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)
114ax.set_xscale("log")
115ax.set_yscale("log")
116plt.minorticks_on()
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:])
121 #, loc='center left'
122 #, bbox_to_anchor=(1, 0.5)
123 , fontsize = fontsize
124 )
125
126plt.tight_layout()
127plt.savefig("benchmark." + dirname + ".root.count.png")


Visualization of the benchmark output


Benchmark moral

  1. Among all summation algorithms, jenkins_type appears to offer the fastest root finding algorithms.
  2. The eigen_type method also tends to offer excellent performance.
  3. 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.

  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 3750 of file pm_polynomial.F90.


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