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

Return the roots of a polynomial of arbitrary degree specified by its coefficients coef.
More...

Detailed Description

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 method.

Parameters
[out]root: The output contiguous vector of type complex of the same kind as the input argument coef and the same size as the degree of the polynomial (i.e., size(coef) - 1), containing the roots of the polynomial determined from the polynomial coefficients coef.
[out]count: The output scalar of type integer of default kind IK, containing the total number of roots computed and stored in the output vector slice root(1 : count).
A value of count = size(root) implies a successful computation of all polynomial roots.
A value of count = 0 implies a total failure of the algorithm in finding any roots.
The condition count < size(root) occurs if either,
  • the algorithm fails to converge, or
  • the algorithm fails to identify all roots of the polynomial, or
  • the condition coef(size(coef)) == 0. occurs (i.e., when the coefficient of the highest power of the polynomial is zero), or
  • the condition size(coef) < 2 occurs (i.e., when the input polynomial is a constant).
[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.
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.


Possible calling interfaces

call setPolyRoot(root(1 : degree), count, coef(0 : degree), method) ! `degree` is the degree of the polynomial.
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...
Warning
The condition 1 < size(coef) must hold for the corresponding input arguments.
The condition coef(size(coef)) /= 0. must hold for the corresponding input arguments (i.e., the coefficient of the highest-degree term must be non-zero).
The condition all(shape(workspace) == size(coef) - 1) must hold for the corresponding input arguments.
The condition size(root) == size(coef) - 1 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.
As such polynomials of high degree can be problematic for root-finding algorithms.
Remarks
The procedures under discussion are impure.
This generic interface combines, significantly extends, and modernizes the SPOLZ, DPOLZ, CPOLZ, and ZPOLZ routines of the MATH77 netlib public-domain library.
See also
getRoot
setRoot
getPolyRoot
setPolyRoot


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK
4 use pm_kind, only: RKS, RKD, RKH ! all processor real/complex kinds are supported.
5 use pm_io, only: display_type
6 use pm_polynomial, only: setPolyRoot, cmplx_roots_gen
10
11 implicit none
12
13 integer(IK) :: count
14
15 type(display_type) :: disp
16 disp = display_type(file = "main.out.F90")
17
18! Template to avoid code duplication in the
19! example for various types and kind parameters.
20! See the output below for the actual code.
21! developer guideline:
22! TYPE - type of coefficient vector: `real`, `complex`
23! RKG - kind of coefficient vector: any supported by the processor
24! CREP - Number of coefficient elements to display per line: use `1` for `real` and `2` for `complex` coefficients.
25#define GET_ROOT_EIG(CREP, TYPE, RKG) \
26block; \
27use pm_val2str, only: getStr; \
28TYPE(RKG) , allocatable :: coef(:); \
29complex(RKG), allocatable :: root(:); \
30call disp%skip(); \
31coef = COEF; \
32call disp%show("coef"); \
33call disp%show( coef , format = "(sp,"//getStr(CREP)//"(g0,:,', '))"); \
34call disp%show("call setResized(root, size(coef, 1, IK) - 1_IK)"); \
35 call setResized(root, size(coef, 1, IK) - 1_IK); \
36call disp%show("call setPolyRoot(root, count, coef, eigen)"); \
37 call setPolyRoot(root, count, coef, eigen); \
38call disp%show("count"); \
39call disp%show( count ); \
40call disp%show("root(1:count)"); \
41call disp%show( root(1:count) , format = "(sp,2(g0,:,', '))"); \
42call disp%show("getPolyVal(coef, root(1:count))"); \
43call disp%show( getPolyVal(coef, root(1:count)) , format = "(sp,2(g0,:,', '))"); \
44call disp%skip(); \
45end block;
46
47! Template to avoid code duplication in the
48! example for various types and kind parameters.
49! See the output below for the actual code.
50! developer guideline:
51! TYPE - type of coefficient vector: `real`, `complex`
52! RKG - kind of coefficient vector: any supported by the processor
53! CREP - Number of coefficient elements to display per line: use `1` for `real` and `2` for `complex` coefficients.
54#define GET_ROOT_JEN(CREP, TYPE, RKG) \
55block; \
56use pm_val2str, only: getStr; \
57TYPE(RKG) , allocatable :: coef(:); \
58complex(RKG), allocatable :: root(:); \
59call disp%skip(); \
60coef = COEF; \
61call disp%show("coef"); \
62call disp%show( coef , format = "(sp,"//getStr(CREP)//"(g0,:,', '))"); \
63call disp%show("call setResized(root, size(coef, 1, IK) - 1_IK)"); \
64 call setResized(root, size(coef, 1, IK) - 1_IK); \
65call disp%show("call setPolyRoot(root, count, coef, jenkins)"); \
66 call setPolyRoot(root, count, coef, jenkins); \
67call disp%show("count"); \
68call disp%show( count ); \
69call disp%show("root(1:count)"); \
70call disp%show( root(1:count) , format = "(sp,2(g0,:,', '))"); \
71call disp%show("getPolyVal(coef, root(1:count))"); \
72call disp%show( getPolyVal(coef, root(1:count)) , format = "(sp,2(g0,:,', '))"); \
73call disp%skip(); \
74end block;
75
76! Template to avoid code duplication in the
77! example for various types and kind parameters.
78! See the output below for the actual code.
79! developer guideline:
80! TYPE - type of coefficient vector: `real`, `complex`
81! RKG - kind of coefficient vector: any supported by the processor
82! CREP - Number of coefficient elements to display per line: use `1` for `real` and `2` for `complex` coefficients.
83#define GET_ROOT_LAG(CREP, TYPE, RKG) \
84block; \
85use pm_val2str, only: getStr; \
86TYPE(RKG) , allocatable :: coef(:); \
87complex(RKG), allocatable :: root(:); \
88call disp%skip(); \
89coef = COEF; \
90call disp%show("coef"); \
91call disp%show( coef , format = "(sp,"//getStr(CREP)//"(g0,:,', '))"); \
92call disp%show("call setResized(root, size(coef, 1, IK) - 1_IK)"); \
93 call setResized(root, size(coef, 1, IK) - 1_IK); \
94call disp%show("call setPolyRoot(root, count, coef, laguerre)"); \
95 call setPolyRoot(root, count, coef, laguerre); \
96call disp%show("count"); \
97call disp%show( count ); \
98call disp%show("root(1:count)"); \
99call disp%show( root(1:count) , format = "(sp,2(g0,:,', '))"); \
100call disp%show("getPolyVal(coef, root(1:count))"); \
101call disp%show( getPolyVal(coef, root(1:count)) , format = "(sp,2(g0,:,', '))"); \
102call disp%skip(); \
103end block;
104
105! Template to avoid code duplication in the
106! example for various types and kind parameters.
107! See the output below for the actual code.
108! developer guideline:
109! TYPE - type of coefficient vector: `real`, `complex`
110! RKG - kind of coefficient vector: any supported by the processor
111! CREP - Number of coefficient elements to display per line: use `1` for `real` and `2` for `complex` coefficients.
112#define GET_ROOT_SKO(CREP, TYPE, RKG) \
113block; \
114use pm_val2str, only: getStr; \
115TYPE(RKG) , allocatable :: coef(:); \
116complex(RKG), allocatable :: root(:); \
117call disp%skip(); \
118coef = COEF; \
119call disp%show("coef"); \
120call disp%show( coef , format = "(sp,"//getStr(CREP)//"(g0,:,', '))"); \
121call disp%show("call setResized(root, size(coef, 1, IK) - 1_IK)"); \
122 call setResized(root, size(coef, 1, IK) - 1_IK); \
123call disp%show("call cmplx_roots_gen(root, coef, size(coef) - 1, .true., .true.)"); \
124 call cmplx_roots_gen(root, coef, size(coef) - 1, .true., .true.); \
125call disp%show("root(1:count)"); \
126call disp%show( root(1:count) , format = "(sp,2(g0,:,', '))"); \
127call disp%show("getPolyVal(coef, root(1:count))"); \
128call disp%show( getPolyVal(coef, root(1:count)) , format = "(sp,2(g0,:,', '))"); \
129call disp%skip(); \
130end block;
131
132 call disp%skip()
133 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
134 call disp%show("! Compute the roots of polynomials with real coefficients.")
135 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
136 call disp%skip()
137
138#define COEF \
139 [8, -8, 16, -16, 8, -8]
140 GET_ROOT_EIG(1, real, RKS)
141 GET_ROOT_EIG(1, real, RKD)
142 GET_ROOT_EIG(1, real, RKH)
143 GET_ROOT_JEN(1, real, RKS)
144 GET_ROOT_JEN(1, real, RKD)
145 GET_ROOT_JEN(1, real, RKH)
146 GET_ROOT_LAG(1, real, RKS)
147 GET_ROOT_LAG(1, real, RKD)
148 GET_ROOT_LAG(1, real, RKH)
149
150 call disp%skip()
151 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
152 call disp%show("! Compute the roots of polynomials with real coefficients.")
153 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
154 call disp%skip()
155
156#define COEF \
157 [3628800, -10628640, 12753576, -8409500, 3416930, -902055, 157773, -18150, 1320, -55, 1]
158 GET_ROOT_EIG(1, real, RKS)
159 GET_ROOT_EIG(1, real, RKD)
160 GET_ROOT_EIG(1, real, RKH)
161 GET_ROOT_JEN(1, real, RKS)
162 GET_ROOT_JEN(1, real, RKD)
163 GET_ROOT_JEN(1, real, RKH)
164 GET_ROOT_LAG(1, real, RKS)
165 GET_ROOT_LAG(1, real, RKD)
166 GET_ROOT_LAG(1, real, RKH)
167
168 call disp%skip()
169 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
170 call disp%show("! Compute the roots of polynomials with complex coefficients with zeros 1,2,...,10.")
171 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
172 call disp%skip()
173
174#define COEF \
175 [3628800, -10628640, 12753576, -8409500, 3416930, -902055, 157773, -18150, 1320, -55, 1]
176 GET_ROOT_EIG(2, complex, RKS)
177 GET_ROOT_EIG(2, complex, RKD)
178 GET_ROOT_EIG(2, complex, RKH)
179 GET_ROOT_JEN(2, complex, RKS)
180 GET_ROOT_JEN(2, complex, RKD)
181 GET_ROOT_JEN(2, complex, RKH)
182 GET_ROOT_LAG(2, complex, RKS)
183 GET_ROOT_LAG(2, complex, RKD)
184 GET_ROOT_LAG(2, complex, RKH)
185
186 call disp%skip()
187 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
188 call disp%show("! Compute the roots of polynomials with complex coefficients with zeros on imaginary axis.")
189 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
190 call disp%skip()
191
192#define COEF \
193 [(0._RKH, 1._RKH), (-10001.0001_RKH, 0._RKH), (0._RKH, -10001.0001_RKH), (1._RKH, 0._RKH)]
194 GET_ROOT_EIG(2, complex, RKD)
195 GET_ROOT_EIG(2, complex, RKH)
196 GET_ROOT_JEN(2, complex, RKD)
197 GET_ROOT_JEN(2, complex, RKH)
198 GET_ROOT_LAG(2, complex, RKD)
199 GET_ROOT_LAG(2, complex, RKH)
200
201 call disp%skip()
202 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
203 call disp%show("! Compute the roots of polynomials with complex coefficients with zeros at 1+i,1/2*(1+i)....1/(2**-9)*(1+i).")
204 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
205 call disp%skip()
206
207#define COEF \
208[ (0._RKH, 9.094947017729282e-13_RKH) \
209, (-4.652065399568528e-10_RKH, -4.652065399568528e-10_RKH) \
210, (1.584803612786345e-7_RKH, 0._RKH) \
211, (-1.154642632172909e-5_RKH, 1.154642632172909e-5_RKH) \
212, (0._RKH, -7.820779428584501e-4_RKH) \
213, (1.271507365163416e-2_RKH, 1.271507365163416e-2_RKH) \
214, (-.2002119533717632e0_RKH, 0._RKH) \
215, (.7567065954208374e0_RKH, -7.567065954208374e-1_RKH) \
216, (0._RKH, 2.658859252929688e0_RKH) \
217, (-1.998046875_RKH, -1.998046875_RKH) \
218, (1._RKH, 0._RKH) \
219]
220GET_ROOT_EIG(2, complex, RKD)
221GET_ROOT_EIG(2, complex, RKH)
222GET_ROOT_JEN(2, complex, RKD)
223GET_ROOT_JEN(2, complex, RKH)
224GET_ROOT_LAG(2, complex, RKD)
225GET_ROOT_LAG(2, complex, RKH)
226GET_ROOT_SKO(2, complex, RKD)
227
228 call disp%skip()
229 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
230 call disp%show("! Compute the roots of polynomials with complex coefficients with multiple zeros.")
231 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
232 call disp%skip()
233
234#define COEF \
235[ (288._RKH, 0._RKH) \
236, (-1344._RKH, 504._RKH) \
237, (2204._RKH, -2352._RKH) \
238, (-920._RKH, 4334._RKH) \
239, (-1587._RKH, -3836._RKH) \
240, (2374._RKH, 1394._RKH) \
241, (-1293._RKH, 200._RKH) \
242, (284._RKH, -334._RKH) \
243, (3._RKH, 100._RKH) \
244, (-10._RKH, -10._RKH) \
245, (1._RKH, 0._RKH) \
246]
247GET_ROOT_EIG(2, complex, RKS)
248GET_ROOT_EIG(2, complex, RKD)
249GET_ROOT_EIG(2, complex, RKH)
250GET_ROOT_JEN(2, complex, RKS)
251GET_ROOT_JEN(2, complex, RKD)
252GET_ROOT_JEN(2, complex, RKH)
253GET_ROOT_LAG(2, complex, RKS)
254GET_ROOT_LAG(2, complex, RKD)
255GET_ROOT_LAG(2, complex, RKH)
256GET_ROOT_SKO(2, complex, RKD)
257
258 call disp%skip()
259 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
260 call disp%show("! Compute the roots of polynomials with complex coefficients with 12 zeros evenly distribute on a circle of radius 1. centered at 0+2i.")
261 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
262 call disp%skip()
263
264#define COEF \
265[ (4095._RKH, 0._RKH) \
266, (0._RKH, 24576._RKH) \
267, (-67584._RKH, 0._RKH) \
268, (0._RKH, -112640._RKH) \
269, (126720._RKH, 0._RKH) \
270, (0._RKH, 101376._RKH) \
271, (-59136._RKH, 0._RKH) \
272, (0._RKH, -25344._RKH) \
273, (7920._RKH, 0._RKH) \
274, (0._RKH, 1760._RKH) \
275, (-264._RKH, 0._RKH) \
276, (0._RKH, -24._RKH) \
277, (1._RKH, 0._RKH) \
278]
279GET_ROOT_EIG(2, complex, RKS)
280GET_ROOT_EIG(2, complex, RKD)
281GET_ROOT_EIG(2, complex, RKH)
282GET_ROOT_JEN(2, complex, RKS)
283GET_ROOT_JEN(2, complex, RKD)
284GET_ROOT_JEN(2, complex, RKH)
285GET_ROOT_LAG(2, complex, RKS)
286GET_ROOT_LAG(2, complex, RKD)
287GET_ROOT_LAG(2, complex, RKH)
288GET_ROOT_SKO(2, complex, RKD)
289
290end program example
Allocate or resize (shrink or expand) an input allocatable scalar string or array of rank 1....
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 ...
Generate and return the conversion of the input value to an output Fortran string,...
Definition: pm_val2str.F90:167
This module contains procedures and generic interfaces for resizing allocatable arrays of various typ...
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 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 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
integer, parameter RKS
The single-precision real kind in Fortran mode. On most platforms, this is an 32-bit real kind.
Definition: pm_kind.F90:567
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...
This module contains the generic procedures for converting values of different types and kinds to For...
Definition: pm_val2str.F90:58
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
2!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
3! Compute the roots of polynomials with real coefficients.
4!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
5
6
7coef
8+8.00000000,
9-8.00000000,
10+16.0000000,
11-16.0000000,
12+8.00000000,
13-8.00000000
14call setResized(root, size(coef, 1, IK) - 1_IK)
15call setPolyRoot(root, count, coef, eigen)
16count
17+5
18root(1:count)
19+0.999999642, +0.00000000,
20-0.429153442E-4, +1.00021696,
21-0.429153442E-4, -1.00021696,
22+0.427067280E-4, +0.999782741,
23+0.427067280E-4, -0.999782741
24getPolyVal(coef, root(1:count))
25+0.114440918E-4, +0.00000000,
26+0.238418579E-5, -0.804837327E-6,
27+0.238418579E-5, +0.804837327E-6,
28+0.190734863E-5, -0.566709787E-6,
29+0.190734863E-5, +0.566709787E-6
30
31
32coef
33+8.0000000000000000,
34-8.0000000000000000,
35+16.000000000000000,
36-16.000000000000000,
37+8.0000000000000000,
38-8.0000000000000000
39call setResized(root, size(coef, 1, IK) - 1_IK)
40call setPolyRoot(root, count, coef, eigen)
41count
42+5
43root(1:count)
44+0.99999999999999922, +0.0000000000000000,
45+0.95202123961968255E-9, +1.0000000121923602,
46+0.95202123961968255E-9, -1.0000000121923602,
47-0.95202196126464855E-9, +0.99999998780764010,
48-0.95202196126464855E-9, -0.99999998780764010
49getPolyVal(coef, root(1:count))
50+0.23980817331903381E-13, +0.0000000000000000,
51+0.26645352591003757E-14, -0.55147887586065774E-14,
52+0.26645352591003757E-14, +0.55147887586065774E-14,
53+0.35527136788005009E-14, -0.50706995570283210E-14,
54+0.35527136788005009E-14, +0.50706995570283210E-14
55
56
57coef
58+8.00000000000000000000000000000000000,
59-8.00000000000000000000000000000000000,
60+16.0000000000000000000000000000000000,
61-16.0000000000000000000000000000000000,
62+8.00000000000000000000000000000000000,
63-8.00000000000000000000000000000000000
64call setResized(root, size(coef, 1, IK) - 1_IK)
65call setPolyRoot(root, count, coef, eigen)
66count
67+5
68root(1:count)
69+1.00000000000000000000000000000000019, +0.00000000000000000000000000000000000,
70-0.108472542273474050653982053325797972E-16, +1.00000000000000000335355606599927526,
71-0.108472542273474050653982053325797972E-16, -1.00000000000000000335355606599927526,
72+0.108472542273474053061394483809842788E-16, +0.999999999999999996646443934000722236,
73+0.108472542273474053061394483809842788E-16, -0.999999999999999996646443934000722236
74getPolyVal(coef, root(1:count))
75-0.770371977754894341222391177033970927E-32, +0.00000000000000000000000000000000000,
76-0.308148791101957736488956470813588371E-32, +0.635981983657862231619686642857405375E-32,
77-0.308148791101957736488956470813588371E-32, -0.635981983657862231619686642857405375E-32,
78-0.308148791101957736488956470813588371E-32, +0.597463384770117558391908329517715426E-32,
79-0.308148791101957736488956470813588371E-32, -0.597463384770117558391908329517715426E-32
80
81
82coef
83+8.00000000,
84-8.00000000,
85+16.0000000,
86-16.0000000,
87+8.00000000,
88-8.00000000
89call setResized(root, size(coef, 1, IK) - 1_IK)
90call setPolyRoot(root, count, coef, jenkins)
91count
92+5
93root(1:count)
94-0.339653343E-5, +1.00000632,
95-0.339653343E-5, -1.00000632,
96+0.339158305E-5, +0.999993563,
97+0.339158305E-5, -0.999993563,
98+1.00000000, +0.00000000
99getPolyVal(coef, root(1:count))
100+0.00000000, -0.469048246E-6,
101+0.00000000, +0.469048246E-6,
102+0.00000000, +0.430129148E-6,
103+0.00000000, -0.430129148E-6,
104+0.00000000, +0.00000000
105
106
107coef
108+8.0000000000000000,
109-8.0000000000000000,
110+16.000000000000000,
111-16.000000000000000,
112+8.0000000000000000,
113-8.0000000000000000
114call setResized(root, size(coef, 1, IK) - 1_IK)
115call setPolyRoot(root, count, coef, jenkins)
116count
117+5
118root(1:count)
119+0.10909933519325752E-10, +0.99999999999780376,
120+0.10909933519325752E-10, -0.99999999999780376,
121-0.10909949390015614E-10, +1.0000000000021962,
122-0.10909949390015614E-10, -1.0000000000021962,
123+1.0000000000000000, +0.0000000000000000
124getPolyVal(coef, root(1:count))
125+0.0000000000000000, -0.48792102521272193E-16,
126+0.0000000000000000, +0.48792102521272193E-16,
127+0.0000000000000000, -0.78172649634843152E-16,
128+0.0000000000000000, +0.78172649634843152E-16,
129+0.0000000000000000, +0.0000000000000000
130
131
132coef
133+8.00000000000000000000000000000000000,
134-8.00000000000000000000000000000000000,
135+16.0000000000000000000000000000000000,
136-16.0000000000000000000000000000000000,
137+8.00000000000000000000000000000000000,
138-8.00000000000000000000000000000000000
139call setResized(root, size(coef, 1, IK) - 1_IK)
140call setPolyRoot(root, count, coef, jenkins)
141count
142+5
143root(1:count)
144+0.852303366282744871748611013224628411E-21, +1.00000000000000000000007148360124131,
145+0.852303366282744871748611013224628411E-21, -1.00000000000000000000007148360124131,
146-0.852303366282752562353390770677014811E-21, +0.999999999999999999999928516398758596,
147-0.852303366282752562353390770677014811E-21, -0.999999999999999999999928516398758596,
148+1.00000000000000000000000000000000000, +0.00000000000000000000000000000000000
149getPolyVal(coef, root(1:count))
150+0.770371977754894341222391177033970927E-33, -0.686054190319130759885022393739307596E-33,
151+0.770371977754894341222391177033970927E-33, +0.686054190319130759885022393739307596E-33,
152+0.00000000000000000000000000000000000, -0.916214605378340388758133288122927895E-33,
153+0.00000000000000000000000000000000000, +0.916214605378340388758133288122927895E-33,
154+0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000
155
156
157coef
158+8.00000000,
159-8.00000000,
160+16.0000000,
161-16.0000000,
162+8.00000000,
163-8.00000000
164call setResized(root, size(coef, 1, IK) - 1_IK)
165call setPolyRoot(root, count, coef, laguerre)
166count
167+5
168root(1:count)
169+0.111623813E-3, +0.999268234,
170-0.395235606E-3, -0.999610364,
171+0.317461498E-3, +0.999678612,
172+1.00000000, -0.395812094E-8,
173+0.999999940, +0.162981451E-8
174getPolyVal(coef, root(1:count))
175+0.219345093E-4, -0.112091075E-4,
176-0.953674316E-5, +0.959518366E-5,
177+0.667572021E-5, +0.645453110E-5,
178+0.00000000, +0.126659870E-6,
179+0.143051147E-5, -0.521540500E-7
180
181
182coef
183+8.0000000000000000,
184-8.0000000000000000,
185+16.000000000000000,
186-16.000000000000000,
187+8.0000000000000000,
188-8.0000000000000000
189call setResized(root, size(coef, 1, IK) - 1_IK)
190call setPolyRoot(root, count, coef, laguerre)
191count
192+5
193root(1:count)
194+0.32308839559944707E-8, +0.99999998512939647,
195-0.31843677481540279E-7, -0.99999997015495556,
196+0.23958124570595574E-7, +0.99999997556987996,
197+1.0000000000000000, +0.0000000000000000,
198+1.0000000000000000, +0.41359030627651384E-24
199getPolyVal(coef, root(1:count))
200+0.97699626167013776E-14, -0.22860883630403326E-14,
201-0.65725203057809267E-13, +0.57909583254923502E-13,
202+0.38191672047105385E-13, +0.36700494578024068E-13,
203+0.0000000000000000, +0.0000000000000000,
204+0.0000000000000000, -0.13234889800848443E-22
205
206
207coef
208+8.00000000000000000000000000000000000,
209-8.00000000000000000000000000000000000,
210+16.0000000000000000000000000000000000,
211-16.0000000000000000000000000000000000,
212+8.00000000000000000000000000000000000,
213-8.00000000000000000000000000000000000
214call setResized(root, size(coef, 1, IK) - 1_IK)
215call setPolyRoot(root, count, coef, laguerre)
216count
217+5
218root(1:count)
219+0.357671243096782216657689025737382332E-17, +0.999999999999999976594857690260151774,
220-0.132424594232994473131643990474026752E-16, -0.999999999999999987824683273111061743,
221+0.972702385889780612655588772995474917E-17, +0.999999999999999990399424321479437958,
222+1.00000000000000000000000000000000000, +0.398272977783113069257220099441927951E-58,
223+1.00000000000000000000000000000000000, +0.155575381946528542678601601344503106E-60
224getPolyVal(coef, root(1:count))
225+0.238815313104017245778941264880530987E-31, -0.132021106190385611906603591813480807E-31,
226-0.123259516440783094595582588325435348E-31, +0.788995639958939921101917270710366269E-32,
227+0.693334779979404907100152059330573835E-32, +0.642010863814198053577283388867887596E-32,
228+0.00000000000000000000000000000000000, -0.127447352890596182162310431821416944E-56,
229+0.00000000000000000000000000000000000, -0.497841222228891336571525124302409939E-59
230
231
232!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
233! Compute the roots of polynomials with real coefficients.
234!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
235
236
237coef
238+3628800.00,
239-10628640.0,
240+12753576.0,
241-8409500.00,
242+3416930.00,
243-902055.000,
244+157773.000,
245-18150.0000,
246+1320.00000,
247-55.0000000,
248+1.00000000
249call setResized(root, size(coef, 1, IK) - 1_IK)
250call setPolyRoot(root, count, coef, eigen)
251count
252+10
253root(1:count)
254+9.70840931, +0.304670960,
255+9.70840931, -0.304670960,
256+7.47363091, +0.879458010,
257+7.47363091, -0.879458010,
258+5.34100008, +0.00000000,
259+5.26379681, +0.00000000,
260+4.04746771, +0.00000000,
261+2.98188233, +0.00000000,
262+2.00178838, +0.00000000,
263+0.999977231, +0.00000000
264getPolyVal(coef, root(1:count))
265-68181.5000, -14448.8984,
266-68181.5000, +14448.8984,
267-8617.75000, -8708.84375,
268-8617.75000, +8708.84375,
269-810.250000, +0.00000000,
270-735.250000, +0.00000000,
271+180.500000, +0.00000000,
272+180.250000, +0.00000000,
273+72.2500000, +0.00000000,
274+8.25000000, +0.00000000
275
276
277coef
278+3628800.0000000000,
279-10628640.000000000,
280+12753576.000000000,
281-8409500.0000000000,
282+3416930.0000000000,
283-902055.00000000000,
284+157773.00000000000,
285-18150.000000000000,
286+1320.0000000000000,
287-55.000000000000000,
288+1.0000000000000000
289call setResized(root, size(coef, 1, IK) - 1_IK)
290call setPolyRoot(root, count, coef, eigen)
291count
292+10
293root(1:count)
294+10.000000000095149, +0.0000000000000000,
295+8.9999999996955644, +0.0000000000000000,
296+8.0000000003491696, +0.0000000000000000,
297+6.9999999998340323, +0.0000000000000000,
298+6.0000000000428484, +0.0000000000000000,
299+4.9999999999586162, +0.0000000000000000,
300+4.0000000000353024, +0.0000000000000000,
301+2.9999999999877933, +0.0000000000000000,
302+2.0000000000016263, +0.0000000000000000,
303+0.99999999999994482, +0.0000000000000000
304getPolyVal(coef, root(1:count))
305+0.19303057342767715E-4, +0.0000000000000000,
306+0.10346993803977966E-4, +0.0000000000000000,
307+0.66170468926429749E-5, +0.0000000000000000,
308+0.89593231678009033E-6, +0.0000000000000000,
309-0.75763091444969177E-6, +0.0000000000000000,
310+0.10197982192039490E-6, +0.0000000000000000,
311+0.14528632164001465E-6, +0.0000000000000000,
312+0.85681676864624023E-7, +0.0000000000000000,
313+0.70314854383468628E-7, +0.0000000000000000,
314+0.19557774066925049E-7, +0.0000000000000000
315
316
317coef
318+3628800.00000000000000000000000000000,
319-10628640.0000000000000000000000000000,
320+12753576.0000000000000000000000000000,
321-8409500.00000000000000000000000000000,
322+3416930.00000000000000000000000000000,
323-902055.000000000000000000000000000000,
324+157773.000000000000000000000000000000,
325-18150.0000000000000000000000000000000,
326+1320.00000000000000000000000000000000,
327-55.0000000000000000000000000000000000,
328+1.00000000000000000000000000000000000
329call setResized(root, size(coef, 1, IK) - 1_IK)
330call setPolyRoot(root, count, coef, eigen)
331count
332+10
333root(1:count)
334+10.0000000000000000000000000000517967, +0.00000000000000000000000000000000000,
335+8.99999999999999999999999999958200849, +0.00000000000000000000000000000000000,
336+8.00000000000000000000000000123870575, +0.00000000000000000000000000000000000,
337+6.99999999999999999999999999814682938, +0.00000000000000000000000000000000000,
338+6.00000000000000000000000000155546576, +0.00000000000000000000000000000000000,
339+4.99999999999999999999999999924843974, +0.00000000000000000000000000000000000,
340+4.00000000000000000000000000020351610, +0.00000000000000000000000000000000000,
341+2.99999999999999999999999999997156942, +0.00000000000000000000000000000000000,
342+2.00000000000000000000000000000172640, +0.00000000000000000000000000000000000,
343+0.999999999999999999999999999999968030, +0.00000000000000000000000000000000000
344getPolyVal(coef, root(1:count))
345+0.190570619346140160075476935146177038E-22, +0.00000000000000000000000000000000000,
346+0.106701452257938892161879271616720111E-22, +0.00000000000000000000000000000000000,
347+0.119020304153870212515188824892176100E-22, +0.00000000000000000000000000000000000,
348+0.954246540633683195630463953068600702E-23, +0.00000000000000000000000000000000000,
349+0.456645703394752484965632161745263673E-23, +0.00000000000000000000000000000000000,
350+0.212247759715144552316858040860725332E-23, +0.00000000000000000000000000000000000,
351+0.846163761376266102956836528264927821E-24, +0.00000000000000000000000000000000000,
352+0.307365452223073271766182624348262564E-24, +0.00000000000000000000000000000000000,
353+0.666429692730710773211828291950897807E-25, +0.00000000000000000000000000000000000,
354+0.117130067207215832867533457373188099E-25, +0.00000000000000000000000000000000000
355
356
357coef
358+3628800.00,
359-10628640.0,
360+12753576.0,
361-8409500.00,
362+3416930.00,
363-902055.000,
364+157773.000,
365-18150.0000,
366+1320.00000,
367-55.0000000,
368+1.00000000
369call setResized(root, size(coef, 1, IK) - 1_IK)
370call setPolyRoot(root, count, coef, jenkins)
371count
372+10
373root(1:count)
374+0.999999166, +0.00000000,
375+2.00011921, +0.00000000,
376+2.99888420, +0.00000000,
377+4.13503027, +0.00000000,
378+4.58229208, +0.00000000,
379+6.40520525, +0.821048856,
380+6.40520525, -0.821048856,
381+8.68674374, +0.593306661,
382+8.68674374, -0.593306661,
383+10.0997829, +0.00000000
384getPolyVal(coef, root(1:count))
385+0.00000000, +0.00000000,
386+7.75000000, +0.00000000,
387+9.25000000, +0.00000000,
388+521.750000, +0.00000000,
389+994.750000, +0.00000000,
390+4793.75000, +3270.65625,
391+4793.75000, -3270.65625,
392+22855.0000, +7855.48438,
393+22855.0000, -7855.48438,
394+57870.2500, +0.00000000
395
396
397coef
398+3628800.0000000000,
399-10628640.000000000,
400+12753576.000000000,
401-8409500.0000000000,
402+3416930.0000000000,
403-902055.00000000000,
404+157773.00000000000,
405-18150.000000000000,
406+1320.0000000000000,
407-55.000000000000000,
408+1.0000000000000000
409call setResized(root, size(coef, 1, IK) - 1_IK)
410call setPolyRoot(root, count, coef, jenkins)
411count
412+10
413root(1:count)
414+0.99999999999999689, +0.0000000000000000,
415+1.9999999999996774, +0.0000000000000000,
416+3.0000000000025926, +0.0000000000000000,
417+3.9999999999885669, +0.0000000000000000,
418+5.0000000000394653, +0.0000000000000000,
419+5.9999999999004254, +0.0000000000000000,
420+7.0000000001635714, +0.0000000000000000,
421+7.9999999998370157, +0.0000000000000000,
422+9.0000000000894520, +0.0000000000000000,
423+9.9999999999792326, +0.0000000000000000
424getPolyVal(coef, root(1:count))
425+0.93132257461547852E-9, +0.0000000000000000,
426-0.14435499906539917E-7, +0.0000000000000000,
427-0.69849193096160889E-8, +0.0000000000000000,
428-0.13317912817001343E-6, +0.0000000000000000,
429-0.15320256352424622E-6, +0.0000000000000000,
430-0.32829120755195618E-6, +0.0000000000000000,
431+0.20023435354232788E-7, +0.0000000000000000,
432+0.18826685845851898E-5, +0.0000000000000000,
433-0.71660615503787994E-5, +0.0000000000000000,
434-0.11947005987167358E-4, +0.0000000000000000
435
436
437coef
438+3628800.00000000000000000000000000000,
439-10628640.0000000000000000000000000000,
440+12753576.0000000000000000000000000000,
441-8409500.00000000000000000000000000000,
442+3416930.00000000000000000000000000000,
443-902055.000000000000000000000000000000,
444+157773.000000000000000000000000000000,
445-18150.0000000000000000000000000000000,
446+1320.00000000000000000000000000000000,
447-55.0000000000000000000000000000000000,
448+1.00000000000000000000000000000000000
449call setResized(root, size(coef, 1, IK) - 1_IK)
450call setPolyRoot(root, count, coef, jenkins)
451count
452+10
453root(1:count)
454+0.999999999999999999999999999999999037, +0.00000000000000000000000000000000000,
455+1.99999999999999999999999999999988136, +0.00000000000000000000000000000000000,
456+3.00000000000000000000000000000386380, +0.00000000000000000000000000000000000,
457+3.99999999999999999999999999995968220, +0.00000000000000000000000000000000000,
458+5.00000000000000000000000000019942003, +0.00000000000000000000000000000000000,
459+5.99999999999999999999999999945674139, +0.00000000000000000000000000000000000,
460+7.00000000000000000000000000086131439, +0.00000000000000000000000000000000000,
461+7.99999999999999999999999999920843509, +0.00000000000000000000000000000000000,
462+9.00000000000000000000000000039093296, +0.00000000000000000000000000000000000,
463+9.99999999999999999999999999991973032, +0.00000000000000000000000000000000000
464getPolyVal(coef, root(1:count))
465+0.121169035041947413311241507627435965E-26, +0.00000000000000000000000000000000000,
466-0.767403888599000284304529548307094444E-26, +0.00000000000000000000000000000000000,
467-0.496793043671984394576090181272487456E-25, +0.00000000000000000000000000000000000,
468-0.226586095528441662892021619263305254E-24, +0.00000000000000000000000000000000000,
469-0.483464449817370179111853615433469500E-24, +0.00000000000000000000000000000000000,
470-0.949157441161921404271391809748248392E-24, +0.00000000000000000000000000000000000,
471-0.437581775214819425271330164545213748E-23, +0.00000000000000000000000000000000000,
472-0.847294672369990945481074782336117224E-23, +0.00000000000000000000000000000000000,
473-0.200878065260375153367819879395017590E-22, +0.00000000000000000000000000000000000,
474-0.137656101743321724682457768765275752E-22, +0.00000000000000000000000000000000000
475
476
477coef
478+3628800.00,
479-10628640.0,
480+12753576.0,
481-8409500.00,
482+3416930.00,
483-902055.000,
484+157773.000,
485-18150.0000,
486+1320.00000,
487-55.0000000,
488+1.00000000
489call setResized(root, size(coef, 1, IK) - 1_IK)
490call setPolyRoot(root, count, coef, laguerre)
491count
492+10
493root(1:count)
494+1.00000036, +0.00000000,
495+1.99991179, +0.00000000,
496+3.00164771, -0.279217958E-2,
497+3.00050926, +0.989809632E-4,
498+2.00066686, -0.105152279E-2,
499+2.99986410, +0.315787271E-3,
500+2.00000215, +0.157766044E-5,
501+1.99925506, +0.763451681E-3,
502+2.99917603, -0.601705164E-3,
503+1.00000143, -0.204890966E-7
504getPolyVal(coef, root(1:count))
505+0.250000000, +0.00000000,
506-6.00000000, +0.00000000,
507-15.0000000, +28.0310059,
508-4.00000000, -0.995964050,
509+29.5000000, -42.3060303,
510-0.500000000, -3.18444824,
511-2.00000000, +0.636050701E-1,
512-22.2500000, +30.8698730,
513-8.25000000, +6.09149170,
514-0.250000000, +0.743500888E-2
515
516
517coef
518+3628800.0000000000,
519-10628640.000000000,
520+12753576.000000000,
521-8409500.0000000000,
522+3416930.0000000000,
523-902055.00000000000,
524+157773.00000000000,
525-18150.000000000000,
526+1320.0000000000000,
527-55.000000000000000,
528+1.0000000000000000
529call setResized(root, size(coef, 1, IK) - 1_IK)
530call setPolyRoot(root, count, coef, laguerre)
531count
532+10
533root(1:count)
534+0.99999999999999911, +0.0000000000000000,
535+2.0000000000000457, +0.0000000000000000,
536+2.9999999999974793, -0.28736254480092181E-19,
537+3.0000000000026730, -0.25939867021780102E-22,
538+2.0000000000001097, +0.66691436887087856E-23,
539+3.0000000000403091, -0.27597116866034943E-10,
540+1.9999999999999909, -0.77452692696933223E-17,
541+2.0000000000001705, -0.11580528575742387E-22,
542+3.0000000000005125, -0.38246763572920617E-21,
543+1.0000000000000016, +0.26469779601696886E-22
544getPolyVal(coef, root(1:count))
545-0.46566128730773926E-9, +0.0000000000000000,
546+0.27939677238464355E-8, +0.0000000000000000,
547+0.26077032089233398E-7, +0.28966144516097633E-15,
548-0.31199306249618530E-7, +0.26147385957797803E-18,
549+0.23283064365386963E-8, +0.26889987352861344E-18,
550-0.42235478758811951E-6, +0.27817893798584754E-6,
551-0.97788870334625244E-8, -0.31228925695390712E-12,
552+0.74505805969238281E-8, -0.46692691217367730E-18,
553+0.17229467630386353E-7, +0.38552737681339307E-17,
554+0.0000000000000000, -0.96053536218636303E-17
555
556
557coef
558+3628800.00000000000000000000000000000,
559-10628640.0000000000000000000000000000,
560+12753576.0000000000000000000000000000,
561-8409500.00000000000000000000000000000,
562+3416930.00000000000000000000000000000,
563-902055.000000000000000000000000000000,
564+157773.000000000000000000000000000000,
565-18150.0000000000000000000000000000000,
566+1320.00000000000000000000000000000000,
567-55.0000000000000000000000000000000000,
568+1.00000000000000000000000000000000000
569call setResized(root, size(coef, 1, IK) - 1_IK)
570call setPolyRoot(root, count, coef, laguerre)
571count
572+10
573root(1:count)
574+1.00000000000000000000000000000000077, +0.00000000000000000000000000000000000,
575+1.99999999999999999999999999999994280, +0.00000000000000000000000000000000000,
576+3.00000000000000000000000000000117328, -0.184186914397489602460971036068311768E-52,
577+2.99999999999999999999999999999705718, +0.116740925043573719397266623324820008E-29,
578+1.99999999999999999999999999999993086, +0.787211432649434425953724102803185716E-58,
579+2.99999999999999999999999999999986673, -0.114087307158325352687449640351479862E-30,
580+2.00000000000000000000000000000020530, -0.200457353256914673460540235036361141E-49,
581+2.00000000000000000000000000000010169, -0.700089218759378442053707206050263977E-59,
582+3.00000000000000000000000000000239316, +0.998793952096713243996622280631709941E-58,
583+1.00000000000000000000000000000000135, -0.509789411562384728649241727285667778E-56
584getPolyVal(coef, root(1:count))
585+0.403896783473158044370805025424786550E-27, +0.00000000000000000000000000000000000,
586+0.403896783473158044370805025424786550E-27, +0.00000000000000000000000000000000000,
587+0.565455496862421262119127035594701169E-26, +0.185660409712669519280658804356040188E-48,
588+0.214065295240773763516526663475136871E-25, -0.117674852443922309152444756312348934E-25,
589-0.363507105125842239933724522882307895E-26, +0.317403649644251960544541558250312043E-53,
590-0.363507105125842239933724522882307895E-26, +0.115000005615591955508949237474271788E-26,
591+0.565455496862421262119127035594701169E-26, -0.808244048331879963392898227665960408E-45,
592+0.282727748431210631059563517797350585E-26, -0.282275973003781387836054745479349811E-54,
593-0.254454973588089567953607166017615526E-25, -0.100678430371348694994859525887176084E-53,
594+0.00000000000000000000000000000000000, +0.184992381667758170332236837997420858E-50
595
596
597!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
598! Compute the roots of polynomials with complex coefficients with zeros 1,2,...,10.
599!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
600
601
602coef
603+3628800.00, +0.00000000,
604-10628640.0, +0.00000000,
605+12753576.0, +0.00000000,
606-8409500.00, +0.00000000,
607+3416930.00, +0.00000000,
608-902055.000, +0.00000000,
609+157773.000, +0.00000000,
610-18150.0000, +0.00000000,
611+1320.00000, +0.00000000,
612-55.0000000, +0.00000000,
613+1.00000000, +0.00000000
614call setResized(root, size(coef, 1, IK) - 1_IK)
615call setPolyRoot(root, count, coef, eigen)
616count
617+10
618root(1:count)
619+9.99373722, +0.196921639E-2,
620+9.07808781, -0.768133765E-2,
621+7.58786488, +0.259252965,
622+7.53517246, -0.258893400,
623+5.48630047, +0.496490821E-1,
624+5.36296797, -0.451893099E-1,
625+3.94950104, +0.102348044E-2,
626+3.00655532, -0.139494587E-3,
627+1.99981821, +0.898552662E-5,
628+0.999996662, -0.528134251E-6
629getPolyVal(coef, root(1:count))
630-4153.50000, +685.126953,
631+1293.75000, +363.354248,
632-3629.75000, -342.171875,
633-2187.25000, +294.710938,
634-1123.75000, -18.1875000,
635-459.750000, +35.3359375,
636-246.750000, +4.63513184,
637-61.2500000, +1.38377380,
638-6.00000000, +0.362543106,
639+0.500000000, +0.191653848
640
641
642coef
643+3628800.0000000000, +0.0000000000000000,
644-10628640.000000000, +0.0000000000000000,
645+12753576.000000000, +0.0000000000000000,
646-8409500.0000000000, +0.0000000000000000,
647+3416930.0000000000, +0.0000000000000000,
648-902055.00000000000, +0.0000000000000000,
649+157773.00000000000, +0.0000000000000000,
650-18150.000000000000, +0.0000000000000000,
651+1320.0000000000000, +0.0000000000000000,
652-55.000000000000000, +0.0000000000000000,
653+1.0000000000000000, +0.0000000000000000
654call setResized(root, size(coef, 1, IK) - 1_IK)
655call setPolyRoot(root, count, coef, eigen)
656count
657+10
658root(1:count)
659+10.000000000092832, +0.14674683172918083E-11,
660+8.9999999995494910, -0.78513163716371057E-11,
661+8.0000000009219665, +0.18126945817416554E-10,
662+6.9999999989793631, -0.23471195338655758E-10,
663+6.0000000006482752, +0.18575051340220250E-10,
664+4.9999999997713829, -0.92418512141289470E-11,
665+4.0000000000382663, +0.28926396100133261E-11,
666+2.9999999999985825, -0.54892097595564787E-12,
667+1.9999999999998195, +0.51688262166740718E-13,
668+1.0000000000000040, -0.69823506140973722E-15
669getPolyVal(coef, root(1:count))
670+0.21335668861865997E-4, +0.53251490324658425E-6,
671+0.15214085578918457E-4, +0.31656507562992751E-6,
672+0.12007541954517365E-4, +0.18271961426412876E-6,
673+0.43502077460289001E-5, +0.10139556374739953E-6,
674+0.16265548765659332E-5, +0.53496147868086052E-7,
675+0.96531584858894348E-6, +0.26616531495547122E-7,
676+0.16111880540847778E-6, +0.12496203114631699E-7,
677-0.93132257461547852E-9, +0.55331234376599884E-8,
678-0.79162418842315674E-8, +0.20840707305641930E-8,
679-0.23283064365386963E-8, +0.25337553908436255E-9
680
681
682coef
683+3628800.00000000000000000000000000000, +0.00000000000000000000000000000000000,
684-10628640.0000000000000000000000000000, +0.00000000000000000000000000000000000,
685+12753576.0000000000000000000000000000, +0.00000000000000000000000000000000000,
686-8409500.00000000000000000000000000000, +0.00000000000000000000000000000000000,
687+3416930.00000000000000000000000000000, +0.00000000000000000000000000000000000,
688-902055.000000000000000000000000000000, +0.00000000000000000000000000000000000,
689+157773.000000000000000000000000000000, +0.00000000000000000000000000000000000,
690-18150.0000000000000000000000000000000, +0.00000000000000000000000000000000000,
691+1320.00000000000000000000000000000000, +0.00000000000000000000000000000000000,
692-55.0000000000000000000000000000000000, +0.00000000000000000000000000000000000,
693+1.00000000000000000000000000000000000, +0.00000000000000000000000000000000000
694call setResized(root, size(coef, 1, IK) - 1_IK)
695call setPolyRoot(root, count, coef, eigen)
696count
697+10
698root(1:count)
699+10.0000000000000000000000000001225200, +0.509435088823196545261926980817026425E-30,
700+8.99999999999999999999999999940435917, -0.209995537749767662358226781607926446E-29,
701+8.00000000000000000000000000121674245, +0.343592927705395245870133814335134133E-29,
702+6.99999999999999999999999999865599133, -0.302445729136684007599014761793556296E-29,
703+6.00000000000000000000000000086041537, +0.192719830774026860653879227602749694E-29,
704+4.99999999999999999999999999968191341, -0.123024474772112043974718002090415763E-29,
705+4.00000000000000000000000000006421590, +0.652456748571305541599111900496149722E-30,
706+2.99999999999999999999999999999352079, -0.193766393931633696713588827264264386E-30,
707+2.00000000000000000000000000000033781, +0.237223756726557254262497294342951649E-31,
708+0.999999999999999999999999999999990659, -0.651188212631313945515795127577401738E-33
709getPolyVal(coef, root(1:count))
710+0.382902228668223289224410580203206145E-22, +0.184863805032161562344648062924845379E-24,
711+0.258530292133333732621308588724151623E-22, +0.846702008207063214628370381684629456E-25,
712+0.141743537192070084091490315622574592E-22, +0.346341671127038407837094885825177763E-25,
713+0.442428536616497321803779824850311186E-23, +0.130656554987047491282774376921933986E-25,
714+0.197545916796721599501760737935263101E-23, +0.555033112629197358683172175557306576E-26,
715+0.891804097908732961970737496137928702E-24, +0.354310487343682686647187846065246503E-26,
716+0.246780934702099565110561870534544582E-24, +0.281861315382803993970816340989114052E-26,
717+0.727014210251684479867449045764615789E-25, +0.195316525083086766287297537885110653E-26,
718+0.137324906380873735086073708644427427E-25, +0.956486187121478849186389090789675841E-27,
719+0.323117426778526435496644020339829240E-26, +0.236303178599651204548771735895300922E-27
720
721
722coef
723+3628800.00, +0.00000000,
724-10628640.0, +0.00000000,
725+12753576.0, +0.00000000,
726-8409500.00, +0.00000000,
727+3416930.00, +0.00000000,
728-902055.000, +0.00000000,
729+157773.000, +0.00000000,
730-18150.0000, +0.00000000,
731+1320.00000, +0.00000000,
732-55.0000000, +0.00000000,
733+1.00000000, +0.00000000
734call setResized(root, size(coef, 1, IK) - 1_IK)
735call setPolyRoot(root, count, coef, jenkins)
736count
737+10
738root(1:count)
739+0.999999464, -0.670552254E-7,
740+2.00012732, -0.640330836E-4,
741+3.00413203, +0.198328402E-2,
742+3.98040628, -0.952805020E-2,
743+5.38614893, -0.139350891E-1,
744+5.33709526, +0.364396721E-1,
745+7.57813692, +0.669925332,
746+7.48616934, -0.654946744,
747+9.30202389, -0.568577349E-1,
748+9.92576122, +0.269834101E-1
749getPolyVal(coef, root(1:count))
750+0.00000000, +0.243330747E-1,
751+1.00000000, -2.58021545,
752-36.7500000, -19.8110352,
753-117.500000, -41.7656250,
754-589.250000, +9.86035156,
755-887.000000, -49.8808594,
756-7647.50000, -2649.03125,
757-7193.50000, +4439.71875,
758-22343.5000, +4553.98828,
759-17839.5000, +6221.86230
760
761
762coef
763+3628800.0000000000, +0.0000000000000000,
764-10628640.000000000, +0.0000000000000000,
765+12753576.000000000, +0.0000000000000000,
766-8409500.0000000000, +0.0000000000000000,
767+3416930.0000000000, +0.0000000000000000,
768-902055.00000000000, +0.0000000000000000,
769+157773.00000000000, +0.0000000000000000,
770-18150.000000000000, +0.0000000000000000,
771+1320.0000000000000, +0.0000000000000000,
772-55.000000000000000, +0.0000000000000000,
773+1.0000000000000000, +0.0000000000000000
774call setResized(root, size(coef, 1, IK) - 1_IK)
775call setPolyRoot(root, count, coef, jenkins)
776count
777+10
778root(1:count)
779+0.99999999999999767, +0.78615245417039750E-19,
780+2.0000000000000497, -0.70753769989184646E-18,
781+3.0000000000002629, +0.38406961987792537E-16,
782+3.9999999999952496, +0.45047880627580722E-12,
783+5.0000000000141975, -0.27036496836422169E-11,
784+6.9999999999270663, -0.90134373226557985E-11,
785+6.0000000000022871, +0.67579228388490423E-11,
786+8.0000000001275087, +0.68963099848637832E-11,
787+8.9999999999098499, -0.29688875676492067E-11,
788+10.000000000023537, +0.58122516591905612E-12
789getPolyVal(coef, root(1:count))
790+0.0000000000000000, -0.28527900256935788E-13,
791+0.27939677238464355E-8, -0.28527920059632989E-13,
792-0.21886080503463745E-7, -0.38714217683733651E-12,
793-0.26542693376541138E-7, +0.19460684431239280E-8,
794+0.13597309589385986E-6, +0.77865110883786122E-8,
795+0.14705583453178406E-5, +0.38938049220170988E-7,
796+0.13969838619232178E-7, +0.19462817776330781E-7,
797-0.69197267293930054E-6, +0.69514804654671607E-7,
798+0.34747645258903503E-5, +0.11970554668839881E-6,
799+0.23906119167804718E-4, +0.21091498824245076E-6
800
801
802coef
803+3628800.00000000000000000000000000000, +0.00000000000000000000000000000000000,
804-10628640.0000000000000000000000000000, +0.00000000000000000000000000000000000,
805+12753576.0000000000000000000000000000, +0.00000000000000000000000000000000000,
806-8409500.00000000000000000000000000000, +0.00000000000000000000000000000000000,
807+3416930.00000000000000000000000000000, +0.00000000000000000000000000000000000,
808-902055.000000000000000000000000000000, +0.00000000000000000000000000000000000,
809+157773.000000000000000000000000000000, +0.00000000000000000000000000000000000,
810-18150.0000000000000000000000000000000, +0.00000000000000000000000000000000000,
811+1320.00000000000000000000000000000000, +0.00000000000000000000000000000000000,
812-55.0000000000000000000000000000000000, +0.00000000000000000000000000000000000,
813+1.00000000000000000000000000000000000, +0.00000000000000000000000000000000000
814call setResized(root, size(coef, 1, IK) - 1_IK)
815call setPolyRoot(root, count, coef, jenkins)
816count
817+10
818root(1:count)
819+0.999999999999999999999999999999999904, -0.184379003047766309026726070351153008E-49,
820+1.99999999999999999999999999999999827, +0.165941102996888701460788044967515520E-48,
821+3.00000000000000000000000000000014252, -0.818772228890815727574784569242905836E-43,
822+4.00000000000000000000000000000214317, -0.393879829577633251033567677421169831E-32,
823+4.99999999999999999999999999997773394, +0.236327897763774097689509421466453076E-31,
824+6.99999999999999999999999999987130628, +0.787759659241237182996329756600879293E-31,
825+6.00000000000000000000000000007577533, -0.590819612978407215382047079706979614E-31,
826+8.00000000000000000000000000011919812, -0.605430221624447494817117793220965173E-31,
827+8.99999999999999999999999999994202335, +0.265548326372253158791183154847662911E-31,
828+10.0000000000000000000000000000116804, -0.539980658158276334206419034384794310E-32
829getPolyVal(coef, root(1:count))
830+0.403896783473158044370805025424786550E-27, +0.669074526259734382196183564090262307E-44,
831-0.323117426778526435496644020339829240E-26, +0.669074527283455244289897397308953489E-44,
832+0.686624531904368675430368543222137134E-26, +0.825322406721942253395382845795219061E-39,
833+0.201948391736579022185402512712393275E-25, -0.170156086377537564446501236645746230E-28,
834+0.238299102249163246178774965000624064E-24, -0.680624345559669401345787133771733833E-28,
835-0.563032116161582313852902205442152450E-24, -0.340312172792214463054414454856269519E-27,
836+0.892611891475679278059479106188778275E-25, -0.170156048537781278030029558955455798E-27,
837+0.926135324503951395742255923299035558E-24, -0.610273663397443074775654735723737081E-27,
838-0.190922009547761807574079535518296602E-23, -0.107069085193292473624605048018225966E-26,
839+0.532578298687706197307343506525123544E-23, -0.195948181232475316156825339210492741E-26
840
841
842coef
843+3628800.00, +0.00000000,
844-10628640.0, +0.00000000,
845+12753576.0, +0.00000000,
846-8409500.00, +0.00000000,
847+3416930.00, +0.00000000,
848-902055.000, +0.00000000,
849+157773.000, +0.00000000,
850-18150.0000, +0.00000000,
851+1320.00000, +0.00000000,
852-55.0000000, +0.00000000,
853+1.00000000, +0.00000000
854call setResized(root, size(coef, 1, IK) - 1_IK)
855call setPolyRoot(root, count, coef, laguerre)
856count
857+10
858root(1:count)
859+1.00000036, +0.00000000,
860+1.99991179, +0.00000000,
861+3.00164771, -0.279217958E-2,
862+3.00050926, +0.989809632E-4,
863+2.00066686, -0.105152279E-2,
864+2.99986410, +0.315787271E-3,
865+2.00000215, +0.157766044E-5,
866+1.99925506, +0.763451681E-3,
867+2.99917603, -0.601705164E-3,
868+1.00000143, -0.204890966E-7
869getPolyVal(coef, root(1:count))
870+0.250000000, +0.00000000,
871-6.00000000, +0.00000000,
872-15.0000000, +28.0310059,
873-4.00000000, -0.995964050,
874+29.5000000, -42.3060303,
875-0.500000000, -3.18444824,
876-2.00000000, +0.636050701E-1,
877-22.2500000, +30.8698730,
878-8.25000000, +6.09149170,
879-0.250000000, +0.743500888E-2
880
881
882coef
883+3628800.0000000000, +0.0000000000000000,
884-10628640.000000000, +0.0000000000000000,
885+12753576.000000000, +0.0000000000000000,
886-8409500.0000000000, +0.0000000000000000,
887+3416930.0000000000, +0.0000000000000000,
888-902055.00000000000, +0.0000000000000000,
889+157773.00000000000, +0.0000000000000000,
890-18150.000000000000, +0.0000000000000000,
891+1320.0000000000000, +0.0000000000000000,
892-55.000000000000000, +0.0000000000000000,
893+1.0000000000000000, +0.0000000000000000
894call setResized(root, size(coef, 1, IK) - 1_IK)
895call setPolyRoot(root, count, coef, laguerre)
896count
897+10
898root(1:count)
899+0.99999999999999911, +0.0000000000000000,
900+2.0000000000000457, +0.0000000000000000,
901+2.9999999999974793, -0.28736254480092181E-19,
902+3.0000000000026730, -0.25939867021780102E-22,
903+2.0000000000001097, +0.66691436887087856E-23,
904+3.0000000000403091, -0.27597116866034943E-10,
905+1.9999999999999909, -0.77452692696933223E-17,
906+2.0000000000001705, -0.11580528575742387E-22,
907+3.0000000000005125, -0.38246763572920617E-21,
908+1.0000000000000016, +0.26469779601696886E-22
909getPolyVal(coef, root(1:count))
910-0.46566128730773926E-9, +0.0000000000000000,
911+0.27939677238464355E-8, +0.0000000000000000,
912+0.26077032089233398E-7, +0.28966144516097633E-15,
913-0.31199306249618530E-7, +0.26147385957797803E-18,
914+0.23283064365386963E-8, +0.26889987352861344E-18,
915-0.42235478758811951E-6, +0.27817893798584754E-6,
916-0.97788870334625244E-8, -0.31228925695390712E-12,
917+0.74505805969238281E-8, -0.46692691217367730E-18,
918+0.17229467630386353E-7, +0.38552737681339307E-17,
919+0.0000000000000000, -0.96053536218636303E-17
920
921
922coef
923+3628800.00000000000000000000000000000, +0.00000000000000000000000000000000000,
924-10628640.0000000000000000000000000000, +0.00000000000000000000000000000000000,
925+12753576.0000000000000000000000000000, +0.00000000000000000000000000000000000,
926-8409500.00000000000000000000000000000, +0.00000000000000000000000000000000000,
927+3416930.00000000000000000000000000000, +0.00000000000000000000000000000000000,
928-902055.000000000000000000000000000000, +0.00000000000000000000000000000000000,
929+157773.000000000000000000000000000000, +0.00000000000000000000000000000000000,
930-18150.0000000000000000000000000000000, +0.00000000000000000000000000000000000,
931+1320.00000000000000000000000000000000, +0.00000000000000000000000000000000000,
932-55.0000000000000000000000000000000000, +0.00000000000000000000000000000000000,
933+1.00000000000000000000000000000000000, +0.00000000000000000000000000000000000
934call setResized(root, size(coef, 1, IK) - 1_IK)
935call setPolyRoot(root, count, coef, laguerre)
936count
937+10
938root(1:count)
939+1.00000000000000000000000000000000077, +0.00000000000000000000000000000000000,
940+1.99999999999999999999999999999994280, +0.00000000000000000000000000000000000,
941+3.00000000000000000000000000000117328, -0.184186914397489602460971036068311768E-52,
942+2.99999999999999999999999999999705718, +0.116740925043573719397266623324820008E-29,
943+1.99999999999999999999999999999993086, +0.787211432649434425953724102803185716E-58,
944+2.99999999999999999999999999999986673, -0.114087307158325352687449640351479862E-30,
945+2.00000000000000000000000000000020530, -0.200457353256914673460540235036361141E-49,
946+2.00000000000000000000000000000010169, -0.700089218759378442053707206050263977E-59,
947+3.00000000000000000000000000000239316, +0.998793952096713243996622280631709941E-58,
948+1.00000000000000000000000000000000135, -0.509789411562384728649241727285667778E-56
949getPolyVal(coef, root(1:count))
950+0.403896783473158044370805025424786550E-27, +0.00000000000000000000000000000000000,
951+0.403896783473158044370805025424786550E-27, +0.00000000000000000000000000000000000,
952+0.565455496862421262119127035594701169E-26, +0.185660409712669519280658804356040188E-48,
953+0.214065295240773763516526663475136871E-25, -0.117674852443922309152444756312348934E-25,
954-0.363507105125842239933724522882307895E-26, +0.317403649644251960544541558250312043E-53,
955-0.363507105125842239933724522882307895E-26, +0.115000005615591955508949237474271788E-26,
956+0.565455496862421262119127035594701169E-26, -0.808244048331879963392898227665960408E-45,
957+0.282727748431210631059563517797350585E-26, -0.282275973003781387836054745479349811E-54,
958-0.254454973588089567953607166017615526E-25, -0.100678430371348694994859525887176084E-53,
959+0.00000000000000000000000000000000000, +0.184992381667758170332236837997420858E-50
960
961
962!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
963! Compute the roots of polynomials with complex coefficients with zeros on imaginary axis.
964!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
965
966
967coef
968+0.0000000000000000, +1.0000000000000000,
969-10001.000099999999, +0.0000000000000000,
970+0.0000000000000000, -10001.000099999999,
971+1.0000000000000000, +0.0000000000000000
972call setResized(root, size(coef, 1, IK) - 1_IK)
973call setPolyRoot(root, count, coef, eigen)
974count
975+3
976root(1:count)
977+0.0000000000000000, +10000.000000000000,
978+0.0000000000000000, +1.0000000000000004,
979+0.0000000000000000, +0.10000000000000255E-3
980getPolyVal(coef, root(1:count))
981+0.0000000000000000, -0.70715235779061913E-4,
982+0.0000000000000000, +0.36375347178818629E-11,
983+0.0000000000000000, -0.25535129566378600E-13
984
985
986coef
987+0.00000000000000000000000000000000000, +1.00000000000000000000000000000000000,
988-10001.0000999999999999999999999999996, +0.00000000000000000000000000000000000,
989+0.00000000000000000000000000000000000, -10001.0000999999999999999999999999996,
990+1.00000000000000000000000000000000000, +0.00000000000000000000000000000000000
991call setResized(root, size(coef, 1, IK) - 1_IK)
992call setPolyRoot(root, count, coef, eigen)
993count
994+3
995root(1:count)
996+0.00000000000000000000000000000000000, +9999.99999999999999999999999999999842,
997+0.00000000000000000000000000000000000, +0.999999999999999999999999999999999904,
998+0.00000000000000000000000000000000000, +0.100000000000000000000000000000000008E-3
999getPolyVal(coef, root(1:count))
1000+0.00000000000000000000000000000000000, +0.115351482477835407341207567853721561E-21,
1001+0.00000000000000000000000000000000000, -0.157752921744758488723815153277131397E-29,
1002+0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000
1003
1004
1005coef
1006+0.0000000000000000, +1.0000000000000000,
1007-10001.000099999999, +0.0000000000000000,
1008+0.0000000000000000, -10001.000099999999,
1009+1.0000000000000000, +0.0000000000000000
1010call setResized(root, size(coef, 1, IK) - 1_IK)
1011call setPolyRoot(root, count, coef, jenkins)
1012count
1013+3
1014root(1:count)
1015+0.0000000000000000, +0.10000000000000000E-3,
1016+0.0000000000000000, +1.0000000000000000,
1017-0.0000000000000000, +10000.000000000000
1018getPolyVal(coef, root(1:count))
1019+0.0000000000000000, +0.0000000000000000,
1020+0.0000000000000000, +0.0000000000000000,
1021+0.0000000000000000, -0.70715235779061913E-4
1022
1023
1024coef
1025+0.00000000000000000000000000000000000, +1.00000000000000000000000000000000000,
1026-10001.0000999999999999999999999999996, +0.00000000000000000000000000000000000,
1027+0.00000000000000000000000000000000000, -10001.0000999999999999999999999999996,
1028+1.00000000000000000000000000000000000, +0.00000000000000000000000000000000000
1029call setResized(root, size(coef, 1, IK) - 1_IK)
1030call setPolyRoot(root, count, coef, jenkins)
1031count
1032+3
1033root(1:count)
1034-0.296736492054993710858538820923811161E-66, +0.100000000000000000000000000000000008E-3,
1035+0.00000000000000000000000000000000000, +1.00000000000000000000000000000000019,
1036+0.296736492054993710858538820923811161E-66, +10000.0000000000000000000000000000000
1037getPolyVal(coef, root(1:count))
1038+0.296706815438720027429570852167189098E-62, +0.00000000000000000000000000000000000,
1039+0.00000000000000000000000000000000000, +0.157752921744758488723815153277131397E-29,
1040-0.296706815438720027429570852167189128E-58, -0.424049213484551264994686342168356274E-22
1041
1042
1043coef
1044+0.0000000000000000, +1.0000000000000000,
1045-10001.000099999999, +0.0000000000000000,
1046+0.0000000000000000, -10001.000099999999,
1047+1.0000000000000000, +0.0000000000000000
1048call setResized(root, size(coef, 1, IK) - 1_IK)
1049call setPolyRoot(root, count, coef, laguerre)
1050count
1051+3
1052root(1:count)
1053+0.0000000000000000, +0.10000000000000000E-3,
1054+0.0000000000000000, +1.0000000000000000,
1055-0.10339757656912846E-24, +1.0000000000000000
1056getPolyVal(coef, root(1:count))
1057+0.0000000000000000, +0.0000000000000000,
1058+0.0000000000000000, +0.0000000000000000,
1059-0.10337689808779039E-20, +0.0000000000000000
1060
1061
1062coef
1063+0.00000000000000000000000000000000000, +1.00000000000000000000000000000000000,
1064-10001.0000999999999999999999999999996, +0.00000000000000000000000000000000000,
1065+0.00000000000000000000000000000000000, -10001.0000999999999999999999999999996,
1066+1.00000000000000000000000000000000000, +0.00000000000000000000000000000000000
1067call setResized(root, size(coef, 1, IK) - 1_IK)
1068call setPolyRoot(root, count, coef, laguerre)
1069count
1070+3
1071root(1:count)
1072+0.00000000000000000000000000000000000, +0.100000000000000000000000000000000008E-3,
1073-0.391601596454741356516116118813642834E-38, +0.999999999999999999999999999999999904,
1074-0.364629801437176271902972503151179155E-62, +0.999999999999999999999999999999999904
1075getPolyVal(coef, root(1:count))
1076+0.00000000000000000000000000000000000, +0.00000000000000000000000000000000000,
1077-0.391523280051466372792226460751041189E-34, -0.157752921744758488723815153277131397E-29,
1078-0.364556879123186851020354627680273858E-58, -0.157752921744758488723815153277131397E-29
1079
1080
1081!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1082! Compute the roots of polynomials with complex coefficients with zeros at 1+i,1/2*(1+i)....1/(2**-9)*(1+i).
1083!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1084
1085
1086coef
1087+0.0000000000000000, +0.90949470177292824E-12,
1088-0.46520653995685279E-9, -0.46520653995685279E-9,
1089+0.15848036127863449E-6, +0.0000000000000000,
1090-0.11546426321729090E-4, +0.11546426321729090E-4,
1091+0.0000000000000000, -0.78207794285845011E-3,
1092+0.12715073651634160E-1, +0.12715073651634160E-1,
1093-0.20021195337176320, +0.0000000000000000,
1094+0.75670659542083740, -0.75670659542083740,
1095+0.0000000000000000, +2.6588592529296879,
1096-1.9980468750000000, -1.9980468750000000,
1097+1.0000000000000000, +0.0000000000000000
1098call setResized(root, size(coef, 1, IK) - 1_IK)
1099call setPolyRoot(root, count, coef, eigen)
1100count
1101+10
1102root(1:count)
1103+0.99999999999999845, +0.99999999999999911,
1104+0.50000000000000133, +0.50000000000000067,
1105+0.25000000000000017, +0.25000000000000028,
1106+0.12499999999999911, +0.12499999999999924,
1107+0.62500000000001193E-1, +0.62500000000001124E-1,
1108+0.31249999999999212E-1, +0.31249999999999153E-1,
1109+0.15625000000000406E-1, +0.15625000000000503E-1,
1110+0.78124999999998820E-2, +0.78124999999998473E-2,
1111+0.39062500000000139E-2, +0.39062500000000173E-2,
1112+0.19531249999999998E-2, +0.19531250000000004E-2
1113getPolyVal(coef, root(1:count))
1114-0.30312292689246898E-14, -0.39170764430220062E-14,
1115-0.68538197226804901E-17, +0.80868759391036751E-17,
1116-0.17941598983115063E-19, +0.37617071523821058E-19,
1117-0.96831376073107395E-22, +0.76998478904539376E-22,
1118+0.18929632499428235E-23, +0.45584800724739300E-23,
1119+0.76689901761965884E-25, +0.68298946085311025E-24,
1120-0.20598735957131060E-25, +0.72903369416905027E-25,
1121-0.16155871338926322E-26, +0.24233807008389483E-26,
1122+0.0000000000000000, -0.80779356694631609E-27,
1123+0.0000000000000000, -0.20194839173657902E-27
1124
1125
1126coef
1127+0.00000000000000000000000000000000000, +0.909494701772928199999999999999999977E-12,
1128-0.465206539956852800000000000000000002E-9, -0.465206539956852800000000000000000002E-9,
1129+0.158480361278634500000000000000000004E-6, +0.00000000000000000000000000000000000,
1130-0.115464263217290900000000000000000001E-4, +0.115464263217290900000000000000000001E-4,
1131+0.00000000000000000000000000000000000, -0.782077942858450099999999999999999970E-3,
1132+0.127150736516341599999999999999999998E-1, +0.127150736516341599999999999999999998E-1,
1133-0.200211953371763200000000000000000000, +0.00000000000000000000000000000000000,
1134+0.756706595420837399999999999999999953, -0.756706595420837399999999999999999953,
1135+0.00000000000000000000000000000000000, +2.65885925292968799999999999999999995,
1136-1.99804687500000000000000000000000000, -1.99804687500000000000000000000000000,
1137+1.00000000000000000000000000000000000, +0.00000000000000000000000000000000000
1138call setResized(root, size(coef, 1, IK) - 1_IK)
1139call setPolyRoot(root, count, coef, eigen)
1140count
1141+10
1142root(1:count)
1143+0.999999999999999160184470589450448173, +0.999999999999999160184470589450448944,
1144+0.500000000000001492398035863052857764, +0.500000000000001492398035863052858149,
1145+0.249999999999999617367224612848345132, +0.249999999999999617367224612848344434,
1146+0.124999999999999149032692721659057207, +0.124999999999999149032692721659057592,
1147+0.625000000000009555071172310128834721E-1, +0.625000000000009555071172310128834842E-1,
1148+0.312499999999994379028135239093267093E-1, +0.312499999999994379028135239093266431E-1,
1149+0.156250000000002691698099555840415078E-1, +0.156250000000002691698099555840415649E-1,
1150+0.781249999999989766148737734007061361E-2, +0.781249999999989766148737734007058653E-2,
1151+0.390625000000002310324757362451072950E-2, +0.390625000000002310324757362451074380E-2,
1152+0.195312499999999767310055151845923676E-2, +0.195312499999999767310055151845923318E-2
1153getPolyVal(coef, root(1:count))
1154-0.330908113531804644550698914860698518E-32, +0.737655648425543017921748959548643300E-33,
1155+0.646304339719146005950884128119047103E-35, -0.199728685471399985252652738182225539E-34,
1156+0.612273361057752245328214721921763582E-37, -0.491287832226650747595597082681272950E-37,
1157+0.380437118781080261768012420837796041E-39, -0.701672705528289773730040710639353438E-39,
1158+0.223655993021642834626120514802966177E-41, -0.822711086520502132446014911659405323E-41,
1159-0.173323103806175811459878802833171501E-42, +0.819759601630017986490381806224600937E-43,
1160+0.516728808719776294903125283838156573E-44, +0.182168800362226219220084845827689097E-43,
1161-0.192678538844662347252012817702363468E-44, +0.262743462060903200798199296866859275E-45,
1162-0.437905770101505334663665494778098791E-45, -0.175162308040602133865466197911239516E-45,
1163-0.218952885050752667331832747389049396E-45, -0.350324616081204267730932395822479033E-45
1164
1165
1166coef
1167+0.0000000000000000, +0.90949470177292824E-12,
1168-0.46520653995685279E-9, -0.46520653995685279E-9,
1169+0.15848036127863449E-6, +0.0000000000000000,
1170-0.11546426321729090E-4, +0.11546426321729090E-4,
1171+0.0000000000000000, -0.78207794285845011E-3,
1172+0.12715073651634160E-1, +0.12715073651634160E-1,
1173-0.20021195337176320, +0.0000000000000000,
1174+0.75670659542083740, -0.75670659542083740,
1175+0.0000000000000000, +2.6588592529296879,
1176-1.9980468750000000, -1.9980468750000000,
1177+1.0000000000000000, +0.0000000000000000
1178call setResized(root, size(coef, 1, IK) - 1_IK)
1179call setPolyRoot(root, count, coef, jenkins)
1180count
1181+10
1182root(1:count)
1183+0.19531249999999967E-2, +0.19531249999999970E-2,
1184+0.39062500000000278E-2, +0.39062500000000243E-2,
1185+0.78124999999999133E-2, +0.78124999999999055E-2,
1186+0.15625000000000160E-1, +0.15625000000000177E-1,
1187+0.31249999999999677E-1, +0.31249999999999771E-1,
1188+0.62500000000000680E-1, +0.62500000000000319E-1,
1189+0.12499999999999963, +0.12499999999999994,
1190+0.24999999999999850, +0.24999999999999886,
1191+0.50000000000000322, +0.50000000000000255,
1192+0.99999999999999811, +0.99999999999999833
1193getPolyVal(coef, root(1:count))
1194-0.50487097934144756E-28, +0.20194839173657902E-27,
1195+0.40389678347315804E-27, +0.0000000000000000,
1196-0.50487097934144756E-27, -0.16155871338926322E-26,
1197-0.75730646901217133E-26, -0.30696155543960011E-25,
1198+0.16620352639920454E-24, -0.70601157751108026E-24,
1199+0.91577537200786489E-23, -0.11752588605501953E-22,
1200-0.99822686138607538E-22, -0.69420325114945901E-21,
1201-0.55985258009756409E-20, -0.10546344237752117E-18,
1202-0.20718055003338878E-16, -0.22954629892174609E-16,
1203-0.67286261814159909E-15, -0.92999469475984681E-14
1204
1205
1206coef
1207+0.00000000000000000000000000000000000, +0.909494701772928199999999999999999977E-12,
1208-0.465206539956852800000000000000000002E-9, -0.465206539956852800000000000000000002E-9,
1209+0.158480361278634500000000000000000004E-6, +0.00000000000000000000000000000000000,
1210-0.115464263217290900000000000000000001E-4, +0.115464263217290900000000000000000001E-4,
1211+0.00000000000000000000000000000000000, -0.782077942858450099999999999999999970E-3,
1212+0.127150736516341599999999999999999998E-1, +0.127150736516341599999999999999999998E-1,
1213-0.200211953371763200000000000000000000, +0.00000000000000000000000000000000000,
1214+0.756706595420837399999999999999999953, -0.756706595420837399999999999999999953,
1215+0.00000000000000000000000000000000000, +2.65885925292968799999999999999999995,
1216-1.99804687500000000000000000000000000, -1.99804687500000000000000000000000000,
1217+1.00000000000000000000000000000000000, +0.00000000000000000000000000000000000
1218call setResized(root, size(coef, 1, IK) - 1_IK)
1219call setPolyRoot(root, count, coef, jenkins)
1220count
1221+10
1222root(1:count)
1223+0.195312499999999767310055151845923262E-2, +0.195312499999999767310055151845923243E-2,
1224+0.390625000000002310324757362451074004E-2, +0.390625000000002310324757362451074154E-2,
1225+0.781249999999989766148737734007059782E-2, +0.781249999999989766148737734007059782E-2,
1226+0.156250000000002691698099555840414927E-1, +0.156250000000002691698099555840415017E-1,
1227+0.312499999999994379028135239093267455E-1, +0.312499999999994379028135239093267033E-1,
1228+0.625000000000009555071172310128835564E-1, +0.625000000000009555071172310128836888E-1,
1229+0.124999999999999149032692721659056846, +0.124999999999999149032692721659056677,
1230+0.249999999999999617367224612848344578, +0.249999999999999617367224612848344265,
1231+0.500000000000001492398035863052858631, +0.500000000000001492398035863052859882,
1232+0.999999999999999160184470589450447692, +0.999999999999999160184470589450446440
1233getPolyVal(coef, root(1:count))
1234-0.437905770101505334663665494778098791E-46, +0.00000000000000000000000000000000000,
1235+0.00000000000000000000000000000000000, -0.262743462060903200798199296866859275E-45,
1236+0.00000000000000000000000000000000000, -0.122613615628421493705826338537867661E-44,
1237+0.726923578368498855541684721331643993E-44, -0.262743462060903200798199296866859275E-44,
1238-0.160361093011171253553834304187739777E-42, -0.553512893408302743014873185399516872E-43,
1239+0.191233449803327379647622721569595742E-41, +0.643397431779337728007937164857669430E-41,
1240+0.146081949240855387902188432467114488E-39, +0.177360857749973751548678596894322439E-39,
1241+0.467314983435330978189620967762502824E-37, -0.114762799293412392402668283566530803E-36,
1242+0.175314932323943594248143528239381890E-34, -0.550582742852535418458451601761261904E-34,
1243+0.612482496979459595638357289641034424E-32, -0.124986688928018887944018337896766708E-31
1244
1245
1246coef
1247+0.0000000000000000, +0.90949470177292824E-12,
1248-0.46520653995685279E-9, -0.46520653995685279E-9,
1249+0.15848036127863449E-6, +0.0000000000000000,
1250-0.11546426321729090E-4, +0.11546426321729090E-4,
1251+0.0000000000000000, -0.78207794285845011E-3,
1252+0.12715073651634160E-1, +0.12715073651634160E-1,
1253-0.20021195337176320, +0.0000000000000000,
1254+0.75670659542083740, -0.75670659542083740,
1255+0.0000000000000000, +2.6588592529296879,
1256-1.9980468750000000, -1.9980468750000000,
1257+1.0000000000000000, +0.0000000000000000
1258call setResized(root, size(coef, 1, IK) - 1_IK)
1259call setPolyRoot(root, count, coef, laguerre)
1260count
1261+10
1262root(1:count)
1263+0.19531249999999967E-2, +0.19531249999999967E-2,
1264+0.39062500000000286E-2, +0.39062500000000191E-2,
1265+0.78124999999998933E-2, +0.78124999999998829E-2,
1266+0.15625000000000246E-1, +0.15625000000000243E-1,
1267+0.31249999999999514E-1, +0.31249999999999493E-1,
1268+0.62500000000000985E-1, +0.62500000000000805E-1,
1269+0.12499999999999883, +0.12499999999999893,
1270+0.24999999999999958, +0.24999999999999978,
1271+0.78124999999996548E-2, +0.78124999999991500E-2,
1272+0.78124999999998829E-2, +0.78124999999999037E-2
1273getPolyVal(coef, root(1:count))
1274+0.0000000000000000, +0.30292258760486853E-27,
1275+0.10097419586828951E-27, -0.40389678347315804E-27,
1276+0.20194839173657902E-27, +0.10097419586828951E-26,
1277-0.49477355975461860E-26, -0.20598735957131060E-25,
1278+0.35239994358033039E-25, -0.27626539989564010E-24,
1279-0.74448274613689857E-24, -0.63977250502148234E-24,
1280-0.18988283361098331E-21, +0.73599849030326141E-21,
1281-0.29902669212774665E-19, -0.21417656476565278E-19,
1282-0.25748419946413825E-25, +0.49780278563066729E-25,
1283+0.75730646901217133E-27, +0.0000000000000000
1284
1285
1286coef
1287+0.00000000000000000000000000000000000, +0.909494701772928199999999999999999977E-12,
1288-0.465206539956852800000000000000000002E-9, -0.465206539956852800000000000000000002E-9,
1289+0.158480361278634500000000000000000004E-6, +0.00000000000000000000000000000000000,
1290-0.115464263217290900000000000000000001E-4, +0.115464263217290900000000000000000001E-4,
1291+0.00000000000000000000000000000000000, -0.782077942858450099999999999999999970E-3,
1292+0.127150736516341599999999999999999998E-1, +0.127150736516341599999999999999999998E-1,
1293-0.200211953371763200000000000000000000, +0.00000000000000000000000000000000000,
1294+0.756706595420837399999999999999999953, -0.756706595420837399999999999999999953,
1295+0.00000000000000000000000000000000000, +2.65885925292968799999999999999999995,
1296-1.99804687500000000000000000000000000, -1.99804687500000000000000000000000000,
1297+1.00000000000000000000000000000000000, +0.00000000000000000000000000000000000
1298call setResized(root, size(coef, 1, IK) - 1_IK)
1299call setPolyRoot(root, count, coef, laguerre)
1300count
1301+10
1302root(1:count)
1303+0.195312499999999767310055151845921456E-2, +0.195312499999999767310055151845921456E-2,
1304+0.390625000000002310324757362451074154E-2, +0.390625000000002310324757362451073703E-2,
1305+0.781249999999989766148737734007059330E-2, +0.781249999999989766148737734007060308E-2,
1306+0.156250000000002691698099555840415108E-1, +0.156250000000002691698099555840414867E-1,
1307+0.312499999999994379028135239093266973E-1, +0.312499999999994379028135239093266823E-1,
1308+0.625000000000009555071172310128833999E-1, +0.625000000000009555071172310128833999E-1,
1309+0.124999999999999149032692721659056738, +0.124999999999999149032692721659056966,
1310+0.249999999999999617367224612848345493, +0.249999999999999617367224612848345108,
1311+0.781249999999989766148737734007061437E-2, +0.781249999999989766148737734007059180E-2,
1312+0.781249999999989766148737734007060308E-2, +0.781249999999989766148737734007060308E-2
1313getPolyVal(coef, root(1:count))
1314+0.00000000000000000000000000000000000, +0.245227231256842987411652677075735323E-44,
1315+0.00000000000000000000000000000000000, -0.262743462060903200798199296866859275E-45,
1316+0.919602117213161202793697539034007461E-45, -0.875811540203010669327330989556197582E-45,
1317-0.700649232162408535461864791644958066E-44, +0.105973196364564290988607049736299907E-43,
1318+0.292083148657704058220664885016991894E-43, +0.185321721906957057629663237390091408E-42,
1319+0.00000000000000000000000000000000000, -0.165668510944801498209957929984450335E-41,
1320-0.351657820198278300722315925421330955E-39, +0.140130547081713869500908420193783258E-39,
1321+0.328615020642129004131902579109303757E-37, +0.110756044348755321795305062654822812E-37,
1322+0.875811540203010669327330989556197582E-46, -0.183920423442632240558739507806801492E-44,
1323+0.00000000000000000000000000000000000, -0.183920423442632240558739507806801492E-44
1324
1325
1326coef
1327+0.0000000000000000, +0.90949470177292824E-12,
1328-0.46520653995685279E-9, -0.46520653995685279E-9,
1329+0.15848036127863449E-6, +0.0000000000000000,
1330-0.11546426321729090E-4, +0.11546426321729090E-4,
1331+0.0000000000000000, -0.78207794285845011E-3,
1332+0.12715073651634160E-1, +0.12715073651634160E-1,
1333-0.20021195337176320, +0.0000000000000000,
1334+0.75670659542083740, -0.75670659542083740,
1335+0.0000000000000000, +2.6588592529296879,
1336-1.9980468750000000, -1.9980468750000000,
1337+1.0000000000000000, +0.0000000000000000
1338call setResized(root, size(coef, 1, IK) - 1_IK)
1339call cmplx_roots_gen(root, coef, size(coef) - 1, .true., .true.)
1340root(1:count)
1341+0.99999999999999989, +0.99999999999999967,
1342+0.50000000000000078, +0.50000000000000122,
1343+0.19531249999999978E-2, +0.19531249999999963E-2,
1344+0.15625000000000257E-1, +0.15625000000000239E-1,
1345+0.31249999999999434E-1, +0.31249999999999410E-1,
1346+0.62500000000000583E-1, +0.62500000000000819E-1,
1347+0.12499999999999913, +0.12499999999999915,
1348+0.24999999999999958, +0.24999999999999978,
1349+0.39062500000000373E-2, +0.39062500000000104E-2,
1350+0.78124999999998829E-2, +0.78124999999999037E-2
1351getPolyVal(coef, root(1:count))
1352+0.78388153247232572E-15, +0.46152081436467382E-14,
1353+0.22472179073438050E-17, +0.18711633712987991E-16,
1354-0.50487097934144756E-28, +0.0000000000000000,
1355+0.14035413225692242E-25, -0.17771458472818954E-25,
1356+0.88554369776489901E-25, -0.14297946134949795E-24,
1357-0.36204297928575204E-23, -0.86728756315191227E-23,
1358-0.11209574623671259E-21, +0.51330232469644973E-22,
1359-0.29902669212774665E-19, -0.21417656476565278E-19,
1360+0.50487097934144756E-27, -0.40389678347315804E-27,
1361+0.75730646901217133E-27, +0.0000000000000000
1362
1363
1364!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1365! Compute the roots of polynomials with complex coefficients with multiple zeros.
1366!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1367
1368
1369coef
1370+288.000000, +0.00000000,
1371-1344.00000, +504.000000,
1372+2204.00000, -2352.00000,
1373-920.000000, +4334.00000,
1374-1587.00000, -3836.00000,
1375+2374.00000, +1394.00000,
1376-1293.00000, +200.000000,
1377+284.000000, -334.000000,
1378+3.00000000, +100.000000,
1379-10.0000000, -10.0000000,
1380+1.00000000, +0.00000000
1381call setResized(root, size(coef, 1, IK) - 1_IK)
1382call setPolyRoot(root, count, coef, eigen)
1383count
1384+10
1385root(1:count)
1386-0.262260437E-5, +4.00000191,
1387+2.99939489, +0.243153726E-2,
1388+3.00061989, -0.242062425E-2,
1389+0.158268567E-1, +2.02321959,
1390+0.118745584E-1, +1.97478998,
1391-0.277059320E-1, +2.00198984,
1392+1.02006865, +0.411952510E-1,
1393+1.04016447, -0.220653433E-1,
1394+0.960722089, +0.191663913E-1,
1395+0.979039192, -0.383066721E-1
1396getPolyVal(coef, root(1:count))
1397+0.239379883, -0.129583254,
1398+0.220031738E-1, +0.146927238E-1,
1399-0.135498047E-1, +0.105561465E-1,
1400-0.311279297E-2, +0.171809196E-1,
1401-0.695800781E-2, +0.117588043E-1,
1402-0.726318359E-2, +0.152904987E-1,
1403+0.823974609E-3, +0.151634216E-3,
1404+0.732421875E-3, -0.166893005E-3,
1405+0.823974609E-3, -0.217437744E-3,
1406+0.457763672E-3, -0.325202942E-3
1407
1408
1409coef
1410+288.00000000000000, +0.0000000000000000,
1411-1344.0000000000000, +504.00000000000000,
1412+2204.0000000000000, -2352.0000000000000,
1413-920.00000000000000, +4334.0000000000000,
1414-1587.0000000000000, -3836.0000000000000,
1415+2374.0000000000000, +1394.0000000000000,
1416-1293.0000000000000, +200.00000000000000,
1417+284.00000000000000, -334.00000000000000,
1418+3.0000000000000000, +100.00000000000000,
1419-10.000000000000000, -10.000000000000000,
1420+1.0000000000000000, +0.0000000000000000
1421call setResized(root, size(coef, 1, IK) - 1_IK)
1422call setPolyRoot(root, count, coef, eigen)
1423count
1424+10
1425root(1:count)
1426-0.26645352591003757E-13, +3.9999999999999880,
1427+3.0000001970524699, -0.23536775645374987E-6,
1428+2.9999998029474919, +0.23536782983636444E-6,
1429+0.10185201113391936E-4, +1.9999738276389909,
1430+0.17573524658964788E-4, +2.0000219075497694,
1431-0.27758725751883115E-4, +2.0000042648112437,
1432+0.99972283047521127, +0.13862002582648678E-4,
1433+1.0000137952404726, +0.27729146568424598E-3,
1434+0.99998601603900350, -0.27723617160111893E-3,
1435+1.0002773582453568, -0.13917296733077833E-4
1436getPolyVal(coef, root(1:count))
1437+0.33077185435104184E-9, +0.20417161294972133E-8,
1438-0.56900262279668823E-10, +0.24780200131034645E-9,
1439-0.16075318853836507E-9, +0.27178263366041855E-9,
1440+0.40927261579781771E-11, +0.18487610964693091E-10,
1441-0.34674485505092889E-11, +0.80048992608106051E-11,
1442-0.26147972675971687E-11, +0.25204519894794775E-10,
1443+0.22737367544323206E-12, +0.10120151244796816E-11,
1444+0.11368683772161603E-12, +0.10641210135275969E-11,
1445+0.56843418860808015E-13, +0.13136297605242930E-11,
1446+0.0000000000000000, +0.14300774106579262E-11
1447
1448
1449coef
1450+288.000000000000000000000000000000000, +0.00000000000000000000000000000000000,
1451-1344.00000000000000000000000000000000, +504.000000000000000000000000000000000,
1452+2204.00000000000000000000000000000000, -2352.00000000000000000000000000000000,
1453-920.000000000000000000000000000000000, +4334.00000000000000000000000000000000,
1454-1587.00000000000000000000000000000000, -3836.00000000000000000000000000000000,
1455+2374.00000000000000000000000000000000, +1394.00000000000000000000000000000000,
1456-1293.00000000000000000000000000000000, +200.000000000000000000000000000000000,
1457+284.000000000000000000000000000000000, -334.000000000000000000000000000000000,
1458+3.00000000000000000000000000000000000, +100.000000000000000000000000000000000,
1459-10.0000000000000000000000000000000000, -10.0000000000000000000000000000000000,
1460+1.00000000000000000000000000000000000, +0.00000000000000000000000000000000000
1461call setResized(root, size(coef, 1, IK) - 1_IK)
1462call setPolyRoot(root, count, coef, eigen)
1463count
1464+10
1465root(1:count)
1466-0.130963236218332038007806500095775058E-31, +4.00000000000000000000000000000001233,
1467+3.00000000000000017551097726216032856, -0.110201355633942977982250385201797701E-15,
1468+2.99999999999999982448902273783963986, +0.110201355633943002360704498217562890E-15,
1469+0.198827846113882935119497667613680559E-10, +2.00000000001021969660975542649721576,
1470-0.109087542230478321562619250537222123E-11, +1.99999999997767115512427343820588314,
1471-0.187919091890835102963115145452168761E-10, +2.00000000001210914826597113529687375,
1472+0.999999998924471798982144453844539050, +0.569327565931999073274615515126895450E-8,
1473+1.00000000569327572067379369961664497, +0.107552815246768798396814854480103414E-8,
1474+0.999999994306724389230175509482966430, -0.107552813966405382053522281260010247E-8,
1475+1.00000000107552809111388633705588402, -0.569327567212362489617908738729673540E-8
1476getPolyVal(coef, root(1:count))
1477+0.119226465062840671740415126035427104E-26, -0.104092661634241323385969495840791379E-28,
1478-0.109996792471754833617097901821618505E-27, +0.130434128372842219928401259706545043E-28,
1479-0.143769899976529401536287531022787790E-27, +0.135376768311312987970879984352679944E-27,
1480+0.833234331139693719466138297079942955E-29, +0.739591563885907217517079181128459199E-29,
1481+0.493038065763132378382330353301741394E-29, +0.103976104123069846169611920461958392E-28,
1482+0.576854536942864882707326513363037430E-29, +0.541027064395165263276723173580121787E-29,
1483-0.345126646034192664867631247311218975E-30, +0.283326920330550838279326406730858078E-30,
1484-0.295822839457879427029398211981044836E-30, +0.781981766077889895393086715455354831E-32,
1485+0.147911419728939713514699105990522418E-30, -0.297212059298750093875342660053652258E-30,
1486-0.542341872339445616220563388631915533E-30, +0.206541271182968673360478096961938272E-30
1487
1488
1489coef
1490+288.000000, +0.00000000,
1491-1344.00000, +504.000000,
1492+2204.00000, -2352.00000,
1493-920.000000, +4334.00000,
1494-1587.00000, -3836.00000,
1495+2374.00000, +1394.00000,
1496-1293.00000, +200.000000,
1497+284.000000, -334.000000,
1498+3.00000000, +100.000000,
1499-10.0000000, -10.0000000,
1500+1.00000000, +0.00000000
1501call setResized(root, size(coef, 1, IK) - 1_IK)
1502call setPolyRoot(root, count, coef, jenkins)
1503count
1504+10
1505root(1:count)
1506+1.00033581, +0.582871586E-2,
1507+0.994213104, +0.120317861E-1,
1508+0.990250826, -0.132678589E-1,
1509+1.01519012, -0.458721817E-2,
1510+0.607423335E-1, +2.00579453,
1511-0.251659974E-1, +1.94538498,
1512-0.355712250E-1, +2.04882479,
1513+3.00005126, -0.752519816E-3,
1514+0.155940652E-4, +4.00000763,
1515+2.99993849, +0.734806061E-3
1516getPolyVal(coef, root(1:count))
1517-0.122070312E-3, +0.337362289E-4,
1518-0.122070312E-3, +0.119447708E-3,
1519+0.00000000, +0.147819519E-3,
1520-0.305175781E-4, +0.857114792E-4,
1521+0.108489990, -0.866241455E-1,
1522+0.114166260, -0.640816689E-1,
1523+0.140686035, -0.817871094E-1,
1524-0.221252441E-1, +0.344342925E-1,
1525-0.312194824, -1.24936664,
1526-0.167846680E-1, -0.233982503E-1
1527
1528
1529coef
1530+288.00000000000000, +0.0000000000000000,
1531-1344.0000000000000, +504.00000000000000,
1532+2204.0000000000000, -2352.0000000000000,
1533-920.00000000000000, +4334.0000000000000,
1534-1587.0000000000000, -3836.0000000000000,
1535+2374.0000000000000, +1394.0000000000000,
1536-1293.0000000000000, +200.00000000000000,
1537+284.00000000000000, -334.00000000000000,
1538+3.0000000000000000, +100.00000000000000,
1539-10.000000000000000, -10.000000000000000,
1540+1.0000000000000000, +0.0000000000000000
1541call setResized(root, size(coef, 1, IK) - 1_IK)
1542call setPolyRoot(root, count, coef, jenkins)
1543count
1544+10
1545root(1:count)
1546+0.99999996474291075, +0.58311142231599888E-7,
1547+0.99999885335962246, -0.37544756228448384E-6,
1548+1.0000008584351527, -0.90929279042689637E-6,
1549+1.0000003234623087, +0.12264294506366602E-5,
1550+0.30410984488284110E-4, +1.9999860648750201,
1551-0.31429571698970132E-5, +2.0000332974285207,
1552-0.27268027223764306E-4, +1.9999806376961058,
1553+2.9999999877058445, +0.89451929269704773E-8,
1554+0.28471985494504760E-14, +4.0000000000001199,
1555+3.0000000122940649, -0.89452001361678413E-8
1556getPolyVal(coef, root(1:count))
1557-0.45474735088646412E-12, +0.50809761264905898E-12,
1558+0.22737367544323206E-12, +0.43439520914807100E-12,
1559+0.56843418860808015E-13, +0.32898705461248401E-13,
1560-0.22737367544323206E-12, +0.10524220384056093E-12,
1561-0.21032064978498966E-10, -0.16478911117612238E-10,
1562-0.18133050616597757E-10, -0.20641887287994876E-10,
1563-0.21088908397359774E-10, -0.21719921868001979E-10,
1564-0.47350567911053076E-10, +0.18249660848882264E-10,
1565+0.45945398596813902E-8, -0.54895050216060023E-8,
1566-0.26147972675971687E-11, +0.23185076734232421E-10
1567
1568
1569coef
1570+288.000000000000000000000000000000000, +0.00000000000000000000000000000000000,
1571-1344.00000000000000000000000000000000, +504.000000000000000000000000000000000,
1572+2204.00000000000000000000000000000000, -2352.00000000000000000000000000000000,
1573-920.000000000000000000000000000000000, +4334.00000000000000000000000000000000,
1574-1587.00000000000000000000000000000000, -3836.00000000000000000000000000000000,
1575+2374.00000000000000000000000000000000, +1394.00000000000000000000000000000000,
1576-1293.00000000000000000000000000000000, +200.000000000000000000000000000000000,
1577+284.000000000000000000000000000000000, -334.000000000000000000000000000000000,
1578+3.00000000000000000000000000000000000, +100.000000000000000000000000000000000,
1579-10.0000000000000000000000000000000000, -10.0000000000000000000000000000000000,
1580+1.00000000000000000000000000000000000, +0.00000000000000000000000000000000000
1581call setResized(root, size(coef, 1, IK) - 1_IK)
1582call setPolyRoot(root, count, coef, jenkins)
1583count
1584+10
1585root(1:count)
1586+1.00000000000063932099124545386125310, -0.711394493616352148533958438749886318E-12,
1587+0.999999999999786843588066074367120034, +0.237128343061650265624618815848318051E-12,
1588+0.999999999999206066535495788955911578, -0.284870330352794067539200364609923454E-12,
1589+1.00000000000036776888519268281573128, +0.759136480907495950452407887491330271E-12,
1590-0.909967953032966924767266632214659318E-13, +2.00000000000007857550779560734224166,
1591+0.113546837508897304279729212539452490E-12, +2.00000000000003951774411646449286381,
1592-0.225500422056006117836351268266463266E-13, +1.99999999999988190674808792816488144,
1593+2.99999999999999999644516009118258827, -0.523331030277165134427624041854420804E-18,
1594-0.175323129324981850099595482517513651E-31, +4.00000000000000000000000000000000077,
1595+3.00000000000000000355483990881739401, +0.523331030277173745165022263987145683E-18
1596getPolyVal(coef, root(1:count))
1597+0.295822839457879427029398211981044836E-30, -0.926844994206600978781222106602910152E-31,
1598-0.246519032881566189191165176650870697E-30, +0.325191254933908150503187817840604140E-30,
1599-0.147911419728939713514699105990522418E-30, +0.922898632587447611097423109983880845E-31,
1600+0.147911419728939713514699105990522418E-30, -0.100369871677403908759506693186486432E-30,
1601-0.404291213925768550273510889707427943E-29, +0.462437080235304506514587087558793324E-29,
1602+0.162702561701833684866169016589574660E-29, +0.196899554385239260912931279623930771E-29,
1603-0.157772181044202361082345713056557246E-29, -0.348877994637031041982717205382988168E-29,
1604+0.230741814777145953082930605345214972E-28, +0.267589680521319767083594237224497672E-28,
1605+0.828846292354401841298535556935557457E-27, +0.102031254573007041868747140456708580E-26,
1606+0.531988072958419836274534451212578964E-28, +0.735860503780615524865880233954677569E-28
1607
1608
1609coef
1610+288.000000, +0.00000000,
1611-1344.00000, +504.000000,
1612+2204.00000, -2352.00000,
1613-920.000000, +4334.00000,
1614-1587.00000, -3836.00000,
1615+2374.00000, +1394.00000,
1616-1293.00000, +200.000000,
1617+284.000000, -334.000000,
1618+3.00000000, +100.000000,
1619-10.0000000, -10.0000000,
1620+1.00000000, +0.00000000
1621call setResized(root, size(coef, 1, IK) - 1_IK)
1622call setPolyRoot(root, count, coef, laguerre)
1623count
1624+10
1625root(1:count)
1626+0.957954764, -0.170177817E-1,
1627+0.951423347, +0.225301757E-1,
1628+0.997259438, -0.460497029E-1,
1629+1.04122233, -0.311625414E-1,
1630-0.289512724E-1, +1.96756065,
1631-0.625183806E-1, +1.99496257,
1632+1.00274229, -0.732975230E-1,
1633+1.00888300, +0.437541790E-1,
1634+0.228498206E-1, +1.97342443,
1635+0.995942533, +0.540917069E-1
1636getPolyVal(coef, root(1:count))
1637-0.823974609E-3, -0.505447388E-4,
1638+0.143432617E-2, -0.185966492E-3,
1639+0.274658203E-3, +0.950813293E-3,
1640+0.946044922E-3, -0.114631653E-2,
1641+0.537109375E-2, -0.559072495E-1,
1642-0.148315430, +0.914115906E-1,
1643-0.143432617E-2, +0.575637817E-2,
1644+0.457763672E-3, +0.347137451E-3,
1645-0.262145996E-1, +0.920677185E-2,
1646-0.488281250E-3, +0.135707855E-2
1647
1648
1649coef
1650+288.00000000000000, +0.0000000000000000,
1651-1344.0000000000000, +504.00000000000000,
1652+2204.0000000000000, -2352.0000000000000,
1653-920.00000000000000, +4334.0000000000000,
1654-1587.0000000000000, -3836.0000000000000,
1655+2374.0000000000000, +1394.0000000000000,
1656-1293.0000000000000, +200.00000000000000,
1657+284.00000000000000, -334.00000000000000,
1658+3.0000000000000000, +100.00000000000000,
1659-10.000000000000000, -10.000000000000000,
1660+1.0000000000000000, +0.0000000000000000
1661call setResized(root, size(coef, 1, IK) - 1_IK)
1662call setPolyRoot(root, count, coef, laguerre)
1663count
1664+10
1665root(1:count)
1666+0.99968405225123580, -0.12620620904501057E-3,
1667+0.99946334030583894, +0.17188842164129054E-3,
1668+0.99989518445399872, -0.44586336082322157E-3,
1669+1.0004086013708342, -0.31976299826058218E-3,
1670+0.99962029840094790, -0.68455853304231278E-5,
1671-0.75712482998442344E-4, +1.9999952681228135,
1672+0.99986307467443547, -0.52164639372122117E-3,
1673-0.70846106203818572E-4, +1.9999952685239863,
1674-0.34064736976397726E-4, +1.9999630961850055,
1675+0.99999854691237156, +0.38853542066632834E-3
1676getPolyVal(coef, root(1:count))
1677-0.23874235921539366E-11, -0.31107061371216105E-12,
1678+0.17394086171407253E-10, +0.73605496697659589E-11,
1679+0.62527760746888816E-11, +0.50842663412709044E-11,
1680+0.68212102632969618E-11, -0.11450271286683744E-10,
1681-0.62527760746888816E-12, +0.38818479680680085E-11,
1682-0.21242385628283955E-9, +0.18444223925939696E-9,
1683+0.12732925824820995E-10, +0.87824469918729164E-11,
1684-0.17075763025786728E-9, +0.13913402201193303E-9,
1685+0.12448708730516955E-10, -0.80609438395107258E-10,
1686-0.56843418860808015E-12, +0.46141285237055740E-11
1687
1688
1689coef
1690+288.000000000000000000000000000000000, +0.00000000000000000000000000000000000,
1691-1344.00000000000000000000000000000000, +504.000000000000000000000000000000000,
1692+2204.00000000000000000000000000000000, -2352.00000000000000000000000000000000,
1693-920.000000000000000000000000000000000, +4334.00000000000000000000000000000000,
1694-1587.00000000000000000000000000000000, -3836.00000000000000000000000000000000,
1695+2374.00000000000000000000000000000000, +1394.00000000000000000000000000000000,
1696-1293.00000000000000000000000000000000, +200.000000000000000000000000000000000,
1697+284.000000000000000000000000000000000, -334.000000000000000000000000000000000,
1698+3.00000000000000000000000000000000000, +100.000000000000000000000000000000000,
1699-10.0000000000000000000000000000000000, -10.0000000000000000000000000000000000,
1700+1.00000000000000000000000000000000000, +0.00000000000000000000000000000000000
1701call setResized(root, size(coef, 1, IK) - 1_IK)
1702call setPolyRoot(root, count, coef, laguerre)
1703count
1704+10
1705root(1:count)
1706+0.999999987186163931177263416074192316, -0.523210113724573224158980506799552893E-8,
1707+0.999999983645692704429053992803497744, +0.516403242552248341950802865482184609E-8,
1708+0.999999996834913808707740552774975021, -0.136112912657541668901042030848548143E-7,
1709+1.00000001250310236797332654499325579, -0.970719584308792676978317369510932310E-8,
1710+0.999999984385238514809745194471578835, -0.280548424142042709182285719596106521E-9,
1711-0.771195183876041062811447140946995581E-10, +1.99999999999511724532936826592821654,
1712+0.999999997457890361873775552226706741, -0.103669219710888875942277199987969308E-7,
1713-0.719073015128987083181422783095156867E-10, +1.99999999999476261447634037342730026,
1714-0.346099727524717880518841359170124199E-10, +1.99999999996294152658095642717453627,
1715+0.999999999809654610234193591766066401, +0.162757841295939920142426163656213815E-7
1716getPolyVal(coef, root(1:count))
1717-0.690253292068385329735262494622437951E-29, -0.103992693465723488793006962745517685E-30,
1718+0.143967115202834654487640463164108487E-28, +0.670912877140922800002062790107526930E-29,
1719+0.483177304447869730814683746235706566E-29, +0.464961255579136580109909497254227733E-29,
1720+0.621227962861546796761736245160194156E-29, -0.955002328589667861027763317238520064E-29,
1721-0.167632942359465008649992320122592074E-29, +0.111187254067668374511683156036595891E-28,
1722-0.226649598831311954342357263412810519E-27, +0.191225446604550158529689992651317627E-27,
1723+0.182424084332358980001462230721644316E-29, +0.189803913024991356676497499453417923E-29,
1724-0.192383453260774254044785303858339492E-27, +0.150775130007754964478287418087745154E-27,
1725+0.202145606962884275136755444853713971E-29, -0.784384562543142124431676932654714688E-28,
1726-0.128189897098414418379405891858452762E-29, +0.125831339954351554987430471299399555E-28
1727
1728
1729coef
1730+288.00000000000000, +0.0000000000000000,
1731-1344.0000000000000, +504.00000000000000,
1732+2204.0000000000000, -2352.0000000000000,
1733-920.00000000000000, +4334.0000000000000,
1734-1587.0000000000000, -3836.0000000000000,
1735+2374.0000000000000, +1394.0000000000000,
1736-1293.0000000000000, +200.00000000000000,
1737+284.00000000000000, -334.00000000000000,
1738+3.0000000000000000, +100.00000000000000,
1739-10.000000000000000, -10.000000000000000,
1740+1.0000000000000000, +0.0000000000000000
1741call setResized(root, size(coef, 1, IK) - 1_IK)
1742call cmplx_roots_gen(root, coef, size(coef) - 1, .true., .true.)
1743root(1:count)
1744+0.79936057773011271E-14, +3.9999999999999916,
1745+3.0000000712765424, -0.55974001247601282E-7,
1746+2.9999999287234407, +0.55974006011966191E-7,
1747+1.0001005681730508, -0.60420290115003961E-4,
1748+0.99989854158045077, +0.61709313521636971E-4,
1749+0.22930793966121343E-4, +1.9999959028873173,
1750+0.99993937508641983, -0.10165852505892826E-3,
1751-0.79168180886902579E-5, +2.0000219071663787,
1752-0.15013975879406113E-4, +1.9999821899463428,
1753+1.0000615151600896, +0.10036950161902832E-3
1754getPolyVal(coef, root(1:count))
1755-0.21429968910524622E-10, -0.62879479401089932E-9,
1756-0.14665602066088468E-10, -0.12111342185633324E-11,
1757-0.17735146684572101E-10, -0.23460623147713926E-10,
1758-0.34106051316484809E-12, -0.91614216213287136E-13,
1759+0.0000000000000000, -0.61336352663587945E-13,
1760+0.13073986337985843E-11, -0.12248734578695908E-10,
1761+0.11368683772161603E-12, -0.56958604499612875E-12,
1762+0.17053025658242404E-11, -0.84408322345574938E-11,
1763+0.47180037654470652E-11, -0.33409741986845098E-11,
1764+0.11368683772161603E-12, -0.65083008427002653E-12
1765
1766
1767!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1768! Compute the roots of polynomials with complex coefficients with 12 zeros evenly distribute on a circle of radius 1. centered at 0+2i.
1769!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1770
1771
1772coef
1773+4095.00000, +0.00000000,
1774+0.00000000, +24576.0000,
1775-67584.0000, +0.00000000,
1776+0.00000000, -112640.000,
1777+126720.000, +0.00000000,
1778+0.00000000, +101376.000,
1779-59136.0000, +0.00000000,
1780+0.00000000, -25344.0000,
1781+7920.00000, +0.00000000,
1782+0.00000000, +1760.00000,
1783-264.000000, +0.00000000,
1784+0.00000000, -24.0000000,
1785+1.00000000, +0.00000000
1786call setResized(root, size(coef, 1, IK) - 1_IK)
1787call setPolyRoot(root, count, coef, eigen)
1788count
1789+12
1790root(1:count)
1791-0.102758408E-2, +3.26374292,
1792+0.743275821, +2.94730639,
1793-0.743724287, +2.94551206,
1794+0.972983360, +2.28771806,
1795-0.971354008, +2.28610396,
1796+0.878654838, +1.98345017,
1797+0.860098600, +1.51398087,
1798-0.878169596, +1.98575449,
1799-0.860523164, +1.51400054,
1800+0.499328941, +1.13599277,
1801-0.566091039E-4, +1.00045466,
1802-0.499485791, +1.13598478
1803getPolyVal(coef, root(1:count))
1804+12.4841309, +0.153654099,
1805-2.42480469, -4.67864990,
1806-2.90429688, +5.80480957,
1807-1.59033203, -0.416015625,
1808-1.83007812, +0.395629883,
1809-0.713867188, -0.938720703E-1,
1810-0.152343750, +0.140258789,
1811-0.572753906, -0.150268555,
1812-0.998535156E-1, -0.616455078E-1,
1813-0.258789062E-1, -0.122070312E-2,
1814-0.170898438E-2, -0.675007701E-3,
1815-0.190429688E-1, -0.146484375E-1
1816
1817
1818coef
1819+4095.0000000000000, +0.0000000000000000,
1820+0.0000000000000000, +24576.000000000000,
1821-67584.000000000000, +0.0000000000000000,
1822+0.0000000000000000, -112640.00000000000,
1823+126720.00000000000, +0.0000000000000000,
1824+0.0000000000000000, +101376.00000000000,
1825-59136.000000000000, +0.0000000000000000,
1826+0.0000000000000000, -25344.000000000000,
1827+7920.0000000000000, +0.0000000000000000,
1828+0.0000000000000000, +1760.0000000000000,
1829-264.00000000000000, +0.0000000000000000,
1830+0.0000000000000000, -24.000000000000000,
1831+1.0000000000000000, +0.0000000000000000
1832call setResized(root, size(coef, 1, IK) - 1_IK)
1833call setPolyRoot(root, count, coef, eigen)
1834count
1835+12
1836root(1:count)
1837-0.22127855103803995E-10, +3.0000000012083348,
1838+0.50000000092178387, +2.8660254036044090,
1839+0.86602540366466041, +2.4999999995817723,
1840-0.50000000092225982, +2.8660254035707347,
1841-0.86602540364973568, +2.4999999995766795,
1842+0.99999999988689936, +2.0000000000036473,
1843-0.99999999988111377, +2.0000000000065561,
1844+0.86602540376976866, +1.5000000000087796,
1845-0.86602540376847770, +1.5000000000108240,
1846+0.49999999999651829, +1.1339745962147472,
1847+0.18475302723801401E-12, +0.99999999999813760,
1848-0.49999999999610956, +1.1339745962153789
1849getPolyVal(coef, root(1:count))
1850+0.14059878594707698E-7, +0.26553426475682450E-9,
1851+0.46070454118307680E-8, -0.11811152944574133E-7,
1852-0.32819116313476115E-8, -0.20859260985162109E-8,
1853+0.47079993237275630E-8, +0.13283624866744503E-7,
1854-0.23110260372050107E-8, +0.18424088921165094E-8,
1855-0.16898411558941007E-8, +0.26602720026858151E-10,
1856-0.92268237494863570E-9, -0.49340087571181357E-10,
1857-0.29012880986556411E-9, +0.29558577807620168E-10,
1858-0.10641088010743260E-9, -0.33196556614711881E-10,
1859-0.22737367544323206E-11, -0.93223206931725144E-11,
1860+0.19554136088117957E-10, +0.22170363269036277E-11,
1861-0.90949470177292824E-12, +0.60936145018786192E-10
1862
1863
1864coef
1865+4095.00000000000000000000000000000000, +0.00000000000000000000000000000000000,
1866+0.00000000000000000000000000000000000, +24576.0000000000000000000000000000000,
1867-67584.0000000000000000000000000000000, +0.00000000000000000000000000000000000,
1868+0.00000000000000000000000000000000000, -112640.000000000000000000000000000000,
1869+126720.000000000000000000000000000000, +0.00000000000000000000000000000000000,
1870+0.00000000000000000000000000000000000, +101376.000000000000000000000000000000,
1871-59136.0000000000000000000000000000000, +0.00000000000000000000000000000000000,
1872+0.00000000000000000000000000000000000, -25344.0000000000000000000000000000000,
1873+7920.00000000000000000000000000000000, +0.00000000000000000000000000000000000,
1874+0.00000000000000000000000000000000000, +1760.00000000000000000000000000000000,
1875-264.000000000000000000000000000000000, +0.00000000000000000000000000000000000,
1876+0.00000000000000000000000000000000000, -24.0000000000000000000000000000000000,
1877+1.00000000000000000000000000000000000, +0.00000000000000000000000000000000000
1878call setResized(root, size(coef, 1, IK) - 1_IK)
1879call setPolyRoot(root, count, coef, eigen)
1880count
1881+12
1882root(1:count)
1883+0.866025403784438646763723170731842991, +2.49999999999999999999999999930164047,
1884+0.500000000000000000000000001261300765, +2.86602540378443864676372317070946253,
1885-0.347698086029325849838930767012779505E-28, +3.00000000000000000000000000152269067,
1886-0.500000000000000000000000001265360529, +2.86602540378443864676372317065155482,
1887-0.866025403784438646763723170698404803, +2.49999999999999999999999999929392211,
1888+0.999999999999999999999999999757820376, +1.99999999999999999999999999995932648,
1889-0.999999999999999999999999999750830695, +1.99999999999999999999999999997135064,
1890+0.866025403784438646763723170716402811, +1.50000000000000000000000000003897716,
1891-0.866025403784438646763723170717795451, +1.50000000000000000000000000004211450,
1892+0.499999999999999999999999999996946149, +1.13397459621556135323627682925377416,
1893+0.400031557749517858232708991895941675E-31, +1.00000000000000000000000000000101670,
1894-0.499999999999999999999999999997192042, +1.13397459621556135323627682925417245
1895getPolyVal(coef, root(1:count))
1896-0.494852445845140705534777329001891802E-26, -0.688083924579027547270380241067910289E-26,
1897+0.716877347619594478167908333700731986E-26, -0.143695944266664931679530181469792529E-25,
1898+0.144759920412581771352079250372217687E-25, +0.417237703235191019806716927084944130E-27,
1899+0.751863328766146351737918495571023556E-26, +0.134045217167417378505074447134264243E-25,
1900-0.646037638330747618041935108538337783E-26, +0.660592122032075285851781500567805189E-26,
1901-0.272866987115947983491916910731315757E-26, -0.464639073175175953387508124951561089E-27,
1902-0.310732310566556550151679881864889496E-26, +0.336054745624151029105396368810466934E-27,
1903-0.737190515929035532157260344256763732E-27, +0.264268403249038954812929069369733387E-27,
1904-0.575474030358728112047855988373792555E-27, -0.318305375256678263483632476091604244E-27,
1905-0.686308987542280270708203851796024020E-28, -0.315544362088404722164691426113114492E-29,
1906-0.138050658413677065947052498924487590E-28, +0.480037869299421429879250790270093257E-30,
1907-0.958465999843529343575250206818585269E-28, -0.457539325028186847138802567864016013E-28
1908
1909
1910coef
1911+4095.00000, +0.00000000,
1912+0.00000000, +24576.0000,
1913-67584.0000, +0.00000000,
1914+0.00000000, -112640.000,
1915+126720.000, +0.00000000,
1916+0.00000000, +101376.000,
1917-59136.0000, +0.00000000,
1918+0.00000000, -25344.0000,
1919+7920.00000, +0.00000000,
1920+0.00000000, +1760.00000,
1921-264.000000, +0.00000000,
1922+0.00000000, -24.0000000,
1923+1.00000000, +0.00000000
1924call setResized(root, size(coef, 1, IK) - 1_IK)
1925call setPolyRoot(root, count, coef, jenkins)
1926count
1927+12
1928root(1:count)
1929+0.232841168E-2, +0.984195650,
1930-0.764802098E-1, +1.36573625,
1931-0.352237761, +2.45377731,
1932+0.665973306, +1.15707934,
1933-0.679132044, +1.19625020,
1934-1.05274248, +1.75213408,
1935+1.06248283, +1.70380104,
1936+0.549214125, +2.57492280,
1937+0.314331472, +3.00236869,
1938-0.975896001, +2.43641162,
1939+0.977696478, +2.40769362,
1940-0.435538173, +2.96562791
1941getPolyVal(coef, root(1:count))
1942+0.211914062, +0.331363678E-1,
1943-1.03515625, -0.161132812E-1,
1944-0.893066406, +0.442382812,
1945-1.37841797, +2.32958984,
1946-2.00976562, -1.55371094,
1947-3.45703125, +0.648193359,
1948-4.32031250, +0.306762695,
1949-0.815429688, -1.86846924,
1950-4.67919922, +2.79812622,
1951-0.570800781, +2.77392578,
1952-1.57812500, -2.06347656,
1953-0.435058594, -2.53369141
1954
1955
1956coef
1957+4095.0000000000000, +0.0000000000000000,
1958+0.0000000000000000, +24576.000000000000,
1959-67584.000000000000, +0.0000000000000000,
1960+0.0000000000000000, -112640.00000000000,
1961+126720.00000000000, +0.0000000000000000,
1962+0.0000000000000000, +101376.00000000000,
1963-59136.000000000000, +0.0000000000000000,
1964+0.0000000000000000, -25344.000000000000,
1965+7920.0000000000000, +0.0000000000000000,
1966+0.0000000000000000, +1760.0000000000000,
1967-264.00000000000000, +0.0000000000000000,
1968+0.0000000000000000, -24.000000000000000,
1969+1.0000000000000000, +0.0000000000000000
1970call setResized(root, size(coef, 1, IK) - 1_IK)
1971call setPolyRoot(root, count, coef, jenkins)
1972count
1973+12
1974root(1:count)
1975-0.48864020079139982E-11, +1.0000000000059499,
1976-0.49999999999999817, +1.1339745962218499,
1977+0.49999999999181693, +1.1339745962167307,
1978+0.86602540377731263, +1.4999999999873355,
1979-0.86602540378127224, +1.4999999999936720,
1980-0.86602540416987528, +2.5000000000049512,
1981-0.99999999998716704, +1.9999999999123383,
1982+1.0000000001130849, +2.0000000000422529,
1983+0.86602540353486657, +2.5000000003320020,
1984+0.83678087047241487E-9, +2.9999999995664282,
1985+0.49999999949958140, +2.8660254031889592,
1986-0.49999999981024412, +2.8660254045275293
1987getPolyVal(coef, root(1:count))
1988-0.74123818194493651E-10, -0.58636824091162948E-10,
1989-0.70940586738288403E-10, -0.56843418860808015E-10,
1990-0.65028871176764369E-10, -0.10049916454590857E-9,
1991+0.15461409930139780E-10, -0.29353941499721259E-9,
1992+0.86401996668428183E-11, -0.10254552762489766E-9,
1993+0.49394657253287733E-8, +0.22826043277746066E-8,
1994-0.71258909883908927E-9, +0.77216100180521607E-9,
1995+0.11218617146369070E-8, +0.69576344685629010E-9,
1996-0.18981154426001012E-8, +0.55988493841141462E-8,
1997-0.22200765670277178E-8, -0.10041370403022984E-7,
1998-0.91895344667136669E-8, -0.18036416804534383E-8,
1999+0.87352418631780893E-8, -0.38570533433812670E-8
2000
2001
2002coef
2003+4095.00000000000000000000000000000000, +0.00000000000000000000000000000000000,
2004+0.00000000000000000000000000000000000, +24576.0000000000000000000000000000000,
2005-67584.0000000000000000000000000000000, +0.00000000000000000000000000000000000,
2006+0.00000000000000000000000000000000000, -112640.000000000000000000000000000000,
2007+126720.000000000000000000000000000000, +0.00000000000000000000000000000000000,
2008+0.00000000000000000000000000000000000, +101376.000000000000000000000000000000,
2009-59136.0000000000000000000000000000000, +0.00000000000000000000000000000000000,
2010+0.00000000000000000000000000000000000, -25344.0000000000000000000000000000000,
2011+7920.00000000000000000000000000000000, +0.00000000000000000000000000000000000,
2012+0.00000000000000000000000000000000000, +1760.00000000000000000000000000000000,
2013-264.000000000000000000000000000000000, +0.00000000000000000000000000000000000,
2014+0.00000000000000000000000000000000000, -24.0000000000000000000000000000000000,
2015+1.00000000000000000000000000000000000, +0.00000000000000000000000000000000000
2016call setResized(root, size(coef, 1, IK) - 1_IK)
2017call setPolyRoot(root, count, coef, jenkins)
2018count
2019+12
2020root(1:count)
2021-0.132101502355678455782266806111711463E-29, +1.00000000000000000000000000000065771,
2022-0.500000000000000000000000000004124283, +1.13397459621556135323627682924768014,
2023+0.500000000000000000000000000001422492, +1.13397459621556135323627682924685237,
2024+0.866025403784438646763723170765615521, +1.50000000000000000000000000000095854,
2025-0.866025403784438646763723170767965348, +1.50000000000000000000000000000229629,
2026-0.866025403784438646763723170674364768, +2.50000000000000000000000000007722748,
2027-1.00000000000000000000000000002751056, +2.00000000000000000000000000004396821,
2028+1.00000000000000000000000000002423398, +2.00000000000000000000000000004125072,
2029+0.866025403784438646763723170671850659, +2.50000000000000000000000000007286602,
2030+0.121557848053039041451101210458919685E-28, +2.99999999999999999999999999993981739,
2031+0.499999999999999999999999999916534385, +2.86602540378443864676372317065553842,
2032-0.499999999999999999999999999916526778, +2.86602540378443864676372317067088808
2033getPolyVal(coef, root(1:count))
2034+0.473316543132607083247037139169671738E-29, -0.158521802826814146938720167332137067E-28,
2035+0.358931711875560371462336497203667735E-28, -0.603478592494074031139972352441331466E-28,
2036+0.299767143983984486056456854807458767E-28, +0.356959559612507841948807175790460769E-28,
2037+0.117934705330541264909053420509776541E-27, +0.154025091744402555006640002371464011E-27,
2038+0.165266359643801973233757134426743715E-27, -0.133711923434961501017287991815432266E-27,
2039-0.119827971503071693242041569066455228E-26, -0.486529963295059030987683592638158407E-27,
2040+0.546286176865550675247622031458329464E-27, -0.611564416772589402145442570235480025E-27,
2041+0.224430927535377858639636776822952682E-27, +0.827515089576841383876903264981642755E-27,
2042+0.805821414683263559228080729436366134E-27, +0.212243026549713226246025570489333635E-26,
2043-0.444523120092040152349509046536850040E-27, -0.145869417663646849741321452469565459E-27,
2044-0.358576724468210916149901219349290481E-26, +0.356989141896453629891510115611658873E-26,
2045-0.622016823766767808567147973725476942E-27, +0.272837404832002195549213970910117652E-26
2046
2047
2048coef
2049+4095.00000, +0.00000000,
2050+0.00000000, +24576.0000,
2051-67584.0000, +0.00000000,
2052+0.00000000, -112640.000,
2053+126720.000, +0.00000000,
2054+0.00000000, +101376.000,
2055-59136.0000, +0.00000000,
2056+0.00000000, -25344.0000,
2057+7920.00000, +0.00000000,
2058+0.00000000, +1760.00000,
2059-264.000000, +0.00000000,
2060+0.00000000, -24.0000000,
2061+1.00000000, +0.00000000
2062call setResized(root, size(coef, 1, IK) - 1_IK)
2063call setPolyRoot(root, count, coef, laguerre)
2064count
2065+12
2066root(1:count)
2067+0.00000000, +1.70626163,
2068+0.150868297E-2, +1.00111639,
2069+0.481668174, +1.13538647,
2070-0.491802990, +1.14549077,
2071-0.498571545, +1.13414574,
2072+0.502268970, +1.13399494,
2073+0.971725762, +1.93425488,
2074-0.366622210, +2.13671899,
2075+0.347302943, +1.49599338,
2076+0.130578488, +1.48512542,
2077+0.588253796, +1.22825062,
2078+0.138924181, +1.76272321
2079getPolyVal(coef, root(1:count))
2080-0.990966797, +0.00000000,
2081-0.153808594E-1, +0.178875923E-1,
2082-0.124267578, -0.143554688,
2083-0.158691406, +0.217285156E-1,
2084-0.634765625E-2, -0.561523438E-2,
2085+0.175781250E-1, +0.372314453E-1,
2086-0.376953125, -0.617187500,
2087-1.20263672, -0.155517578,
2088-1.01171875, +0.823974609E-2,
2089-1.00097656, +0.469970703E-2,
2090-0.968750000, +0.664550781,
2091-1.07080078, +0.395202637E-1
2092
2093
2094coef
2095+4095.0000000000000, +0.0000000000000000,
2096+0.0000000000000000, +24576.000000000000,
2097-67584.000000000000, +0.0000000000000000,
2098+0.0000000000000000, -112640.00000000000,
2099+126720.00000000000, +0.0000000000000000,
2100+0.0000000000000000, +101376.00000000000,
2101-59136.000000000000, +0.0000000000000000,
2102+0.0000000000000000, -25344.000000000000,
2103+7920.0000000000000, +0.0000000000000000,
2104+0.0000000000000000, +1760.0000000000000,
2105-264.00000000000000, +0.0000000000000000,
2106+0.0000000000000000, -24.000000000000000,
2107+1.0000000000000000, +0.0000000000000000
2108call setResized(root, size(coef, 1, IK) - 1_IK)
2109call setPolyRoot(root, count, coef, laguerre)
2110count
2111+12
2112root(1:count)
2113-0.99999999995068578, +2.0000000000185447,
2114+0.73453638394708858E-21, +1.0000000000003504,
2115-0.10311355583400225E-14, +1.0000000000012070,
2116-0.50000000000139555, +1.1339745962161143,
2117+0.49999999999915301, +1.1339745962166152,
2118+0.49999999999514877, +1.1339745962201375,
2119-0.99999999999267275, +2.0000000000129092,
2120-0.59777740889872463E-15, +2.9999999998944782,
2121-0.86602540377771808, +1.4999999999975995,
2122-0.86602540378634441, +1.4999999999945266,
2123+0.86602540374995229, +2.4999999999971356,
2124+0.50000000000276834, +1.1339745962161540
2125getPolyVal(coef, root(1:count))
2126-0.37744030123576522E-10, -0.67802830017171800E-9,
2127+0.95496943686157465E-11, +0.88144366072918776E-20,
2128-0.13187673175707459E-10, -0.12373626699900132E-13,
2129+0.63664629124104977E-11, -0.19326762412674725E-10,
2130+0.40927261579781771E-11, +0.16370904631912708E-10,
2131-0.72759576141834259E-10, -0.95496943686157465E-11,
2132+0.47020876081660390E-9, +0.47748471843078732E-11,
2133-0.22314452507998794E-8, +0.71733288970084749E-14,
2134-0.41836756281554699E-10, +0.10459189070388675E-10,
2135+0.17143975128419697E-9, +0.0000000000000000,
2136+0.10063558875117451E-8, +0.28076101443730295E-8,
2137+0.0000000000000000, +0.10004441719502211E-10
2138
2139
2140coef
2141+4095.00000000000000000000000000000000, +0.00000000000000000000000000000000000,
2142+0.00000000000000000000000000000000000, +24576.0000000000000000000000000000000,
2143-67584.0000000000000000000000000000000, +0.00000000000000000000000000000000000,
2144+0.00000000000000000000000000000000000, -112640.000000000000000000000000000000,
2145+126720.000000000000000000000000000000, +0.00000000000000000000000000000000000,
2146+0.00000000000000000000000000000000000, +101376.000000000000000000000000000000,
2147-59136.0000000000000000000000000000000, +0.00000000000000000000000000000000000,
2148+0.00000000000000000000000000000000000, -25344.0000000000000000000000000000000,
2149+7920.00000000000000000000000000000000, +0.00000000000000000000000000000000000,
2150+0.00000000000000000000000000000000000, +1760.00000000000000000000000000000000,
2151-264.000000000000000000000000000000000, +0.00000000000000000000000000000000000,
2152+0.00000000000000000000000000000000000, -24.0000000000000000000000000000000000,
2153+1.00000000000000000000000000000000000, +0.00000000000000000000000000000000000
2154call setResized(root, size(coef, 1, IK) - 1_IK)
2155call setPolyRoot(root, count, coef, laguerre)
2156count
2157+12
2158root(1:count)
2159-0.999999999999999999999999999952989687, +1.99999999999999999999999999997059375,
2160+0.176794967929835023895557031022669585E-52, +1.00000000000000000000000000000005913,
2161+0.502360025638320651109473759790750207E-44, +0.999999999999999999999999999999272865,
2162-0.499999999999999999999999999998095063, +1.13397459621556135323627682924526175,
2163+0.500000000000000000000000000000621016, +1.13397459621556135323627682924703110,
2164+0.499999999999999999999999999999574899, +1.13397459621556135323627682924621913,
2165-1.00000000000000000000000000001824068, +2.00000000000000000000000000000100842,
2166-0.499999999999999999999999999623053139, +2.86602540378443864676372317086667542,
2167-0.866025403784438646763723170753606674, +1.50000000000000000000000000000857694,
2168-0.866025403784438646763723170751245195, +1.50000000000000000000000000000166343,
2169-0.260808262955316027176952067679347635E-52, +0.999999999999999999999999999999835429,
2170+0.500000000000000000000000000000089652, +1.13397459621556135323627682924599052
2171getPolyVal(coef, root(1:count))
2172-0.206681557167905093017872884104089992E-27, +0.220092192556662293709872269713897358E-27,
2173+0.394430452610505902705864282641393115E-29, +0.212153961515802028674668437226824642E-51,
2174-0.749417859959961215141142137018646918E-29, +0.602832030765984781331368511756793015E-43,
2175+0.197215226305252951352932141320696557E-29, +0.201159530831358010379990784147110489E-28,
2176-0.315544362088404722164691426113114492E-29, -0.116356983520099241298229963379210969E-28,
2177-0.106496222204836593730583356313176141E-28, -0.130162049361466947892935213271659728E-28,
2178+0.240602576092408600650577212411249800E-28, +0.242771943581766383115459465965777462E-27,
2179-0.295822839457879427029398211981044836E-26, -0.745966593499619288492465824545534728E-26,
2180+0.157772181044202361082345713056557246E-28, -0.743501403170803626600554172779026021E-28,
2181+0.560091242706918381842327281350778223E-28, -0.997909045104579933845836635082724581E-28,
2182-0.749417859959961215141142137018646918E-29, -0.312969915546379232612342481216386687E-51,
2183+0.591645678915758854058796423962089672E-29, +0.136078506150624536433523177511280625E-28
2184
2185
2186coef
2187+4095.0000000000000, +0.0000000000000000,
2188+0.0000000000000000, +24576.000000000000,
2189-67584.000000000000, +0.0000000000000000,
2190+0.0000000000000000, -112640.00000000000,
2191+126720.00000000000, +0.0000000000000000,
2192+0.0000000000000000, +101376.00000000000,
2193-59136.000000000000, +0.0000000000000000,
2194+0.0000000000000000, -25344.000000000000,
2195+7920.0000000000000, +0.0000000000000000,
2196+0.0000000000000000, +1760.0000000000000,
2197-264.00000000000000, +0.0000000000000000,
2198+0.0000000000000000, -24.000000000000000,
2199+1.0000000000000000, +0.0000000000000000
2200call setResized(root, size(coef, 1, IK) - 1_IK)
2201call cmplx_roots_gen(root, coef, size(coef) - 1, .true., .true.)
2202root(1:count)
2203-0.49999999970081721, +2.8660254041280107,
2204+0.49999999986772581, +2.8660254035722428,
2205+0.99999999997852063, +1.9999999999519149,
2206-0.50000000000066047, +1.1339745962149845,
2207-0.60552205778295854E-23, +1.0000000000005400,
2208+0.86602540374488346, +1.4999999999921991,
2209-0.86602540388409677, +2.5000000002639675,
2210+0.38753530091931363E-9, +2.9999999997259361,
2211-1.0000000001682818, +2.0000000000982645,
2212-0.86602540379604520, +1.4999999999997076,
2213+0.86602540374995229, +2.4999999999971356,
2214+0.50000000000276834, +1.1339745962161540
2215getPolyVal(coef, root(1:count))
2216+0.25720510166138411E-8, -0.23226220946526155E-9,
2217-0.20977495296392590E-8, -0.32093794288812205E-9,
2218-0.28785507311113179E-9, -0.63982952269725502E-9,
2219+0.21373125491663814E-10, +0.50022208597511053E-11,
2220-0.14097167877480388E-10, -0.72662646933768450E-22,
2221-0.28740032576024532E-9, -0.37857716961298138E-9,
2222+0.17521415429655463E-8, -0.10079475032398477E-8,
2223-0.61882019508630037E-8, -0.46504235927026007E-8,
2224+0.26066118152812123E-8, -0.82241058407817036E-9,
2225+0.98225427791476250E-10, +0.49340087571181357E-10,
2226+0.10063558875117451E-8, +0.28076101443730295E-8,
2227+0.0000000000000000, +0.10004441719502211E-10
2228
2229
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#if 0
28 , bench_type(name = SK_"Skowron", exec = setPolyRootSkowron, overhead = setOverhead) &
29#endif
30 ]
31 allocate(rootCount(size(bench)))
32
33 write(*,"(*(g0,:,' '))")
34 write(*,"(*(g0,:,' '))") "Benchmarking setPolyRoot()"
35 write(*,"(*(g0,:,' '))")
36
37 open(newunit = fileUnit, file = "main.out", status = "replace")
38 open(newunit = rootUnit, file = "root.count", status = "replace")
39
40 write(fileUnit, "(*(g0,:,','))") "Polynomial Degree", (bench(ibench)%name, ibench = 1, size(bench))
41 write(rootUnit, "(*(g0,:,','))") "Polynomial Degree", (bench(ibench)%name, ibench = 1, size(bench))
42 loopOverArraySize: do ideg = 0, degmax
43
44 degree = 2**ideg
45 call setUnifRand(coef(0 : degree), (1._RKG, 1._RKG), (2._RKG, 2._RKG))
46 write(*,"(*(g0,:,' '))") "Benchmarking with coef size", degree
47 do ibench = 1, size(bench)
48 bench(ibench)%timing = bench(ibench)%getTiming(minsec = 0.07_RK)
49 end do
50 write(rootUnit,"(*(g0,:,','))") degree, rootCount
51 write(fileUnit,"(*(g0,:,','))") degree, (bench(ibench)%timing%mean, ibench = 1, size(bench))
52
53 end do loopOverArraySize
54 write(*,"(*(g0,:,' '))") dumsum
55 write(*,"(*(g0,:,' '))")
56
57 close(fileUnit)
58
59contains
60
61 !%%%%%%%%%%%%%%%%%%%%
62 ! procedure wrappers.
63 !%%%%%%%%%%%%%%%%%%%%
64
65 subroutine setOverhead()
66 call getDummy()
67 end subroutine
68
69 subroutine getDummy()
70 dumsum = dumsum + rootCount(ibench)
71 end subroutine
72
73 subroutine setPolyRootEigen()
74 use pm_polynomial, only: setPolyRoot, method => eigen
75 call setPolyRoot(root(1 : degree), rootCount(ibench), coef(0 : degree), method)
76 call getDummy()
77 end subroutine
78
79 subroutine setPolyRootJenkins()
80 use pm_polynomial, only: setPolyRoot, method => jenkins
81 call setPolyRoot(root(1 : degree), rootCount(ibench), coef(0 : degree), method)
82 call getDummy()
83 end subroutine
84
85 subroutine setPolyRootLaguerre()
86 use pm_polynomial, only: setPolyRoot, method => laguerre
87 call setPolyRoot(root(1 : degree), rootCount(ibench), coef(0 : degree), method)
88 call getDummy()
89 end subroutine
90
91#if 0
92 subroutine setPolyRootSkowron()
93 use pm_polynomial, only: cmplx_roots_gen
94 call cmplx_roots_gen(root(1 : degree), coef(0 : degree), degree, polish_roots_after = .true., use_roots_as_starting_points = .false.)
95 rootCount(ibench) = root(1)%re
96 call getDummy()
97 end subroutine
98#endif
99
100end program benchmark
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...
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
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
Todo:
Critical Priority: The current implementation relies on several local allocations, the most important of which are called workspace in the current implementation.
Although this generic interface used to accept a workspace argument, it was subsequently removed in favor of simplicity.
A new set of interfaces can be added in the future to allow specification of scratch space arguments to avoid repeated local allocations.
This could boost performance when this generic interface is called many times repeatedly.
Todo:
Critical Priority: The generic interface for the Eigenvalue method internally uses the eigenvalue computing routine of EISPACK for upper Hessenberg matrices to compute the roots of the polynomial.
This internal implementation should be substituted with the corresponding external routine once it is implemented in the ParaMonte library.
Todo:
Low Priority: The existing implementations of the algorithms for the Jenkins-Traub method contain relics of F77 coding style in the form of a few goto statements.
These remaining goto statements should be carefully removed in the future.
Todo:
Critical Priority: The method of Skowron-Gould must be fully implemented.
Todo:
Normal Priority: A benchmark comparing the runtime performances of the complex vs. real coefficient routines of this module should be added here.


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


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