comparison src/corefcn/quadcc.cc @ 15039:e753177cde93

maint: Move non-dynamically linked functions from DLD-FUNCTIONS/ to corefcn/ directory * __contourc__.cc, __dispatch__.cc, __lin_interpn__.cc, __pchip_deriv__.cc, __qp__.cc, balance.cc, besselj.cc, betainc.cc, bsxfun.cc, cellfun.cc, colloc.cc, conv2.cc, daspk.cc, dasrt.cc, dassl.cc, det.cc, dlmread.cc, dot.cc, eig.cc, fft.cc, fft2.cc, fftn.cc, filter.cc, find.cc, gammainc.cc, gcd.cc, getgrent.cc, getpwent.cc, getrusage.cc, givens.cc, hess.cc, hex2num.cc, inv.cc, kron.cc, lookup.cc, lsode.cc, lu.cc, luinc.cc, matrix_type.cc, max.cc, md5sum.cc, mgorth.cc, nproc.cc, pinv.cc, quad.cc, quadcc.cc, qz.cc, rand.cc, rcond.cc, regexp.cc, schur.cc, spparms.cc, sqrtm.cc, str2double.cc, strfind.cc, sub2ind.cc, svd.cc, syl.cc, time.cc, tril.cc, typecast.cc: Move functions from DLD-FUNCTIONS/ to corefcn/ directory. Include "defun.h", not "defun-dld.h". Change docstring to refer to these as "Built-in Functions". * build-aux/mk-opts.pl: Generate options code with '#include "defun.h"'. Change option docstrings to refer to these as "Built-in Functions". * corefcn/module.mk: List of functions to build in corefcn/ dir. * DLD-FUNCTIONS/config-module.awk: Update to new build system. * DLD-FUNCTIONS/module-files: Remove functions which are now in corefcn/ directory. * src/Makefile.am: Update to build "convenience library" in corefcn/. Octave program now links against all other libraries + corefcn libary. * src/find-defun-files.sh: Strip $srcdir from filename. * src/link-deps.mk: Add REGEX and FFTW link dependencies for liboctinterp. * type.m, which.m: Change failing tests to use 'amd', still a dynamic function, rather than 'dot', which isn't.
author Rik <rik@octave.org>
date Fri, 27 Jul 2012 15:35:00 -0700
parents src/DLD-FUNCTIONS/quadcc.cc@5ae9f0f77635
children
comparison
equal deleted inserted replaced
15038:ab18578c2ade 15039:e753177cde93
1 /*
2
3 Copyright (C) 2010-2012 Pedro Gonnet
4
5 This file is part of Octave.
6
7 Octave is free software; you can redistribute it and/or modify it
8 under the terms of the GNU General Public License as published by the
9 Free Software Foundation; either version 3 of the License, or (at your
10 option) any later version.
11
12 Octave is distributed in the hope that it will be useful, but WITHOUT
13 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
14 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
15 for more details.
16
17 You should have received a copy of the GNU General Public License
18 along with Octave; see the file COPYING. If not, see
19 <http://www.gnu.org/licenses/>.
20
21 */
22
23 #ifdef HAVE_CONFIG_H
24 #include <config.h>
25 #endif
26
27 #include "lo-ieee.h"
28 #include "parse.h"
29 #include "variables.h"
30
31 #include "defun.h"
32 #include "error.h"
33 #include "oct-obj.h"
34 #include "utils.h"
35
36 //#include "oct.h"
37 //#include "defun.h"
38
39 /* Define the size of the interval heap. */
40 #define cquad_heapsize 200
41
42
43 /* Data of a single interval */
44 typedef struct
45 {
46 double a, b;
47 double c[64];
48 double fx[33];
49 double igral, err;
50 int depth, rdepth, ndiv;
51 } cquad_ival;
52
53 /* Some constants and matrices that we'll need. */
54
55 static const double xi[33] = {
56 -1., -0.99518472667219688624, -0.98078528040323044912,
57 -0.95694033573220886493, -0.92387953251128675612,
58 -0.88192126434835502970, -0.83146961230254523708,
59 -0.77301045336273696082, -0.70710678118654752440,
60 -0.63439328416364549822, -0.55557023301960222475,
61 -0.47139673682599764857, -0.38268343236508977173,
62 -0.29028467725446236764, -0.19509032201612826785,
63 -0.098017140329560601995, 0., 0.098017140329560601995,
64 0.19509032201612826785, 0.29028467725446236764, 0.38268343236508977173,
65 0.47139673682599764857, 0.55557023301960222475, 0.63439328416364549822,
66 0.70710678118654752440, 0.77301045336273696082, 0.83146961230254523708,
67 0.88192126434835502970, 0.92387953251128675612, 0.95694033573220886493,
68 0.98078528040323044912, 0.99518472667219688624, 1.
69 };
70
71 static const double bee[68] = {
72 0.00000000000000e+00, 2.28868854108532e-01, 0.00000000000000e+00,
73 -8.15740215243451e-01, 0.00000000000000e+00, 5.31212715259731e-01,
74 0.00000000000000e+00, 1.38538036812454e-02, 0.00000000000000e+00,
75 3.74405228908818e-02, 0.00000000000000e+00, 2.12224115039342e-01,
76 0.00000000000000e+00, -8.16362644507898e-01, 0.00000000000000e+00,
77 5.35648426691481e-01, 0.00000000000000e+00, 1.52417902753662e-03,
78 0.00000000000000e+00, 2.63058840550873e-03, 0.00000000000000e+00,
79 4.15292106318904e-03, 0.00000000000000e+00, 6.97106011119775e-03,
80 0.00000000000000e+00, 1.35535708431058e-02, 0.00000000000000e+00,
81 3.52132898424856e-02, 0.00000000000000e+00, 2.06946714741884e-01,
82 0.00000000000000e+00, -8.15674251283876e-01, 0.00000000000000e+00,
83 5.38841175520580e-01, 0.00000000000000e+00, 1.84909689577590e-04,
84 0.00000000000000e+00, 2.90936325007499e-04, 0.00000000000000e+00,
85 3.84877750950089e-04, 0.00000000000000e+00, 4.86436656735046e-04,
86 0.00000000000000e+00, 6.08688640346879e-04, 0.00000000000000e+00,
87 7.66732830740331e-04, 0.00000000000000e+00, 9.82753336104205e-04,
88 0.00000000000000e+00, 1.29359957505615e-03, 0.00000000000000e+00,
89 1.76616363801885e-03, 0.00000000000000e+00, 2.53323433039089e-03,
90 0.00000000000000e+00, 3.88872172121956e-03, 0.00000000000000e+00,
91 6.58635106468291e-03, 0.00000000000000e+00, 1.30326736343254e-02,
92 0.00000000000000e+00, 3.44353850696714e-02, 0.00000000000000e+00,
93 2.05025409531915e-01, 0.00000000000000e+00, -8.14985893995401e-01,
94 0.00000000000000e+00, 5.40679930965238e-01
95 };
96
97 static const double Lalpha[33] = {
98 5.77350269189626e-01, 5.16397779494322e-01, 5.07092552837110e-01,
99 5.03952630678970e-01, 5.02518907629606e-01, 5.01745206004255e-01,
100 5.01280411827603e-01, 5.00979432868120e-01, 5.00773395667191e-01,
101 5.00626174321759e-01, 5.00517330712619e-01, 5.00434593736979e-01,
102 5.00370233297676e-01, 5.00319182924304e-01, 5.00278009473803e-01,
103 5.00244319584578e-01, 5.00216403386025e-01, 5.00193012939056e-01,
104 5.00173220168024e-01, 5.00156323280355e-01, 5.00141783641018e-01,
105 5.00129182278347e-01, 5.00118189340972e-01, 5.00108542278496e-01,
106 5.00100030010004e-01, 5.00092481273333e-01, 5.00085755939229e-01,
107 5.00079738458365e-01, 5.00074332862969e-01, 5.00069458915387e-01,
108 5.00065049112355e-01, 5.00061046334395e-01, 5.00057401986298e-01
109 };
110
111 static const double Lgamma[33] = {
112 0.0, 0.0, 5.16397779494322e-01, 5.07092552837110e-01, 5.03952630678970e-01,
113 5.02518907629606e-01, 5.01745206004255e-01, 5.01280411827603e-01,
114 5.00979432868120e-01, 5.00773395667191e-01, 5.00626174321759e-01,
115 5.00517330712619e-01, 5.00434593736979e-01, 5.00370233297676e-01,
116 5.00319182924304e-01, 5.00278009473803e-01, 5.00244319584578e-01,
117 5.00216403386025e-01, 5.00193012939056e-01, 5.00173220168024e-01,
118 5.00156323280355e-01, 5.00141783641018e-01, 5.00129182278347e-01,
119 5.00118189340972e-01, 5.00108542278496e-01, 5.00100030010003e-01,
120 5.00092481273333e-01, 5.00085755939229e-01, 5.00079738458365e-01,
121 5.00074332862969e-01, 5.00069458915387e-01, 5.00065049112355e-01,
122 5.00061046334395e-01
123 };
124
125 static const double V1inv[5 * 5] = {
126 .47140452079103168293e-1, .37712361663282534635, .56568542494923801952,
127 .37712361663282534635, .47140452079103168293e-1,
128 -.81649658092772603273e-1, -.46188021535170061160, 0,
129 .46188021535170061160, .81649658092772603273e-1, .15058465048420853962,
130 .12046772038736683169, -.54210474174315074262, .12046772038736683169,
131 .15058465048420853962, -.21380899352993950775, .30237157840738178177, -0.,
132 -.30237157840738178177, .21380899352993950775, .10774960475223581324,
133 -.21549920950447162648, .21549920950447162648, -.21549920950447162648,
134 .10774960475223581324
135 };
136
137 static const double V2inv[9 * 9] = {
138 .11223917161691230546e-1, .10339219839658349826, .19754094204576565761,
139 .25577315077753587922, .27835314560994251755, .25577315077753587922,
140 .19754094204576565761, .10339219839658349826, .11223917161691230546e-1,
141 -.19440394783993476970e-1, -.16544884625069155470, -.24193725566041460608,
142 -.16953338808305493604, 0.0, .16953338808305493604, .24193725566041460608,
143 .16544884625069155470, .19440394783993476970e-1, .26466393115406349388e-1,
144 .17766815796285469394, .11316664642449611462, -.16306601003711325980,
145 -.30847037493128779631, -.16306601003711325980, .11316664642449611462,
146 .17766815796285469394, .26466393115406349388e-1,
147 -.32395302049990834508e-1, -.15521142532414866547,
148 .88573492664788602740e-1, .29570405784974857322, 0.0,
149 -.29570405784974857322, -.88573492664788602740e-1, .15521142532414866547,
150 .32395302049990834508e-1, .41442155673936851246e-1,
151 .98186757907405608245e-1, -.23056908429499411784,
152 -.68047008326360625520e-1, .31797435808002456774,
153 -.68047008326360625520e-1, -.23056908429499411784,
154 .98186757907405608245e-1, .41442155673936851246e-1,
155 -.49981120317798783134e-1, -.24861810572835756217e-1,
156 .23561326072010832539, -.24472785656448415351, 0.0, .24472785656448415351,
157 -.23561326072010832539, .24861810572835756217e-1,
158 .49981120317798783134e-1, .79691635865674781228e-1,
159 -.95725617891693941833e-1, -.57957553356854386344e-1,
160 .21164072460540271452, -.27529837844505833514, .21164072460540271452,
161 -.57957553356854386344e-1, -.95725617891693941833e-1,
162 .79691635865674781228e-1,
163 -.10894869830716590913, .20131094491947531782, -.15407672674888869038,
164 .83385723639789791384e-1, 0.0, -.83385723639789791384e-1,
165 .15407672674888869038, -.20131094491947531782, .10894869830716590913,
166 .54581057089643838221e-1, -.10916211417928767644, .10916211417928767644,
167 -.10916211417928767644, .10916211417928767644, -.10916211417928767644,
168 .10916211417928767644, -.10916211417928767644, .54581057089643838221e-1
169 };
170
171 static const double V3inv[17 * 17] = {
172 .27729677693590098996e-2, .26423663180333065153e-1,
173 .53374068493933898312e-1, .77007854739523195947e-1,
174 .98257061072911596869e-1, .11538049741786835604, .12832134344120884559,
175 .13612785914022865001, .13888293186236181317, .13612785914022865001,
176 .12832134344120884559, .11538049741786835604, .98257061072911596869e-1,
177 .77007854739523195947e-1, .53374068493933898312e-1,
178 .26423663180333065153e-1, .27729677693590098996e-2,
179 -.48029210642807413690e-2, -.44887724635478800254e-1,
180 -.85409520147301089416e-1, -.11090267822061423050, -.12033983162705862441,
181 -.11102786862182788886, -.85054870109799336515e-1,
182 -.45998467987742225160e-1, 0.0, .45998467987742225160e-1,
183 .85054870109799336515e-1, .11102786862182788886, .12033983162705862441,
184 .11090267822061423050, .85409520147301089416e-1, .44887724635478800254e-1,
185 .48029210642807413690e-2, .62758546879582030087e-2,
186 .55561297093529155869e-1,
187 .93281491021051539742e-1, .92320151237493695139e-1,
188 .55077987469605684531e-1,
189 -.96998141716497488255e-2, -.80285961895427405567e-1,
190 -.13496839655913850224,
191 -.15512521776684524331, -.13496839655913850224, -.80285961895427405567e-1,
192 -.96998141716497488255e-2, .55077987469605684531e-1,
193 .92320151237493695139e-1, .93281491021051539742e-1,
194 .55561297093529155869e-1, .62758546879582030087e-2,
195 -.74850969394858555939e-2, -.61751608943839234096e-1,
196 -.82974150437304275958e-1, -.38437763431942633378e-1,
197 .45745502025779701366e-1, .12369235652734542162, .14720439712852868239,
198 .98768034347019704401e-1, 0.0,
199 -.98768034347019704401e-1, -.14720439712852868239, -.12369235652734542162,
200 -.45745502025779701366e-1, .38437763431942633378e-1,
201 .82974150437304275958e-1, .61751608943839234096e-1,
202 .74850969394858555939e-2, .86710099994384056338e-2,
203 .64006230103659573344e-1, .58517426396091675690e-1,
204 -.29743410528985802680e-1,
205 -.11934127779157114754, -.12686773515361299409, -.30729137153877447035e-1,
206 .97307836256600731568e-1, .15635811574451401023, .97307836256600731568e-1,
207 -.30729137153877447035e-1, -.12686773515361299409, -.11934127779157114754,
208 -.29743410528985802680e-1, .58517426396091675690e-1,
209 .64006230103659573344e-1, .86710099994384056338e-2,
210 -.97486395666294840165e-2, -.62995604908060224672e-1,
211 -.24373234450275529219e-1, .87760984413626872730e-1,
212 .12205204576993351394,
213 .16216004196864002088e-1, -.12422320942156845775, -.13682714580929614678,
214 0.0, .13682714580929614678, .12422320942156845775,
215 -.16216004196864002088e-1, -.12205204576993351394,
216 -.87760984413626872730e-1, .24373234450275529219e-1,
217 .62995604908060224672e-1, .97486395666294840165e-2,
218 .10956271233750488468e-1, .58613204255294358939e-1,
219 -.13306063940736618859e-1, -.11606666444978454399,
220 -.52059598001115805639e-1, .10868540217796151849, .12594452879014618005,
221 -.44678658254872910434e-1, -.15617684362128533405,
222 -.44678658254872910434e-1, .12594452879014618005, .10868540217796151849,
223 -.52059598001115805639e-1, -.11606666444978454399,
224 -.13306063940736618859e-1, .58613204255294358939e-1,
225 .10956271233750488468e-1, -.12098893000863087230e-1,
226 -.51626244709126208453e-1, .48919433304746979330e-1,
227 .10467644465949427090,
228 -.48729879523084673782e-1, -.13668732103524749234, .28190838706814496438e-1,
229 .15434223333238741600, 0.0, -.15434223333238741600,
230 -.28190838706814496438e-1, .13668732103524749234,
231 .48729879523084673782e-1, -.10467644465949427090,
232 -.48919433304746979330e-1, .51626244709126208453e-1,
233 .12098893000863087230e-1, .13542668300437944822e-1,
234 .41712033418258689308e-1,
235 -.76190463272803434388e-1, -.58303943170068132010e-1, .12158068748245606853,
236 .42121099930651007882e-1, -.14684425840766337756,
237 -.16108203535058647043e-1, .15698075850757976092,
238 -.16108203535058647043e-1, -.14684425840766337756,
239 .42121099930651007882e-1, .12158068748245606853,
240 -.58303943170068132010e-1, -.76190463272803434388e-1,
241 .41712033418258689308e-1, .13542668300437944822e-1,
242 -.14939634995117694417e-1, -.30047246373341564039e-1,
243 .91624635082546425678e-1, -.79133374319110026377e-2,
244 -.12292558212072233355, .90013382617762643524e-1,
245 .84013717196539593395e-1, -.14813033309980695856, 0.0,
246 .14813033309980695856, -.84013717196539593395e-1,
247 -.90013382617762643524e-1,
248 .12292558212072233355, .79133374319110026377e-2, -.91624635082546425678e-1,
249 .30047246373341564039e-1, .14939634995117694417e-1,
250 .16986031342807474208e-1,
251 .15760203882617033601e-1, -.91494054040950941996e-1,
252 .70082459207876130806e-1,
253 .53390713710144539104e-1, -.14340746778352039430, .84048122493418898508e-1,
254 .72456667788091316868e-1, -.15564535320096811360,
255 .72456667788091316868e-1, .84048122493418898508e-1,
256 -.14340746778352039430, .53390713710144539104e-1,
257 .70082459207876130806e-1, -.91494054040950941996e-1,
258 .15760203882617033601e-1,
259 .16986031342807474208e-1, -.18994065631858742028e-1,
260 -.82901821370405592927e-3, .77239669773015192888e-1,
261 -.10850735431039424680, .47524484622086496464e-1,
262 .69148184871588737021e-1, -.14829314646228194928, .11992057742398672066,
263 0.0, -.11992057742398672066, .14829314646228194928,
264 -.69148184871588737021e-1, -.47524484622086496464e-1,
265 .10850735431039424680, -.77239669773015192888e-1,
266 .82901821370405592927e-3, .18994065631858742028e-1,
267 .22761703826371535132e-1, -.17728848711449643358e-1,
268 -.47496371572480503788e-1, .10659958402328690063, -.11696013966166296514,
269 .63073750910894244526e-1, .32928881123602721303e-1,
270 -.12280950532497593683, .15926189077282729505, -.12280950532497593683,
271 .32928881123602721303e-1, .63073750910894244526e-1,
272 -.11696013966166296514, .10659958402328690063, -.47496371572480503788e-1,
273 -.17728848711449643358e-1, .22761703826371535132e-1,
274 -.26493215276042203434e-1, .35579780856128386192e-1,
275 .10447309718398935122e-1, -.68616154085314996709e-1,
276 .11775363082763954214, -.13918901977011837274, .12312819418827395690,
277 -.72053565748259077905e-1, 0.0, .72053565748259077905e-1,
278 -.12312819418827395690, .13918901977011837274, -.11775363082763954214,
279 .68616154085314996709e-1, -.10447309718398935122e-1,
280 -.35579780856128386192e-1,
281 .26493215276042203434e-1, .40742523354399706918e-1,
282 -.73124912999529117195e-1, .49317266444153837821e-1,
283 -.13686605413876015320e-1, -.28342624942191100464e-1,
284 .70371855298258216249e-1, -.10600251632853603875, .12981016288391131812,
285 -.13817029659318161476, .12981016288391131812, -.10600251632853603875,
286 .70371855298258216249e-1, -.28342624942191100464e-1,
287 -.13686605413876015320e-1,
288 .49317266444153837821e-1, -.73124912999529117195e-1,
289 .40742523354399706918e-1, -.54944368958699908688e-1,
290 .10777725663147408190, -.10152395581538265428, .91369146312596428468e-1,
291 -.77703071757424700773e-1, .61050911730999815031e-1,
292 -.42052599404498348871e-1, .21438229266251454773e-1, 0.0,
293 -.21438229266251454773e-1, .42052599404498348871e-1,
294 -.61050911730999815031e-1, .77703071757424700773e-1,
295 -.91369146312596428468e-1,
296 .10152395581538265428, -.10777725663147408190, .54944368958699908688e-1,
297 .27485608464748840573e-1, -.54971216929497681146e-1,
298 .54971216929497681146e-1,
299 -.54971216929497681146e-1, .54971216929497681146e-1,
300 -.54971216929497681146e-1, .54971216929497681146e-1,
301 -.54971216929497681146e-1, .54971216929497681146e-1,
302 -.54971216929497681146e-1, .54971216929497681146e-1,
303 -.54971216929497681146e-1, .54971216929497681146e-1,
304 -.54971216929497681146e-1, .54971216929497681146e-1,
305 -.54971216929497681146e-1, .27485608464748840573e-1
306 };
307
308 static const double V4inv[33 * 33] = {
309 .69120897476690862600e-3, .66419939766331555194e-2,
310 .13600665164323186111e-1, .20122785860913684493e-1,
311 .26583214101668429944e-1, .32712713318999268739e-1,
312 .38576221976287138036e-1, .44033030938268925133e-1,
313 .49092709529622799673e-1, .53657949874312515646e-1,
314 .57724533144734311859e-1, .61219564530655179096e-1,
315 .64138907503837875026e-1, .66427905189318792009e-1,
316 .68088956652280022887e-1, .69083051391555695878e-1,
317 .69422738116739271449e-1, .69083051391555695878e-1,
318 .68088956652280022887e-1, .66427905189318792009e-1,
319 .64138907503837875026e-1, .61219564530655179096e-1,
320 .57724533144734311859e-1, .53657949874312515646e-1,
321 .49092709529622799673e-1, .44033030938268925133e-1,
322 .38576221976287138036e-1, .32712713318999268739e-1,
323 .26583214101668429944e-1, .20122785860913684493e-1,
324 .13600665164323186111e-1, .66419939766331555194e-2,
325 .69120897476690862600e-3, -.11972090629438798134e-2,
326 -.11448874821643225573e-1, -.23104401104002905904e-1,
327 -.33352899418646530133e-1, -.42538626424075425908e-1,
328 -.49969730733911825941e-1, -.55555454015360728353e-1,
329 -.58955533624852604918e-1, -.60126044219122513907e-1,
330 -.58959430451175833624e-1, -.55546925396227130606e-1,
331 -.49984739749347973762e-1, -.42513009141170294365e-1,
332 -.33399140950669746346e-1, -.23007690803851790829e-1,
333 -.11728275717520066169e-1, 0.0, .11728275717520066169e-1,
334 .23007690803851790829e-1, .33399140950669746346e-1,
335 .42513009141170294365e-1, .49984739749347973762e-1,
336 .55546925396227130606e-1, .58959430451175833624e-1,
337 .60126044219122513907e-1, .58955533624852604918e-1,
338 .55555454015360728353e-1, .49969730733911825941e-1,
339 .42538626424075425908e-1, .33352899418646530133e-1,
340 .23104401104002905904e-1, .11448874821643225573e-1,
341 .11972090629438798134e-2, .15501585012936019146e-2,
342 .14628781502199620482e-1, .28684915921474815271e-1,
343 .39299396074628048026e-1, .46393418975496284204e-1,
344 .48756902531094699526e-1, .46331333488337494692e-1,
345 .39012645376980228775e-1, .27452795421085791153e-1,
346 .12430953621169863781e-1, -.47682978056024928800e-2,
347 -.22825828045428973853e-1,
348 -.40195512090720278312e-1, -.55503004262826221955e-1,
349 -.67424537752827046308e-1, -.75020199300113606452e-1,
350 -.77607844312483656131e-1, -.75020199300113606452e-1,
351 -.67424537752827046308e-1, -.55503004262826221955e-1,
352 -.40195512090720278312e-1, -.22825828045428973853e-1,
353 -.47682978056024928800e-2, .12430953621169863781e-1,
354 .27452795421085791153e-1, .39012645376980228775e-1,
355 .46331333488337494692e-1, .48756902531094699526e-1,
356 .46393418975496284204e-1, .39299396074628048026e-1,
357 .28684915921474815271e-1, .14628781502199620482e-1,
358 .15501585012936019146e-2, -.18377757558949194214e-2,
359 -.17050470050949761565e-1, -.31952119564923250836e-1,
360 -.40197423449026348155e-1,
361 -.41205649520281371624e-1, -.33909965817492272248e-1,
362 -.19393664422115332144e-1, .56661049630886784692e-3,
363 .22948272173686561721e-1, .44489719570904738207e-1,
364 .61790363672287920596e-1, .72121014727028013894e-1,
365 .73627151185287858579e-1, .65784665375961398923e-1,
366 .49369676372333667559e-1, .26444326317059715065e-1, 0.0,
367 -.26444326317059715065e-1, -.49369676372333667559e-1,
368 -.65784665375961398923e-1, -.73627151185287858579e-1,
369 -.72121014727028013894e-1, -.61790363672287920596e-1,
370 -.44489719570904738207e-1, -.22948272173686561721e-1,
371 -.56661049630886784692e-3, .19393664422115332144e-1,
372 .33909965817492272248e-1, .41205649520281371624e-1,
373 .40197423449026348155e-1, .31952119564923250836e-1,
374 .17050470050949761565e-1, .18377757558949194214e-2,
375 .20942714740729767769e-2, .18935902405146518232e-1,
376 .33335840852491735126e-1, .36770680999102286065e-1,
377 .28873194534132768509e-1, .10267303017729535513e-1,
378 -.14607738306201572890e-1, -.40139568545572305818e-1,
379 -.59808326733858291561e-1, -.68528358823372627506e-1,
380 -.63306535387619244879e-1, -.44508601817574921056e-1,
381 -.15449116105605395357e-1, .17941083795006546367e-1,
382 .48747356011657242123e-1, .70329553984201665523e-1,
383 .78106117292526169663e-1, .70329553984201665523e-1,
384 .48747356011657242123e-1, .17941083795006546367e-1,
385 -.15449116105605395357e-1, -.44508601817574921056e-1,
386 -.63306535387619244879e-1, -.68528358823372627506e-1,
387 -.59808326733858291561e-1,
388 -.40139568545572305818e-1, -.14607738306201572890e-1,
389 .10267303017729535513e-1, .28873194534132768509e-1,
390 .36770680999102286065e-1, .33335840852491735126e-1,
391 .18935902405146518232e-1, .20942714740729767769e-2,
392 -.23245285491878278419e-2, -.20401404737639389919e-1,
393 -.33019548231022514097e-1, -.29709828426463720091e-1,
394 -.11760070922697422156e-1, .15987584743850393793e-1,
395 .43619012891472813485e-1, .61177322409671487721e-1,
396 .61144030218486655594e-1,
397 .41895377620089086167e-1, .80232011820644308033e-2,
398 -.30574701186675900915e-1,
399 -.62072243008844865848e-1, -.76336186183574765586e-1,
400 -.68435466095345537115e-1, -.40237669208466966207e-1, 0.0,
401 .40237669208466966207e-1, .68435466095345537115e-1,
402 .76336186183574765586e-1, .62072243008844865848e-1,
403 .30574701186675900915e-1, -.80232011820644308033e-2,
404 -.41895377620089086167e-1, -.61144030218486655594e-1,
405 -.61177322409671487721e-1, -.43619012891472813485e-1,
406 -.15987584743850393793e-1, .11760070922697422156e-1,
407 .29709828426463720091e-1, .33019548231022514097e-1,
408 .20401404737639389919e-1, .23245285491878278419e-2,
409 .25451717261579269307e-2, .21480418595666878775e-1,
410 .31177212469293007998e-1, .19816333607013379373e-1,
411 -.72439496274458793681e-2, -.38404203906598342397e-1,
412 -.57633632255322221046e-1, -.54070547403585392952e-1,
413 -.26249823354368866005e-1, .15643058212336881516e-1,
414 .54539832735118677194e-1, .73283028002473989724e-1,
415 .62835303524135936213e-1, .26175977027801048141e-1,
416 -.22193636309998606610e-1, -.62597049956093311234e-1,
417 -.78206986173170212505e-1, -.62597049956093311234e-1,
418 -.22193636309998606610e-1, .26175977027801048141e-1,
419 .62835303524135936213e-1,
420 .73283028002473989724e-1, .54539832735118677194e-1,
421 .15643058212336881516e-1,
422 -.26249823354368866005e-1, -.54070547403585392952e-1,
423 -.57633632255322221046e-1, -.38404203906598342397e-1,
424 -.72439496274458793681e-2, .19816333607013379373e-1,
425 .31177212469293007998e-1, .21480418595666878775e-1,
426 .25451717261579269307e-2, -.27506573922483820005e-2,
427 -.22224442095099251870e-1, -.27949927254215773020e-1,
428 -.80918481053370034987e-2, .25121859354449306916e-1,
429 .51563535009373061074e-1, .51936965107145960512e-1,
430 .22146626648171527753e-1,
431 -.24172689882103382748e-1, -.61731229104853568296e-1,
432 -.68477262429344201201e-1, -.38311232728303704742e-1,
433 .14160578713659552679e-1, .61248813427564184033e-1,
434 .77136328841293031805e-1, .52514801765183697988e-1, 0.0,
435 -.52514801765183697988e-1, -.77136328841293031805e-1,
436 -.61248813427564184033e-1, -.14160578713659552679e-1,
437 .38311232728303704742e-1,
438 .68477262429344201201e-1, .61731229104853568296e-1,
439 .24172689882103382748e-1,
440 -.22146626648171527753e-1, -.51936965107145960512e-1,
441 -.51563535009373061074e-1, -.25121859354449306916e-1,
442 .80918481053370034987e-2, .27949927254215773020e-1,
443 .22224442095099251870e-1, .27506573922483820005e-2,
444 .29562461131654311467e-2, .22630271480554450613e-1,
445 .23547399831373800971e-1, -.43964593440902476642e-2,
446 -.39055315767504970597e-1, -.52369643937940066804e-1,
447 -.28506131614971613422e-1, .19906048093338832322e-1,
448 .60408880866392420279e-1, .62493397473656883090e-1,
449 .21391278377641297859e-1, -.37302864786623254746e-1,
450 -.73665127933539496872e-1, -.61706142476854010202e-1,
451 -.78065168882546327888e-2, .52335307373945544428e-1,
452 .78278746279419264777e-1, .52335307373945544428e-1,
453 -.78065168882546327888e-2, -.61706142476854010202e-1,
454 -.73665127933539496872e-1, -.37302864786623254746e-1,
455 .21391278377641297859e-1, .62493397473656883090e-1,
456 .60408880866392420279e-1, .19906048093338832322e-1,
457 -.28506131614971613422e-1, -.52369643937940066804e-1,
458 -.39055315767504970597e-1, -.43964593440902476642e-2,
459 .23547399831373800971e-1, .22630271480554450613e-1,
460 .29562461131654311467e-2, -.31515718415504761303e-2,
461 -.22739451096655080673e-1, -.18157123602272119779e-1,
462 .16496480897167303621e-1, .46921166788569301124e-1,
463 .40644395739978416354e-1, -.46275803430732216900e-2,
464 -.52883375891308909486e-1, -.61116483226324111734e-1,
465 -.17411698764545629853e-1, .44773430013166822765e-1,
466 .73441577962383869198e-1, .42127368371995472815e-1,
467 -.25504645957196772465e-1, -.74126818045972742488e-1,
468 -.62780077864719287317e-1, 0.0, .62780077864719287317e-1,
469 .74126818045972742488e-1, .25504645957196772465e-1,
470 -.42127368371995472815e-1, -.73441577962383869198e-1,
471 -.44773430013166822765e-1, .17411698764545629853e-1,
472 .61116483226324111734e-1, .52883375891308909486e-1,
473 .46275803430732216900e-2, -.40644395739978416354e-1,
474 -.46921166788569301124e-1, -.16496480897167303621e-1,
475 .18157123602272119779e-1, .22739451096655080673e-1,
476 .31515718415504761303e-2, .33536559294882188208e-2,
477 .22535348942792006185e-1,
478 .12048629300953560767e-1, -.27166076791299493403e-1,
479 -.47492745604230978367e-1, -.19246623430993153174e-1,
480 .36231297307556299322e-1, .61713617181636122004e-1,
481 .25928029734266134490e-1, -.40478700752883602818e-1,
482 -.71053889866326412049e-1, -.31870824482961751482e-1,
483 .41515251100219081281e-1, .76481960760098381651e-1,
484 .36726509155999912440e-1, -.40090067032627055969e-1,
485 -.78270742903374539397e-1, -.40090067032627055969e-1,
486 .36726509155999912440e-1, .76481960760098381651e-1,
487 .41515251100219081281e-1, -.31870824482961751482e-1,
488 -.71053889866326412049e-1, -.40478700752883602818e-1,
489 .25928029734266134490e-1, .61713617181636122004e-1,
490 .36231297307556299322e-1, -.19246623430993153174e-1,
491 -.47492745604230978367e-1, -.27166076791299493403e-1,
492 .12048629300953560767e-1, .22535348942792006185e-1,
493 .33536559294882188208e-2,
494 -.35481220456925318865e-2, -.22062913693073191150e-1,
495 -.54487362861834144999e-2, .35438821865804087489e-1,
496 .40733077820527411302e-1, -.67403098138950720914e-2,
497 -.55559584405239171054e-1, -.42417050790865158745e-1,
498 .24499901971884704925e-1, .68721232891705409302e-1,
499 .34086082787461126592e-1, -.43441000373118474002e-1,
500 -.73878085292669148950e-1, -.18846995664706657127e-1,
501 .59827776178286834498e-1, .70644634584085901794e-1, 0.0,
502 -.70644634584085901794e-1, -.59827776178286834498e-1,
503 .18846995664706657127e-1, .73878085292669148950e-1,
504 .43441000373118474002e-1, -.34086082787461126592e-1,
505 -.68721232891705409302e-1, -.24499901971884704925e-1,
506 .42417050790865158745e-1, .55559584405239171054e-1,
507 .67403098138950720914e-2, -.40733077820527411302e-1,
508 -.35438821865804087489e-1, .54487362861834144999e-2,
509 .22062913693073191150e-1, .35481220456925318865e-2,
510 .37554176816665075631e-2, .21297045781589919482e-1,
511 -.13327293083183431816e-2,
512 -.40635299172764596484e-1, -.27659860508374175359e-1,
513 .31089232744083445986e-1, .56113781541334176109e-1,
514 .37577840643257763400e-2, -.60511227350664590865e-1,
515 -.46670556446129053853e-1, .33263195878575888247e-1,
516 .72757324720645228775e-1, .15011712351692283635e-1,
517 -.65601212994924119078e-1, -.60016855838843789772e-1,
518 .26220858553188665966e-1, .78322776605833552980e-1,
519 .26220858553188665966e-1, -.60016855838843789772e-1,
520 -.65601212994924119078e-1,
521 .15011712351692283635e-1, .72757324720645228775e-1,
522 .33263195878575888247e-1,
523 -.46670556446129053853e-1, -.60511227350664590865e-1,
524 .37577840643257763400e-2, .56113781541334176109e-1,
525 .31089232744083445986e-1, -.27659860508374175359e-1,
526 -.40635299172764596484e-1, -.13327293083183431816e-2,
527 .21297045781589919482e-1, .37554176816665075631e-2,
528 -.39566995305720591229e-2, -.20291873414438919995e-1,
529 .80617453830770930551e-2, .42270189157016547906e-1,
530 .10332624526759093004e-1, -.48054759547616142024e-1,
531 -.37678032941171643972e-1,
532 .36617192625732482394e-1, .61009425973424865714e-1,
533 -.95589113168026591466e-2,
534 -.71023202645076922361e-1, -.25097788086808784456e-1,
535 .62406621963267050244e-1, .56907293171100693511e-1,
536 -.36435383083882206257e-1, -.75790105119208756348e-1, 0.0,
537 .75790105119208756348e-1, .36435383083882206257e-1,
538 -.56907293171100693511e-1, -.62406621963267050244e-1,
539 .25097788086808784456e-1, .71023202645076922361e-1,
540 .95589113168026591466e-2,
541 -.61009425973424865714e-1, -.36617192625732482394e-1,
542 .37678032941171643972e-1, .48054759547616142024e-1,
543 -.10332624526759093004e-1, -.42270189157016547906e-1,
544 -.80617453830770930551e-2, .20291873414438919995e-1,
545 .39566995305720591229e-2, .41776092289182138591e-2,
546 .19013221163904414395e-1, -.14420609729849899876e-1,
547 -.40259160586844441220e-1, .86327811113710831649e-2,
548 .53564430703021034399e-1, .65469185402150431933e-2,
549 -.60383116311280629856e-1,
550 -.25657793784058876939e-1, .58745680576829226900e-1,
551 .45649937869034420296e-1,
552 -.49167932056844167772e-1, -.62696614328552187977e-1,
553 .32540234556426699997e-1, .74280410383464269758e-1,
554 -.11425672633410999870e-1, -.78280649404686404903e-1,
555 -.11425672633410999870e-1, .74280410383464269758e-1,
556 .32540234556426699997e-1, -.62696614328552187977e-1,
557 -.49167932056844167772e-1, .45649937869034420296e-1,
558 .58745680576829226900e-1, -.25657793784058876939e-1,
559 -.60383116311280629856e-1, .65469185402150431933e-2,
560 .53564430703021034399e-1,
561 .86327811113710831649e-2, -.40259160586844441220e-1,
562 -.14420609729849899876e-1, .19013221163904414395e-1,
563 .41776092289182138591e-2, -.43935502082478059199e-2,
564 -.17528761237509401631e-1, .20208915249153872535e-1,
565 .34734743119040669109e-1, -.26275910172353637955e-1,
566 -.46368003346018878786e-1,
567 .26800056330709381025e-1, .56681476464606609921e-1,
568 -.24749011438127255898e-1,
569 -.64934612189056658992e-1, .20333742247679279535e-1,
570 .71429299070059318651e-1,
571 -.14452513210428671266e-1, -.75793341281736586582e-1,
572 .74717094137184935270e-2, .78034921554757317374e-1, 0.0,
573 -.78034921554757317374e-1, -.74717094137184935270e-2,
574 .75793341281736586582e-1, .14452513210428671266e-1,
575 -.71429299070059318651e-1, -.20333742247679279535e-1,
576 .64934612189056658992e-1, .24749011438127255898e-1,
577 -.56681476464606609921e-1,
578 -.26800056330709381025e-1, .46368003346018878786e-1,
579 .26275910172353637955e-1,
580 -.34734743119040669109e-1, -.20208915249153872535e-1,
581 .17528761237509401631e-1, .43935502082478059199e-2,
582 .46379089482818671473e-2, .15791188144791287229e-1,
583 -.25134290048737455284e-1, -.26249795071946841205e-1,
584 .39960457575789924651e-1, .28111892450146525404e-1,
585 -.51026476400767918226e-1,
586 -.27266747278681831364e-1, .60708796647861610865e-1,
587 .23532306960642115854e-1,
588 -.68169639871532441111e-1, -.18204924701958312032e-1,
589 .73822890510656128485e-1, .11373392486424717019e-1,
590 -.77133324017644609416e-1, -.39295877480342619961e-2,
591 .78351902829418987960e-1, -.39295877480342619961e-2,
592 -.77133324017644609416e-1, .11373392486424717019e-1,
593 .73822890510656128485e-1, -.18204924701958312032e-1,
594 -.68169639871532441111e-1, .23532306960642115854e-1,
595 .60708796647861610865e-1, -.27266747278681831364e-1,
596 -.51026476400767918226e-1, .28111892450146525404e-1,
597 .39960457575789924651e-1, -.26249795071946841205e-1,
598 -.25134290048737455284e-1, .15791188144791287229e-1,
599 .46379089482818671473e-2, -.48780095920069827068e-2,
600 -.13886961667516983541e-1, .29071311049368895844e-1,
601 .15480559452075811600e-1, -.47527977686242313065e-1,
602 -.31929089844361042178e-2, .58015667638415922967e-1,
603 -.14547915466597622925e-1, -.61067668299848923244e-1,
604 .35093678009090186851e-1, .55378399159800654657e-1,
605 -.54277226474891610385e-1, -.42023830782434076509e-1,
606 .69197384645944912066e-1, .22610783557709586445e-1,
607 -.77269275900637030185e-1, 0.0, .77269275900637030185e-1,
608 -.22610783557709586445e-1,
609 -.69197384645944912066e-1, .42023830782434076509e-1,
610 .54277226474891610385e-1,
611 -.55378399159800654657e-1, -.35093678009090186851e-1,
612 .61067668299848923244e-1, .14547915466597622925e-1,
613 -.58015667638415922967e-1, .31929089844361042178e-2,
614 .47527977686242313065e-1, -.15480559452075811600e-1,
615 -.29071311049368895844e-1, .13886961667516983541e-1,
616 .48780095920069827068e-2, .51591759101720291381e-2,
617 .11747497650231330965e-1, -.31777863364694653331e-1,
618 -.34555825499804605557e-2, .47914131921157015198e-1,
619 -.22573685920142225247e-1, -.45320344390022666738e-1,
620 .49660630547172186418e-1, .25707858143963615736e-1,
621 -.68132707341917233933e-1, .67534860185243140399e-2,
622 .69268150370037450063e-1, -.41585011920451477177e-1,
623 -.51622397460510041271e-1, .68408139576363036148e-1,
624 .18981259024768933323e-1, -.78265472429342305554e-1,
625 .18981259024768933323e-1, .68408139576363036148e-1,
626 -.51622397460510041271e-1,
627 -.41585011920451477177e-1, .69268150370037450063e-1,
628 .67534860185243140399e-2,
629 -.68132707341917233933e-1, .25707858143963615736e-1,
630 .49660630547172186418e-1,
631 -.45320344390022666738e-1, -.22573685920142225247e-1,
632 .47914131921157015198e-1, -.34555825499804605557e-2,
633 -.31777863364694653331e-1, .11747497650231330965e-1,
634 .51591759101720291381e-2, -.54365757412741340377e-2,
635 -.94862516619529080191e-2, .33240472093448190877e-1,
636 -.88698898099681552229e-2,
637 -.40973252097216337576e-1, .42995673349795657065e-1,
638 .17320914507876958783e-1,
639 -.62201292691914856803e-1, .24726274174637346693e-1,
640 .51320859246515407288e-1,
641 -.62882063373810501763e-1, -.11003569131725622672e-1,
642 .73842261324108943465e-1, -.39240120294802923208e-1,
643 -.49293966443941122807e-1, .73552644778818223475e-1, 0.0,
644 -.73552644778818223475e-1, .49293966443941122807e-1,
645 .39240120294802923208e-1, -.73842261324108943465e-1,
646 .11003569131725622672e-1, .62882063373810501763e-1,
647 -.51320859246515407288e-1,
648 -.24726274174637346693e-1, .62201292691914856803e-1,
649 -.17320914507876958783e-1, -.42995673349795657065e-1,
650 .40973252097216337576e-1, .88698898099681552229e-2,
651 -.33240472093448190877e-1, .94862516619529080191e-2,
652 .54365757412741340377e-2, .57750194549356126240e-2,
653 .69981166020044116791e-2, -.33274982140403110792e-1,
654 .20297071020698356116e-1, .27898517839646066582e-1,
655 -.53368678853282030262e-1, .16656482990394548343e-1,
656 .46342901447260614255e-1,
657 -.60536796508149003365e-1, .29109107483842596340e-2,
658 .63224486124385124504e-1,
659 -.59028872851312033411e-1, -.14783105962696191734e-1,
660 .74269399241069253865e-1, -.49053677339382384625e-1,
661 -.33525466624811186739e-1, .78397349622515386647e-1,
662 -.33525466624811186739e-1, -.49053677339382384625e-1,
663 .74269399241069253865e-1, -.14783105962696191734e-1,
664 -.59028872851312033411e-1,
665 .63224486124385124504e-1, .29109107483842596340e-2,
666 -.60536796508149003365e-1,
667 .46342901447260614255e-1, .16656482990394548343e-1,
668 -.53368678853282030262e-1,
669 .27898517839646066582e-1, .20297071020698356116e-1,
670 -.33274982140403110792e-1,
671 .69981166020044116791e-2, .57750194549356126240e-2,
672 -.61100308370519200637e-2, -.44383614355738148616e-2,
673 .32011283412619094811e-1, -.29965011866372897633e-1,
674 -.10560682331349193348e-1, .51110336443392506342e-1,
675 -.45012284729681775492e-1, -.94236825555873320102e-2,
676 .60860695783141264746e-1,
677 -.55014628647083368926e-1, -.73474782382499482121e-2,
678 .66640148475243034781e-1, -.62533116045749887988e-1,
679 -.38650525912400102585e-2, .68429769005837003777e-1,
680 -.66984505412544901945e-1, 0.0, .66984505412544901945e-1,
681 -.68429769005837003777e-1, .38650525912400102585e-2,
682 .62533116045749887988e-1, -.66640148475243034781e-1,
683 .73474782382499482121e-2,
684 .55014628647083368926e-1, -.60860695783141264746e-1,
685 .94236825555873320102e-2,
686 .45012284729681775492e-1, -.51110336443392506342e-1,
687 .10560682331349193348e-1,
688 .29965011866372897633e-1, -.32011283412619094811e-1,
689 .44383614355738148616e-2,
690 .61100308370519200637e-2, .65409373892036191538e-2,
691 .16350101107071157065e-2, -.29301957285983144319e-1,
692 .36838667173388832579e-1, -.81922703976491586393e-2,
693 -.36955670021050133434e-1, .58374851095540469865e-1,
694 -.31977016246946181856e-1, -.25311073698658094646e-1,
695 .66674413950106952577e-1,
696 -.54865713324521039571e-1, -.39797027891537985440e-2,
697 .62830285264808449064e-1, -.72226313251296100676e-1,
698 .22560232697133353980e-1, .46455784709904033738e-1,
699 -.78200930751070349956e-1, .46455784709904033738e-1,
700 .22560232697133353980e-1, -.72226313251296100676e-1,
701 .62830285264808449064e-1, -.39797027891537985440e-2,
702 -.54865713324521039571e-1, .66674413950106952577e-1,
703 -.25311073698658094646e-1, -.31977016246946181856e-1,
704 .58374851095540469865e-1, -.36955670021050133434e-1,
705 -.81922703976491586393e-2, .36838667173388832579e-1,
706 -.29301957285983144319e-1, .16350101107071157065e-2,
707 .65409373892036191538e-2, -.69686180931868703196e-2,
708 .11849538727632789870e-2, .25452286414610537766e-1,
709 -.40522480651713943230e-1, .25694679053362813183e-1,
710 .14057118113748390637e-1, -.52037614725803488893e-1,
711 .58849342223684035589e-1,
712 -.25075229077361409271e-1, -.29559771094034181083e-1,
713 .68296746944165720199e-1, -.62890462146423984955e-1,
714 .14457636466274596445e-1, .45787612031322361496e-1,
715 -.77231759014655809742e-1, .57881203613910543657e-1, 0.0,
716 -.57881203613910543657e-1, .77231759014655809742e-1,
717 -.45787612031322361496e-1, -.14457636466274596445e-1,
718 .62890462146423984955e-1,
719 -.68296746944165720199e-1, .29559771094034181083e-1,
720 .25075229077361409271e-1,
721 -.58849342223684035589e-1, .52037614725803488893e-1,
722 -.14057118113748390637e-1, -.25694679053362813183e-1,
723 .40522480651713943230e-1, -.25452286414610537766e-1,
724 -.11849538727632789870e-2, .69686180931868703196e-2,
725 .75611653617520254845e-2, -.43290610418608409141e-2,
726 -.20277062025115566914e-1,
727 .40362947027704828926e-1, -.38938808024132120254e-1,
728 .11831186195916702262e-1,
729 .28476667401744525357e-1, -.59320969056617684621e-1,
730 .61101629747436200186e-1,
731 -.29514834848355389223e-1, -.20668001885001084821e-1,
732 .62923592802445122793e-1, -.73558456263588833115e-1,
733 .45314556330160999776e-1, .79031645918426015574e-2,
734 -.58136953576334689357e-1, .78538474524006405758e-1,
735 -.58136953576334689357e-1, .79031645918426015574e-2,
736 .45314556330160999776e-1, -.73558456263588833115e-1,
737 .62923592802445122793e-1, -.20668001885001084821e-1,
738 -.29514834848355389223e-1, .61101629747436200186e-1,
739 -.59320969056617684621e-1, .28476667401744525357e-1,
740 .11831186195916702262e-1, -.38938808024132120254e-1,
741 .40362947027704828926e-1, -.20277062025115566914e-1,
742 -.43290610418608409141e-2, .75611653617520254845e-2,
743 -.81505692478987769484e-2, .74297333588288568430e-2,
744 .14314212513540223314e-1, -.36711242251332751607e-1,
745 .46240027755503814626e-1, -.34921532671769023773e-1,
746 .46930051972353714773e-2,
747 .32842770336385381562e-1, -.61317813706529588466e-1,
748 .67000809902468893103e-1,
749 -.45337449655535622885e-1, .35794459576271920867e-2,
750 .41830061526027213385e-1,
751 -.72091371931944711708e-1, .74150028530317793195e-1,
752 -.46487632538609942002e-1, 0.0, .46487632538609942002e-1,
753 -.74150028530317793195e-1, .72091371931944711708e-1,
754 -.41830061526027213385e-1, -.35794459576271920867e-2,
755 .45337449655535622885e-1, -.67000809902468893103e-1,
756 .61317813706529588466e-1, -.32842770336385381562e-1,
757 -.46930051972353714773e-2, .34921532671769023773e-1,
758 -.46240027755503814626e-1, .36711242251332751607e-1,
759 -.14314212513540223314e-1, -.74297333588288568430e-2,
760 .81505692478987769484e-2, .90693182942442189743e-2,
761 -.11121000903959576737e-1, -.71308296141317458546e-2,
762 .29219439765986671645e-1, -.45820286629778129593e-1,
763 .49088381175879124421e-1, -.35614888785023038938e-1,
764 .78906970900092777895e-2,
765 .26262843038404929480e-1, -.56143674270125757857e-1,
766 .71700220472378350694e-1,
767 -.66963544500697307945e-1, .42215091779892228883e-1,
768 -.41338867413966866997e-2, -.36164891772995367321e-1,
769 .66584367783847858225e-1, -.77874712365070098328e-1,
770 .66584367783847858225e-1, -.36164891772995367321e-1,
771 -.41338867413966866997e-2, .42215091779892228883e-1,
772 -.66963544500697307945e-1,
773 .71700220472378350694e-1, -.56143674270125757857e-1,
774 .26262843038404929480e-1,
775 .78906970900092777895e-2, -.35614888785023038938e-1,
776 .49088381175879124421e-1,
777 -.45820286629778129593e-1, .29219439765986671645e-1,
778 -.71308296141317458546e-2, -.11121000903959576737e-1,
779 .90693182942442189743e-2, -.99848472706332791043e-2,
780 .14701271465939718856e-1, -.32917820356048383366e-3,
781 -.19201195309873585230e-1, .38409681836626963278e-1,
782 -.51647324405878909521e-1, .54522171113149311354e-1,
783 -.45040302741689006270e-1, .24183738595685990149e-1,
784 .42204134165479735097e-2, -.34317295181348742251e-1,
785 .59542472465494579941e-1, -.74135115907618101263e-1,
786 .74491937840566532596e-1, -.60042604725161994304e-1,
787 .33437677409000083169e-1, 0.0,
788 -.33437677409000083169e-1, .60042604725161994304e-1,
789 -.74491937840566532596e-1, .74135115907618101263e-1,
790 -.59542472465494579941e-1, .34317295181348742251e-1,
791 -.42204134165479735097e-2, -.24183738595685990149e-1,
792 .45040302741689006270e-1, -.54522171113149311354e-1,
793 .51647324405878909521e-1, -.38409681836626963278e-1,
794 .19201195309873585230e-1, .32917820356048383366e-3,
795 -.14701271465939718856e-1, .99848472706332791043e-2,
796 .11775579274769383373e-1, -.19892153937316935880e-1,
797 .95335114477449041055e-2, .57661528440359081617e-2,
798 -.23382690532380910781e-1, .40237257037170725321e-1,
799 -.53280289903551636474e-1, .59974361806023689068e-1,
800 -.58701684061992853224e-1, .49033407111597129616e-1,
801 -.31818835267847249219e-1, .90800541261162098886e-2,
802 .16272906819312603838e-1, -.40863896581186229487e-1,
803 .61346046297517367703e-1,
804 -.74896047554167268919e-1, .79632642148310325817e-1,
805 -.74896047554167268919e-1, .61346046297517367703e-1,
806 -.40863896581186229487e-1, .16272906819312603838e-1,
807 .90800541261162098886e-2, -.31818835267847249219e-1,
808 .49033407111597129616e-1, -.58701684061992853224e-1,
809 .59974361806023689068e-1, -.53280289903551636474e-1,
810 .40237257037170725321e-1, -.23382690532380910781e-1,
811 .57661528440359081617e-2, .95335114477449041055e-2,
812 -.19892153937316935880e-1,
813 .11775579274769383373e-1, -.13562702617218467450e-1,
814 .24885419969649845849e-1, -.18368693901908875583e-1,
815 .81673147806084084638e-2, .47890591326129587131e-2,
816 -.19313752945227974024e-1, .34065953398362954708e-1,
817 -.47667045133463415672e-1, .58820377816690514309e-1,
818 -.66424139824618415970e-1,
819 .69667606260856092515e-1, -.68102459384364543253e-1,
820 .61683024923302547971e-1,
821 -.50771943476441639136e-1, .36110771847327189215e-1,
822 -.18758028464284563358e-1, 0.0, .18758028464284563358e-1,
823 -.36110771847327189215e-1, .50771943476441639136e-1,
824 -.61683024923302547971e-1, .68102459384364543253e-1,
825 -.69667606260856092515e-1, .66424139824618415970e-1,
826 -.58820377816690514309e-1, .47667045133463415672e-1,
827 -.34065953398362954708e-1, .19313752945227974024e-1,
828 -.47890591326129587131e-2, -.81673147806084084638e-2,
829 .18368693901908875583e-1, -.24885419969649845849e-1,
830 .13562702617218467450e-1, .20576545037980523979e-1,
831 -.40093155172981004337e-1, .36954083167944054826e-1,
832 -.31856506837591907746e-1, .24996323181546255126e-1,
833 -.16637165210473614136e-1, .71002706773325085237e-2,
834 .32478629093205201133e-2,
835 -.14009562579050569518e-1, .24771262248780618922e-1,
836 -.35119395835433647559e-1, .44656290368574753171e-1,
837 -.53015448339647394161e-1, .59875631995693046782e-1,
838 -.64973208326045193862e-1, .68112280331082143373e-1,
839 -.69172215234062186994e-1, .68112280331082143373e-1,
840 -.64973208326045193862e-1, .59875631995693046782e-1,
841 -.53015448339647394161e-1, .44656290368574753171e-1,
842 -.35119395835433647559e-1, .24771262248780618922e-1,
843 -.14009562579050569518e-1, .32478629093205201133e-2,
844 .71002706773325085237e-2, -.16637165210473614136e-1,
845 .24996323181546255126e-1, -.31856506837591907746e-1,
846 .36954083167944054826e-1, -.40093155172981004337e-1,
847 .20576545037980523979e-1, -.27584914609096156163e-1,
848 .54904171411058497973e-1, -.54109756419563083153e-1,
849 .52794234894345577483e-1, -.50970276026831042415e-1,
850 .48655445537990983379e-1,
851 -.45872036510847994332e-1, .42646854695899611372e-1,
852 -.39010960357087507670e-1, .34999369144476467749e-1,
853 -.30650714874402762189e-1, .26006877464703437057e-1,
854 -.21112579608213651273e-1, .16014956068786763273e-1,
855 -.10763099747751940252e-1, .54075888924374485533e-2, 0.0,
856 -.54075888924374485533e-2, .10763099747751940252e-1,
857 -.16014956068786763273e-1,
858 .21112579608213651273e-1, -.26006877464703437057e-1,
859 .30650714874402762189e-1,
860 -.34999369144476467749e-1, .39010960357087507670e-1,
861 -.42646854695899611372e-1, .45872036510847994332e-1,
862 -.48655445537990983379e-1, .50970276026831042415e-1,
863 -.52794234894345577483e-1, .54109756419563083153e-1,
864 -.54904171411058497973e-1, .27584914609096156163e-1,
865 .13794141262469565740e-1, -.27588282524939131481e-1,
866 .27588282524939131481e-1, -.27588282524939131481e-1,
867 .27588282524939131481e-1, -.27588282524939131481e-1,
868 .27588282524939131481e-1,
869 -.27588282524939131481e-1, .27588282524939131481e-1,
870 -.27588282524939131481e-1, .27588282524939131481e-1,
871 -.27588282524939131481e-1, .27588282524939131481e-1,
872 -.27588282524939131481e-1, .27588282524939131481e-1,
873 -.27588282524939131481e-1, .27588282524939131481e-1,
874 -.27588282524939131481e-1, .27588282524939131481e-1,
875 -.27588282524939131481e-1, .27588282524939131481e-1,
876 -.27588282524939131481e-1, .27588282524939131481e-1,
877 -.27588282524939131481e-1, .27588282524939131481e-1,
878 -.27588282524939131481e-1, .27588282524939131481e-1,
879 -.27588282524939131481e-1, .27588282524939131481e-1,
880 -.27588282524939131481e-1, .27588282524939131481e-1,
881 -.27588282524939131481e-1, .13794141262469565740e-1
882 };
883
884 static const double Tleft[33 * 33] = {
885 1., -.86602540378443864678, 0., .33071891388307382381, 0.,
886 -.20728904939721249057, 0., .15128841196122722208, 0.,
887 -.11918864298744029244, 0., .98352013661686631224e-1, 0.,
888 -.83727065404940845733e-1, 0., .72893399403505841203e-1, 0.,
889 -.64544632643375022436e-1, 0., .57913170372415565639e-1, 0.,
890 -.52518242575729562263e-1, 0., .48043311993977520457e-1, 0.,
891 -.44271433659733990243e-1, 0., .41048928022856771981e-1, 0.,
892 -.38263878662008271459e-1, 0., .35832844026365304501e-1, 0., 0.,
893 .50000000000000000000, -.96824583655185422130, .57282196186948000082,
894 .21650635094610966169, -.35903516540862679125, -.97578093724974971969e-1,
895 .26203921611325660506, .55792409597991015609e-1, -.20644078533943456204,
896 -.36172381205961199479e-1, .17035068468874958194,
897 .25371838001497225980e-1, -.14501953125000000000,
898 -.18786835250972344757e-1, .12625507130328301066,
899 .14473795929590520582e-1, -.11179458309419422675,
900 -.11494434254897626155e-1, .10030855351241635862,
901 .93498556820544479096e-2, -.90964264465390582629e-1,
902 -.77546391824364392762e-2, .83213457337452292745e-1,
903 .65358085945588638605e-2, -.76680372422574234569e-1,
904 -.55835321940047427169e-2, .71098828931825789428e-1,
905 .48253327982967591019e-2, -.66274981937248958553e-1,
906 -.42118078245337801387e-2, .62064306433355646267e-1,
907 .37083386598903548973e-2, 0., 0., .25000000000000000000,
908 -.73950997288745200531, .83852549156242113615, -.23175620272173946716,
909 -.37791833195149451496, .25710129174850522325, .21608307321780204633,
910 -.22844049245646009157, -.14009503000335388415, .19897685605518413847,
911 .98264706042471226893e-1, -.17445445004279014046,
912 -.72761100054958328401e-1, .15463589893742108388,
913 .56056770591708784481e-1, -.13855313872640495158,
914 -.44517752443294564781e-1, .12534277657695128850,
915 .36211835346039665762e-1, -.11434398255136139683,
916 -.30033588409423828125e-1, .10506705408753910481,
917 .25313077840725783008e-1, -.97149327637744872155e-1,
918 -.21624927200393328444e-1, .90319582367202122625e-1,
919 .18688433567711780666e-1, -.84372291635345108584e-1,
920 -.16312261561845420752e-1, .79149526894804751586e-1,
921 .14362333871852474757e-1, 0., 0., 0., .12500000000000000000,
922 -.49607837082461073572, .82265291131801144317, -.59621200088559103072,
923 -.80054302859059362371e-1, .42612156697795759420,
924 -.90098145270865592887e-1, -.29769623255090078484, .13630307904779758221,
925 .21638835185708931831, -.14600247270306082052, -.16348801804014290453,
926 .14340708728599057249, .12755243353979286190, -.13661523715071346961,
927 -.10215585947881057394, .12864248070157166547, .83592528025348693602e-1,
928 -.12066728689302565222, -.69633728678718053052e-1, .11314245177331919532,
929 .58882939251410088028e-1, -.10621835858758221487,
930 -.50432266865187597572e-1, .99916834723527771581e-1,
931 .43672094283057258509e-1, -.94206380251950852413e-1,
932 -.38181356812697746418e-1, .89035739656537771225e-1,
933 .33661934598216332678e-1, 0., 0., 0., 0., .62500000000000000000e-1,
934 -.31093357409581873586, .67604086414949799246, -.75644205980613611039,
935 .28990586430124175741, .30648508196770360914, -.35801372616842500052,
936 -.91326869828709014708e-1, .31127929687500000000,
937 -.90915752838698393094e-2, -.25637381283965534330,
938 .57601077850322797594e-1, .21019685709225757945,
939 -.81244992138514014256e-1, -.17375078516720988858,
940 .92289437277967051125e-1, .14527351914265391374,
941 -.96675340792832019889e-1, -.12289485697108543415,
942 .97448175340011084006e-1, .10511755943298339844,
943 -.96242247086378239657e-1, -.90822942272780513537e-1,
944 .93966350452322132384e-1, .79189411876493712558e-1,
945 -.91139307067989309325e-1, -.69613039934383197265e-1,
946 .88062491671135767870e-1, .61646331729340817494e-1, 0., 0., 0., 0., 0.,
947 .31250000000000000000e-1, -.18684782411095934408, .50176689760410660236,
948 -.74784031498626095398, .56472001151566251186, .14842464993721351203e-1,
949 -.41162920273003120936, .20243071230196532282, .23772054897172750436,
950 -.24963810923972235950, -.12116179938394678936, .24330535483519110663,
951 .47903849781124471359e-1, -.22133299683101224293,
952 -.20542915138527200983e-2, .19653465717678146728,
953 -.26818172626509178444e-1, -.17319122357631210944,
954 .45065391411065545445e-1, .15253391395444065941,
955 -.56543897711725408302e-1, -.13469154928743585367,
956 .63632471400208840155e-1, .11941684923913523817,
957 -.67828850207933293098e-1, -.10636309084510652670,
958 .70095786922999181504e-1, .95187373095150709082e-1, 0., 0., 0., 0., 0.,
959 0., .15625000000000000000e-1, -.10909562534194485289,
960 .34842348626527747318, -.64461114561628111443, .69382480527334683659,
961 -.29551102358528827763, -.25527584713978439819, .38878771718544715394,
962 -.82956185835347407489e-2, -.31183177761966943912, .12831420840372374767,
963 .22067618205599434368, -.17569196937129496961, -.14598057000132284135,
964 .18864406621763419484, .89921002550386645767e-1, -.18571835020187122114,
965 -.48967672227195481777e-1, .17584685670380332798,
966 .19267984545067426324e-1, -.16335437520503462738,
967 .22598055455032407594e-2, .15032800884170631129,
968 -.17883358353754640871e-1, -.13774837869432209951,
969 .29227555960587143675e-1, .12604194747513151053, 0., 0., 0., 0., 0., 0.,
970 0., .78125000000000000000e-2, -.62377810244809812496e-1,
971 .23080781467370883845, -.50841310636012325368, .69834547012574056043,
972 -.52572723156526459672, .11464215704954976471e-1, .38698869011491210342,
973 -.26125646622255207507, -.16951698812361607510, .29773875898928782269,
974 .20130501202570367491e-1, -.26332493149159310198,
975 .67734613690401207009e-1, .21207315477103762715, -.11541543390889415193,
976 -.16249634759782417533, .13885887405041735068, .11996491328010275427,
977 -.14810432001630926895, -.85177658352556243411e-1, .14918860659904380587,
978 .57317789510444151564e-1, -.14569827645586660151,
979 -.35213090145965327390e-1, .13975998126844578198, 0., 0., 0., 0., 0., 0.,
980 0., 0., .39062500000000000000e-2, -.35101954600803571207e-1,
981 .14761284084133737720, -.37655033076080192966, .62410290231517322776,
982 -.64335622317683389875, .28188168266139524244, .22488495672137010675,
983 -.39393811089283576186, .75184777995770096714e-1, .28472023119398293003,
984 -.20410910833705899572, -.15590046962908511750, .23814567544617953125,
985 .54442805556829031204e-1, -.22855930338589720954,
986 .16303223615756629897e-1, .20172722433875559213,
987 -.62723406421217419404e-1, -.17012230831020922010,
988 .91754642766136561612e-1, .13927644821381121197, -.10886600968068418181,
989 -.11139075654373395292, .11797455976331702879, 0., 0., 0., 0., 0., 0., 0.,
990 0., 0., .19531250000000000000e-2, -.19506820659607596598e-1,
991 .91865676095362231937e-1, -.26604607809696493849, .51425874205091288223,
992 -.66047561132505329292, .48660109511591303851, -.17575661168678285615e-1,
993 -.36594333408055703366, .29088854695378694533, .11318677346656537927,
994 -.31110645235730182168, .60733219161008787341e-1, .24333848233620420826,
995 -.15254312332655419708, -.15995968483455388613, .19010344455215289289,
996 .86040636766440260000e-1, -.19652589954665259945,
997 -.27633388517205837713e-1, .18660848552712880387,
998 -.15942583868416775867e-1, -.16902042462382064786,
999 .47278526495327740646e-1, 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1000 .97656250000000000000e-3, -.10731084460857378207e-1,
1001 .55939644713816406331e-1, -.18118487371914493668, .39914857299829864263,
1002 -.60812322949933902435, .60011887183061967583, -.26002695805835928795,
1003 -.20883922404786010096, .38988130966114638081, -.11797833550782589082,
1004 -.25231824756239520077, .24817859972953934712, .90516417677868996417e-1,
1005 -.26079073291293066798, .30259468817169480161e-1, .22178195264114178432,
1006 -.10569877864302048175, -.16679648389266977455, .14637718550245050850,
1007 .11219272032739559870, -.16359363640525750353, -.64358194509092101393e-1,
1008 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., .48828125000000000000e-3,
1009 -.58542865274813470967e-2, .33461741635290096452e-1,
1010 -.11979993155896201271, .29580223766987206958, -.51874761979436016742,
1011 .62861483498014306968, -.44868895761051453296, .12567502628371529386e-1,
1012 .35040366183235474275, -.30466868455569500886, -.70903913601490112666e-1,
1013 .30822791893032512740, -.11969443264190207736, -.20764760317621313946,
1014 .20629838355452128532, .95269702915334718507e-1, -.22432624768705133300,
1015 -.33103381593477797101e-2, .20570036048155716333,
1016 -.62208282720094518964e-1, -.17095309330441436348, 0., 0., 0., 0., 0., 0.,
1017 0., 0., 0., 0., 0., 0., .24414062500000000000e-3,
1018 -.31714797501871532475e-2, .19721062526127334100e-1,
1019 -.77311181185536498246e-1, .21124871792841566575, -.41777980401893650886,
1020 .59401977834943551650, -.56132417807488349048, .23433675061367565951,
1021 .20222775295220942126, -.38280372496506190127, .14443804214023095767,
1022 .22268950939178466797, -.27211314150777981984, -.34184876506180717313e-1,
1023 .26006498895669734842, -.97650425186005090107e-1, -.19024527660129101293,
1024 .16789164198044635671, .10875811641651905252, -.19276785058805921298, 0.,
1025 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., .12207031250000000000e-3,
1026 -.17078941137247586143e-2, .11477733754843910060e-1,
1027 -.48887017020924625462e-1, .14634927241421789683, -.32156282683019547854,
1028 .52165811920227223937, -.60001958466396926460, .41208501541480733755,
1029 -.11366945503190350975e-2, -.33968093962672089159, .30955190935923386766,
1030 .40657421856578262210e-1, -.29873400409871531764, .16094481791768257440,
1031 .16876122436206497694, -.23650217045022161255, -.33070260090574765012e-1,
1032 .22985258456375907796, -.68645651043827097771e-1, 0., 0., 0., 0., 0., 0.,
1033 0., 0., 0., 0., 0., 0., 0., 0., .61035156250000000000e-4,
1034 -.91501857608428649078e-3, .66085179496951987952e-2,
1035 -.30383171695850355404e-1, .98840838845366876117e-1,
1036 -.23855447246420318989, .43322017468145613917, -.58049033744876107191,
1037 .52533893203742699346, -.20681056202371946180, -.20180000924562504384,
1038 .37503922291962681797, -.15988102869837429062, -.19823558102762374094,
1039 .28393023878803799622, -.11188133439357510403e-1, -.24730368377168229255,
1040 .14731529061377942839, .14878558042884266021, 0., 0., 0., 0., 0., 0., 0.,
1041 0., 0., 0., 0., 0., 0., 0., 0., .30517578125000000000e-4,
1042 -.48804277318479845551e-3, .37696080990601968396e-2,
1043 -.18603912108994738255e-1, .65325006755649582964e-1,
1044 -.17162960707938819795, .34411527956476971322, -.52289350347082497959,
1045 .57319653625674910592, -.37662253421045430413, -.14099055105384663902e-1,
1046 .33265570610216904208, -.30921265572647566661, -.19911390594166455281e-1,
1047 .28738590811031797718, -.18912130469738472647, -.13235936203215819193,
1048 .25076406142356675279, 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1049 0., 0., 0., .15258789062500000000e-4, -.25928719280954633249e-3,
1050 .21327398937568540428e-2, -.11244626133630732010e-1,
1051 .42375605740664331966e-1, -.12031130345907846211, .26352562258934426830,
1052 -.44590628258512682078, .56682835613700749379, -.49116715128261660395,
1053 .17845943097110339078, .20541650677432497477, -.36739803642257458221,
1054 .16776034069210108273, .17920950989905112908, -.28867732805385066532,
1055 .46473465543376206337e-1, 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1056 0., 0., 0., 0., 0., .76293945312500000000e-5, -.13727610943181290891e-3,
1057 .11979683091449349286e-2, -.67195313034570709806e-2,
1058 .27044920779931968175e-1, -.82472196498517457862e-1,
1059 .19570475044896150093, -.36391620788543817693, .52241392782736588032,
1060 -.54727504974907879912, .34211551468813581183, .31580472732719957762e-1,
1061 -.32830006549176759667, .30563797665254420769, .64905014620683140120e-2,
1062 -.27642986248995073032, 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1063 0., 0., 0., 0., 0., 0., .38146972656250000000e-5,
1064 -.72454147007837596854e-4, .66859847582761390285e-3,
1065 -.39751311980366118437e-2, .17015198650201528366e-1,
1066 -.55443621868993855715e-1, .14157060481641692131, -.28641242619559616836,
1067 .45610665490966615415, -.55262786406029265394, .45818352706035500108,
1068 -.14984403004611673047, -.21163807462970713245, .36007252928843413718,
1069 -.17030961385712954159, 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1070 0., 0., 0., 0., 0., 0., 0., .19073486328125000000e-5,
1071 -.38135049864067468562e-4, .37101393638555730015e-3,
1072 -.23305339886279723213e-2, .10569913448297127219e-1,
1073 -.36640175162216897547e-1, .10010476414320235508, -.21860074212675559892,
1074 .38124757096345313719, -.52020999209879669177, .52172632730659212045,
1075 -.30841620620308814614, -.50322546186721500184e-1, .32577618885114899053,
1076 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1077 0., 0., .95367431640625000000e-6, -.20021483206955925244e-4,
1078 .20481807322420625431e-3, -.13553476938058909882e-2,
1079 .64919676350791905019e-2, -.23848725425069251903e-1,
1080 .69384632678886421292e-1, -.16249711393618776934, .30736618106830314788,
1081 -.46399909601971539157, .53765031034002467225, -.42598991476520183929,
1082 .12130445348350215652, 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1083 0., 0., 0., 0., 0., 0., 0., 0., .47683715820312500000e-6,
1084 -.10487707828484902486e-4, .11254146162337528943e-3,
1085 -.78248929534271987118e-3, .39468337145306794566e-2,
1086 -.15313546659475671763e-1, .47249070825218564146e-1,
1087 -.11804374107101480543, .24031796927792491122, -.39629215049166341285,
1088 .51629108968402548545, -.49622372075429782915, 0., 0., 0., 0., 0., 0., 0.,
1089 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1090 .23841857910156250000e-6, -.54823314130625337326e-5,
1091 .61575377321535518154e-4, -.44877834366497538134e-3,
1092 .23774612048621955857e-2, -.97136347645161687796e-2,
1093 .31671599547606636717e-1, -.84028665767000747480e-1,
1094 .18298487576742964949, -.32647878537696945218, .46970971486488895077, 0.,
1095 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1096 0., 0., 0., 0., .11920928955078125000e-6, -.28604020001177375838e-5,
1097 .33559227978295551013e-4, -.25583821662860610560e-3,
1098 .14201552747787302339e-2, -.60938046986874414969e-2,
1099 .20930869247951926793e-1, -.58745021125678072911e-1,
1100 .13613725780285953720, -.26083988356030237586, 0., 0., 0., 0., 0., 0., 0.,
1101 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1102 .59604644775390625000e-7, -.14898180663526043291e-5,
1103 .18224991282807693921e-4, -.14504433444608833821e-3,
1104 .84184722720281809548e-3, -.37846965430000478789e-2,
1105 .13656355548211376864e-1, -.40409541997718853934e-1,
1106 .99226988101858325902e-1, 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1107 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1108 .29802322387695312500e-7, -.77471708843445529468e-6,
1109 .98649879372606876995e-5, -.81814934772838523887e-4,
1110 .49554483992403011328e-3, -.23290922072351413938e-2,
1111 .88068134250844034186e-2, -.27393666952485719070e-1, 0., 0., 0., 0., 0.,
1112 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1113 0., 0., 0., .14901161193847656250e-7, -.40226235946098233685e-6,
1114 .53236418690561306700e-5, -.45933829691164002269e-4,
1115 .28982005232838857913e-3, -.14212974043211018374e-2,
1116 .56192363087488842264e-2, 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1117 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1118 .74505805969238281250e-8, -.20858299254133430408e-6,
1119 .28648457300134381744e-5, -.25677535898258910850e-4,
1120 .16849420429491355445e-3, -.86062824010315834002e-3, 0., 0., 0., 0., 0.,
1121 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1122 0., 0., 0., 0., 0., .37252902984619140625e-8, -.10801736017613096861e-6,
1123 .15376606719887104015e-5, -.14296523739727437959e-4,
1124 .97419023656050887203e-4, 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1125 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1126 .18626451492309570312e-8, -.55871592916438890146e-7,
1127 .82331193828137454068e-6, -.79302250528382787666e-5, 0., 0., 0., 0., 0.,
1128 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1129 0., 0., 0., 0., 0., 0., 0., .93132257461547851562e-9,
1130 -.28867244235852488244e-7, .43982811713864556957e-6, 0., 0., 0., 0., 0.,
1131 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1132 0., 0., 0., 0., 0., 0., 0., 0., .46566128730773925781e-9,
1133 -.14899342093408253335e-7, 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1134 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1135 0., 0., .23283064365386962891e-9
1136 };
1137
1138 static const double Tright[33 * 33] = {
1139 1., .86602540378443864678, 0., -.33071891388307382381, 0.,
1140 .20728904939721249057, 0., -.15128841196122722208, 0.,
1141 .11918864298744029244, 0., -.98352013661686631224e-1, 0.,
1142 .83727065404940845733e-1, 0., -.72893399403505841203e-1, 0.,
1143 .64544632643375022436e-1, 0., -.57913170372415565639e-1, 0.,
1144 .52518242575729562263e-1, 0., -.48043311993977520457e-1, 0.,
1145 .44271433659733990243e-1, 0., -.41048928022856771981e-1, 0.,
1146 .38263878662008271459e-1, 0., -.35832844026365304501e-1, 0., 0.,
1147 .50000000000000000000, .96824583655185422130, .57282196186948000082,
1148 -.21650635094610966169, -.35903516540862679125, .97578093724974971969e-1,
1149 .26203921611325660506, -.55792409597991015609e-1, -.20644078533943456204,
1150 .36172381205961199479e-1, .17035068468874958194,
1151 -.25371838001497225980e-1, -.14501953125000000000,
1152 .18786835250972344757e-1, .12625507130328301066,
1153 -.14473795929590520582e-1, -.11179458309419422675,
1154 .11494434254897626155e-1, .10030855351241635862,
1155 -.93498556820544479096e-2, -.90964264465390582629e-1,
1156 .77546391824364392762e-2, .83213457337452292745e-1,
1157 -.65358085945588638605e-2, -.76680372422574234569e-1,
1158 .55835321940047427169e-2, .71098828931825789428e-1,
1159 -.48253327982967591019e-2, -.66274981937248958553e-1,
1160 .42118078245337801387e-2, .62064306433355646267e-1,
1161 -.37083386598903548973e-2, 0., 0., .25000000000000000000,
1162 .73950997288745200531, .83852549156242113615, .23175620272173946716,
1163 -.37791833195149451496, -.25710129174850522325, .21608307321780204633,
1164 .22844049245646009157, -.14009503000335388415, -.19897685605518413847,
1165 .98264706042471226893e-1, .17445445004279014046,
1166 -.72761100054958328401e-1, -.15463589893742108388,
1167 .56056770591708784481e-1, .13855313872640495158,
1168 -.44517752443294564781e-1, -.12534277657695128850,
1169 .36211835346039665762e-1, .11434398255136139683,
1170 -.30033588409423828125e-1, -.10506705408753910481,
1171 .25313077840725783008e-1, .97149327637744872155e-1,
1172 -.21624927200393328444e-1, -.90319582367202122625e-1,
1173 .18688433567711780666e-1, .84372291635345108584e-1,
1174 -.16312261561845420752e-1, -.79149526894804751586e-1,
1175 .14362333871852474757e-1, 0., 0., 0., .12500000000000000000,
1176 .49607837082461073572, .82265291131801144317, .59621200088559103072,
1177 -.80054302859059362371e-1, -.42612156697795759420,
1178 -.90098145270865592887e-1, .29769623255090078484, .13630307904779758221,
1179 -.21638835185708931831, -.14600247270306082052, .16348801804014290453,
1180 .14340708728599057249, -.12755243353979286190, -.13661523715071346961,
1181 .10215585947881057394, .12864248070157166547, -.83592528025348693602e-1,
1182 -.12066728689302565222, .69633728678718053052e-1, .11314245177331919532,
1183 -.58882939251410088028e-1, -.10621835858758221487,
1184 .50432266865187597572e-1, .99916834723527771581e-1,
1185 -.43672094283057258509e-1, -.94206380251950852413e-1,
1186 .38181356812697746418e-1, .89035739656537771225e-1,
1187 -.33661934598216332678e-1, 0., 0., 0., 0., .62500000000000000000e-1,
1188 .31093357409581873586, .67604086414949799246, .75644205980613611039,
1189 .28990586430124175741, -.30648508196770360914, -.35801372616842500052,
1190 .91326869828709014708e-1, .31127929687500000000, .90915752838698393094e-2,
1191 -.25637381283965534330, -.57601077850322797594e-1, .21019685709225757945,
1192 .81244992138514014256e-1, -.17375078516720988858,
1193 -.92289437277967051125e-1, .14527351914265391374,
1194 .96675340792832019889e-1, -.12289485697108543415,
1195 -.97448175340011084006e-1, .10511755943298339844,
1196 .96242247086378239657e-1, -.90822942272780513537e-1,
1197 -.93966350452322132384e-1, .79189411876493712558e-1,
1198 .91139307067989309325e-1, -.69613039934383197265e-1,
1199 -.88062491671135767870e-1, .61646331729340817494e-1, 0., 0., 0., 0., 0.,
1200 .31250000000000000000e-1, .18684782411095934408, .50176689760410660236,
1201 .74784031498626095398, .56472001151566251186, -.14842464993721351203e-1,
1202 -.41162920273003120936, -.20243071230196532282, .23772054897172750436,
1203 .24963810923972235950, -.12116179938394678936, -.24330535483519110663,
1204 .47903849781124471359e-1, .22133299683101224293,
1205 -.20542915138527200983e-2, -.19653465717678146728,
1206 -.26818172626509178444e-1, .17319122357631210944,
1207 .45065391411065545445e-1, -.15253391395444065941,
1208 -.56543897711725408302e-1, .13469154928743585367,
1209 .63632471400208840155e-1, -.11941684923913523817,
1210 -.67828850207933293098e-1, .10636309084510652670,
1211 .70095786922999181504e-1, -.95187373095150709082e-1, 0., 0., 0., 0., 0.,
1212 0., .15625000000000000000e-1, .10909562534194485289,
1213 .34842348626527747318, .64461114561628111443, .69382480527334683659,
1214 .29551102358528827763, -.25527584713978439819, -.38878771718544715394,
1215 -.82956185835347407489e-2, .31183177761966943912, .12831420840372374767,
1216 -.22067618205599434368, -.17569196937129496961, .14598057000132284135,
1217 .18864406621763419484, -.89921002550386645767e-1, -.18571835020187122114,
1218 .48967672227195481777e-1, .17584685670380332798,
1219 -.19267984545067426324e-1, -.16335437520503462738,
1220 -.22598055455032407594e-2, .15032800884170631129,
1221 .17883358353754640871e-1, -.13774837869432209951,
1222 -.29227555960587143675e-1, .12604194747513151053, 0., 0., 0., 0., 0., 0.,
1223 0., .78125000000000000000e-2, .62377810244809812496e-1,
1224 .23080781467370883845, .50841310636012325368, .69834547012574056043,
1225 .52572723156526459672, .11464215704954976471e-1, -.38698869011491210342,
1226 -.26125646622255207507, .16951698812361607510, .29773875898928782269,
1227 -.20130501202570367491e-1, -.26332493149159310198,
1228 -.67734613690401207009e-1, .21207315477103762715, .11541543390889415193,
1229 -.16249634759782417533, -.13885887405041735068, .11996491328010275427,
1230 .14810432001630926895, -.85177658352556243411e-1, -.14918860659904380587,
1231 .57317789510444151564e-1, .14569827645586660151,
1232 -.35213090145965327390e-1, -.13975998126844578198, 0., 0., 0., 0., 0., 0.,
1233 0., 0., .39062500000000000000e-2, .35101954600803571207e-1,
1234 .14761284084133737720, .37655033076080192966, .62410290231517322776,
1235 .64335622317683389875, .28188168266139524244, -.22488495672137010675,
1236 -.39393811089283576186, -.75184777995770096714e-1, .28472023119398293003,
1237 .20410910833705899572, -.15590046962908511750, -.23814567544617953125,
1238 .54442805556829031204e-1, .22855930338589720954, .16303223615756629897e-1,
1239 -.20172722433875559213, -.62723406421217419404e-1, .17012230831020922010,
1240 .91754642766136561612e-1, -.13927644821381121197, -.10886600968068418181,
1241 .11139075654373395292, .11797455976331702879, 0., 0., 0., 0., 0., 0., 0.,
1242 0., 0., .19531250000000000000e-2, .19506820659607596598e-1,
1243 .91865676095362231937e-1, .26604607809696493849, .51425874205091288223,
1244 .66047561132505329292, .48660109511591303851, .17575661168678285615e-1,
1245 -.36594333408055703366, -.29088854695378694533, .11318677346656537927,
1246 .31110645235730182168, .60733219161008787341e-1, -.24333848233620420826,
1247 -.15254312332655419708, .15995968483455388613, .19010344455215289289,
1248 -.86040636766440260000e-1, -.19652589954665259945,
1249 .27633388517205837713e-1, .18660848552712880387, .15942583868416775867e-1,
1250 -.16902042462382064786, -.47278526495327740646e-1, 0., 0., 0., 0., 0., 0.,
1251 0., 0., 0., 0., .97656250000000000000e-3, .10731084460857378207e-1,
1252 .55939644713816406331e-1, .18118487371914493668, .39914857299829864263,
1253 .60812322949933902435, .60011887183061967583, .26002695805835928795,
1254 -.20883922404786010096, -.38988130966114638081, -.11797833550782589082,
1255 .25231824756239520077, .24817859972953934712, -.90516417677868996417e-1,
1256 -.26079073291293066798, -.30259468817169480161e-1, .22178195264114178432,
1257 .10569877864302048175, -.16679648389266977455, -.14637718550245050850,
1258 .11219272032739559870, .16359363640525750353, -.64358194509092101393e-1,
1259 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., .48828125000000000000e-3,
1260 .58542865274813470967e-2, .33461741635290096452e-1, .11979993155896201271,
1261 .29580223766987206958, .51874761979436016742, .62861483498014306968,
1262 .44868895761051453296, .12567502628371529386e-1, -.35040366183235474275,
1263 -.30466868455569500886, .70903913601490112666e-1, .30822791893032512740,
1264 .11969443264190207736, -.20764760317621313946, -.20629838355452128532,
1265 .95269702915334718507e-1, .22432624768705133300,
1266 -.33103381593477797101e-2, -.20570036048155716333,
1267 -.62208282720094518964e-1, .17095309330441436348, 0., 0., 0., 0., 0., 0.,
1268 0., 0., 0., 0., 0., 0., .24414062500000000000e-3,
1269 .31714797501871532475e-2, .19721062526127334100e-1,
1270 .77311181185536498246e-1, .21124871792841566575, .41777980401893650886,
1271 .59401977834943551650, .56132417807488349048, .23433675061367565951,
1272 -.20222775295220942126, -.38280372496506190127, -.14443804214023095767,
1273 .22268950939178466797, .27211314150777981984, -.34184876506180717313e-1,
1274 -.26006498895669734842, -.97650425186005090107e-1, .19024527660129101293,
1275 .16789164198044635671, -.10875811641651905252, -.19276785058805921298, 0.,
1276 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., .12207031250000000000e-3,
1277 .17078941137247586143e-2, .11477733754843910060e-1,
1278 .48887017020924625462e-1, .14634927241421789683, .32156282683019547854,
1279 .52165811920227223937, .60001958466396926460, .41208501541480733755,
1280 .11366945503190350975e-2, -.33968093962672089159, -.30955190935923386766,
1281 .40657421856578262210e-1, .29873400409871531764, .16094481791768257440,
1282 -.16876122436206497694, -.23650217045022161255, .33070260090574765012e-1,
1283 .22985258456375907796, .68645651043827097771e-1, 0., 0., 0., 0., 0., 0.,
1284 0., 0., 0., 0., 0., 0., 0., 0., .61035156250000000000e-4,
1285 .91501857608428649078e-3, .66085179496951987952e-2,
1286 .30383171695850355404e-1, .98840838845366876117e-1, .23855447246420318989,
1287 .43322017468145613917, .58049033744876107191, .52533893203742699346,
1288 .20681056202371946180, -.20180000924562504384, -.37503922291962681797,
1289 -.15988102869837429062, .19823558102762374094, .28393023878803799622,
1290 .11188133439357510403e-1, -.24730368377168229255, -.14731529061377942839,
1291 .14878558042884266021, 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1292 0., 0., .30517578125000000000e-4, .48804277318479845551e-3,
1293 .37696080990601968396e-2, .18603912108994738255e-1,
1294 .65325006755649582964e-1, .17162960707938819795, .34411527956476971322,
1295 .52289350347082497959, .57319653625674910592, .37662253421045430413,
1296 -.14099055105384663902e-1, -.33265570610216904208, -.30921265572647566661,
1297 .19911390594166455281e-1, .28738590811031797718, .18912130469738472647,
1298 -.13235936203215819193, -.25076406142356675279, 0., 0., 0., 0., 0., 0.,
1299 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., .15258789062500000000e-4,
1300 .25928719280954633249e-3, .21327398937568540428e-2,
1301 .11244626133630732010e-1, .42375605740664331966e-1, .12031130345907846211,
1302 .26352562258934426830, .44590628258512682078, .56682835613700749379,
1303 .49116715128261660395, .17845943097110339078, -.20541650677432497477,
1304 -.36739803642257458221, -.16776034069210108273, .17920950989905112908,
1305 .28867732805385066532, .46473465543376206337e-1, 0., 0., 0., 0., 0., 0.,
1306 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., .76293945312500000000e-5,
1307 .13727610943181290891e-3, .11979683091449349286e-2,
1308 .67195313034570709806e-2, .27044920779931968175e-1,
1309 .82472196498517457862e-1, .19570475044896150093, .36391620788543817693,
1310 .52241392782736588032, .54727504974907879912, .34211551468813581183,
1311 -.31580472732719957762e-1, -.32830006549176759667, -.30563797665254420769,
1312 .64905014620683140120e-2, .27642986248995073032, 0., 0., 0., 0., 0., 0.,
1313 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., .38146972656250000000e-5,
1314 .72454147007837596854e-4, .66859847582761390285e-3,
1315 .39751311980366118437e-2, .17015198650201528366e-1,
1316 .55443621868993855715e-1, .14157060481641692131, .28641242619559616836,
1317 .45610665490966615415, .55262786406029265394, .45818352706035500108,
1318 .14984403004611673047, -.21163807462970713245, -.36007252928843413718,
1319 -.17030961385712954159, 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1320 0., 0., 0., 0., 0., 0., 0., .19073486328125000000e-5,
1321 .38135049864067468562e-4, .37101393638555730015e-3,
1322 .23305339886279723213e-2, .10569913448297127219e-1,
1323 .36640175162216897547e-1, .10010476414320235508, .21860074212675559892,
1324 .38124757096345313719, .52020999209879669177, .52172632730659212045,
1325 .30841620620308814614, -.50322546186721500184e-1, -.32577618885114899053,
1326 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1327 0., 0., .95367431640625000000e-6, .20021483206955925244e-4,
1328 .20481807322420625431e-3, .13553476938058909882e-2,
1329 .64919676350791905019e-2, .23848725425069251903e-1,
1330 .69384632678886421292e-1, .16249711393618776934, .30736618106830314788,
1331 .46399909601971539157, .53765031034002467225, .42598991476520183929,
1332 .12130445348350215652, 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1333 0., 0., 0., 0., 0., 0., 0., 0., .47683715820312500000e-6,
1334 .10487707828484902486e-4, .11254146162337528943e-3,
1335 .78248929534271987118e-3, .39468337145306794566e-2,
1336 .15313546659475671763e-1, .47249070825218564146e-1, .11804374107101480543,
1337 .24031796927792491122, .39629215049166341285, .51629108968402548545,
1338 .49622372075429782915, 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1339 0., 0., 0., 0., 0., 0., 0., 0., 0., .23841857910156250000e-6,
1340 .54823314130625337326e-5, .61575377321535518154e-4,
1341 .44877834366497538134e-3, .23774612048621955857e-2,
1342 .97136347645161687796e-2, .31671599547606636717e-1,
1343 .84028665767000747480e-1, .18298487576742964949, .32647878537696945218,
1344 .46970971486488895077, 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1345 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., .11920928955078125000e-6,
1346 .28604020001177375838e-5, .33559227978295551013e-4,
1347 .25583821662860610560e-3, .14201552747787302339e-2,
1348 .60938046986874414969e-2, .20930869247951926793e-1,
1349 .58745021125678072911e-1, .13613725780285953720, .26083988356030237586,
1350 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1351 0., 0., 0., 0., 0., 0., .59604644775390625000e-7,
1352 .14898180663526043291e-5, .18224991282807693921e-4,
1353 .14504433444608833821e-3, .84184722720281809548e-3,
1354 .37846965430000478789e-2, .13656355548211376864e-1,
1355 .40409541997718853934e-1, .99226988101858325902e-1, 0., 0., 0., 0., 0.,
1356 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1357 0., 0., .29802322387695312500e-7, .77471708843445529468e-6,
1358 .98649879372606876995e-5, .81814934772838523887e-4,
1359 .49554483992403011328e-3, .23290922072351413938e-2,
1360 .88068134250844034186e-2, .27393666952485719070e-1, 0., 0., 0., 0., 0.,
1361 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1362 0., 0., 0., .14901161193847656250e-7, .40226235946098233685e-6,
1363 .53236418690561306700e-5, .45933829691164002269e-4,
1364 .28982005232838857913e-3, .14212974043211018374e-2,
1365 .56192363087488842264e-2, 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1366 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1367 .74505805969238281250e-8, .20858299254133430408e-6,
1368 .28648457300134381744e-5, .25677535898258910850e-4,
1369 .16849420429491355445e-3, .86062824010315834002e-3, 0., 0., 0., 0., 0.,
1370 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1371 0., 0., 0., 0., 0., .37252902984619140625e-8, .10801736017613096861e-6,
1372 .15376606719887104015e-5, .14296523739727437959e-4,
1373 .97419023656050887203e-4, 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1374 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1375 .18626451492309570312e-8, .55871592916438890146e-7,
1376 .82331193828137454068e-6, .79302250528382787666e-5, 0., 0., 0., 0., 0.,
1377 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1378 0., 0., 0., 0., 0., 0., 0., .93132257461547851562e-9,
1379 .28867244235852488244e-7, .43982811713864556957e-6, 0., 0., 0., 0., 0.,
1380 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1381 0., 0., 0., 0., 0., 0., 0., 0., .46566128730773925781e-9,
1382 .14899342093408253335e-7, 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1383 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
1384 0., 0., .23283064365386962891e-9
1385 };
1386
1387 /* Allocates a workspace for the given maximum number of intervals.
1388 Note that if the workspace gets filled, the intervals with the
1389 lowest error estimates are dropped. The maximum number of
1390 intervals is therefore not the maximum number of intervals
1391 that will be computed, but merely the size of the buffer.
1392 */
1393
1394 /* Compute the product of the fx with one of the inverse
1395 Vandermonde-like matrices. */
1396
1397 void
1398 Vinvfx (const double *fx, double *c, const int d)
1399 {
1400
1401 int i, j;
1402
1403 switch (d)
1404 {
1405 case 0:
1406 for (i = 0; i <= 4; i++)
1407 {
1408 c[i] = 0.0;
1409 for (j = 0; j <= 4; j++)
1410 c[i] += V1inv[i * 5 + j] * fx[j * 8];
1411 }
1412 break;
1413 case 1:
1414 for (i = 0; i <= 8; i++)
1415 {
1416 c[i] = 0.0;
1417 for (j = 0; j <= 8; j++)
1418 c[i] += V2inv[i * 9 + j] * fx[j * 4];
1419 }
1420 break;
1421 case 2:
1422 for (i = 0; i <= 16; i++)
1423 {
1424 c[i] = 0.0;
1425 for (j = 0; j <= 16; j++)
1426 c[i] += V3inv[i * 17 + j] * fx[j * 2];
1427 }
1428 break;
1429 case 3:
1430 for (i = 0; i <= 32; i++)
1431 {
1432 c[i] = 0.0;
1433 for (j = 0; j <= 32; j++)
1434 c[i] += V4inv[i * 33 + j] * fx[j];
1435 }
1436 break;
1437 }
1438
1439 }
1440
1441
1442 /* Downdate the interpolation given by the n coefficients c
1443 by removing the nodes with indices in nans. */
1444
1445 void
1446 downdate (double *c, int n, int d, int *nans, int nnans)
1447 {
1448
1449 static const int bidx[4] = { 0, 6, 16, 34 };
1450 double b_new[34], alpha;
1451 int i, j;
1452
1453 for (i = 0; i <= n + 1; i++)
1454 b_new[i] = bee[bidx[d] + i];
1455 for (i = 0; i < nnans; i++)
1456 {
1457 b_new[n + 1] = b_new[n + 1] / Lalpha[n];
1458 b_new[n] = (b_new[n] + xi[nans[i]] * b_new[n + 1]) / Lalpha[n - 1];
1459 for (j = n - 1; j > 0; j--)
1460 b_new[j] =
1461 (b_new[j] + xi[nans[i]] * b_new[j + 1] -
1462 Lgamma[j + 1] * b_new[j + 2]) / Lalpha[j - 1];
1463 for (j = 0; j <= n; j++)
1464 b_new[j] = b_new[j + 1];
1465 alpha = c[n] / b_new[n];
1466 for (j = 0; j < n; j++)
1467 c[j] -= alpha * b_new[j];
1468 c[n] = 0;
1469 n--;
1470 }
1471
1472 }
1473
1474
1475 /* The actual integration routine. */
1476
1477 DEFUN (quadcc, args, nargout,
1478 "-*- texinfo -*-\n\
1479 @deftypefn {Function File} {@var{q} =} quadcc (@var{f}, @var{a}, @var{b})\n\
1480 @deftypefnx {Function File} {@var{q} =} quadcc (@var{f}, @var{a}, @var{b}, @var{tol})\n\
1481 @deftypefnx {Function File} {@var{q} =} quadcc (@var{f}, @var{a}, @var{b}, @var{tol}, @var{sing})\n\
1482 @deftypefnx {Function File} {[@var{q}, @var{err}, @var{nr_points}] =} quadcc (@dots{})\n\
1483 Numerically evaluate the integral of @var{f} from @var{a} to @var{b}\n\
1484 using the doubly-adaptive Clenshaw-Curtis quadrature described by P. Gonnet\n\
1485 in @cite{Increasing the Reliability of Adaptive Quadrature Using Explicit\n\
1486 Interpolants}.\n\
1487 @var{f} is a function handle, inline function, or string\n\
1488 containing the name of the function to evaluate.\n\
1489 The function @var{f} must be vectorized and must return a vector of output\n\
1490 values if given a vector of input values. For example,\n\
1491 \n\
1492 @example\n\
1493 f = @@(x) x .* sin (1./x) .* sqrt (abs (1 - x));\n\
1494 @end example\n\
1495 \n\
1496 @noindent\n\
1497 which uses the element-by-element `dot' form for all operators.\n\
1498 \n\
1499 @var{a} and @var{b} are the lower and upper limits of integration. Either\n\
1500 or both limits may be infinite. @code{quadcc} handles an inifinite limit\n\
1501 by substituting the variable of integration with @code{x = tan (pi/2*u)}.\n\
1502 \n\
1503 The optional argument @var{tol} defines the relative tolerance used to stop\n\
1504 the integration procedure. The default value is @math{1e^{-6}}.\n\
1505 \n\
1506 The optional argument @var{sing} contains a list of points where the\n\
1507 integrand has known singularities, or discontinuities\n\
1508 in any of its derivatives, inside the integration interval.\n\
1509 For the example above, which has a discontinuity at x=1, the call to\n\
1510 @code{quadcc} would be as follows\n\
1511 \n\
1512 @example\n\
1513 int = quadcc (f, a, b, 1.0e-6, [ 1 ]);\n\
1514 @end example\n\
1515 \n\
1516 The result of the integration is returned in @var{q}.\n\
1517 @var{err} is an estimate of the absolute integration error and\n\
1518 @var{nr_points} is the number of points at which the integrand was evaluated.\n\
1519 If the adaptive integration did not converge, the value of\n\
1520 @var{err} will be larger than the requested tolerance. Therefore, it is\n\
1521 recommended to verify this value for difficult integrands.\n\
1522 \n\
1523 @code{quadcc} is capable of dealing with non-numeric\n\
1524 values of the integrand such as @code{NaN} or @code{Inf}.\n\
1525 If the integral diverges, and @code{quadcc} detects this,\n\
1526 then a warning is issued and @code{Inf} or @code{-Inf} is returned.\n\
1527 \n\
1528 Note: @code{quadcc} is a general purpose quadrature algorithm\n\
1529 and, as such, may be less efficient for a smooth or otherwise\n\
1530 well-behaved integrand than other methods such as @code{quadgk}.\n\
1531 \n\
1532 The algorithm uses Clenshaw-Curtis quadrature rules of increasing\n\
1533 degree in each interval and bisects the interval if either the\n\
1534 function does not appear to be smooth or a rule of maximum\n\
1535 degree has been reached. The error estimate is computed from the\n\
1536 L2-norm of the difference between two successive interpolations\n\
1537 of the integrand over the nodes of the respective quadrature rules.\n\
1538 \n\
1539 Reference: P. Gonnet, @cite{Increasing the Reliability of Adaptive\n\
1540 Quadrature Using Explicit Interpolants}, ACM Transactions on\n\
1541 Mathematical Software, Vol. 37, Issue 3, Article No. 3, 2010.\n\
1542 @seealso{quad, quadv, quadl, quadgk, trapz, dblquad, triplequad}\n\
1543 @end deftypefn")
1544 {
1545 octave_value_list retval;
1546
1547 /* Some constants that we will need. */
1548 static const int n[4] = { 4, 8, 16, 32 };
1549 static const int skip[4] = { 8, 4, 2, 1 };
1550 static const int idx[4] = { 0, 5, 14, 31 };
1551 static const double w = M_SQRT2 / 2;
1552 static const int ndiv_max = 20;
1553
1554 /* The interval heap. */
1555 cquad_ival ivals[cquad_heapsize];
1556 int heap[cquad_heapsize];
1557
1558 /* Arguments left and right */
1559 int nargin = args.length ();
1560 octave_function *fcn;
1561 double a, b, tol, iivals[cquad_heapsize], *sing;
1562
1563 /* Variables needed for transforming the integrand. */
1564 bool wrap = false;
1565 double xw;
1566
1567 /* Stuff we will need to call the integrand. */
1568 octave_value_list fargs, fvals;
1569
1570 /* Actual variables (as opposed to constants above). */
1571 double m, h, ml, hl, mr, hr, temp;
1572 double igral, err, igral_final, err_final;
1573 int nivals, neval = 0;
1574 int i, j, d, split, t;
1575 int nnans, nans[33];
1576 cquad_ival *iv, *ivl, *ivr;
1577 double nc, ncdiff;
1578
1579
1580 /* Parse the input arguments. */
1581 if (nargin < 3)
1582 {
1583 print_usage ();
1584 return retval;
1585 }
1586
1587 if (args(0).is_function_handle () || args(0).is_inline_function ())
1588 fcn = args(0).function_value ();
1589 else
1590 {
1591 std::string fcn_name = unique_symbol_name ("__quadcc_fcn_");
1592 std::string fname = "function y = ";
1593 fname.append (fcn_name);
1594 fname.append ("(x) y = ");
1595 fcn = extract_function (args(0), "quadcc", fcn_name, fname,
1596 "; endfunction");
1597 }
1598
1599 if (!args(1).is_real_scalar ())
1600 {
1601 error ("quadcc: lower limit of integration (A) must be a single real scalar");
1602 return retval;
1603 }
1604 else
1605 a = args(1).double_value ();
1606
1607 if (!args(2).is_real_scalar ())
1608 {
1609 error ("quadcc: upper limit of integration (B) must be a single real scalar");
1610 return retval;
1611 }
1612 else
1613 b = args(2).double_value ();
1614
1615 if (nargin < 4 || args(3).is_empty ())
1616 tol = 1.0e-6;
1617 else if (!args(3).is_real_scalar () || args(3).double_value () <= 0)
1618 {
1619 error ("quadcc: tolerance (TOL) must be a single real scalar > 0");
1620 return retval;
1621 }
1622 else
1623 tol = args(3).double_value ();
1624
1625 if (nargin < 5)
1626 {
1627 nivals = 1;
1628 iivals[0] = a;
1629 iivals[1] = b;
1630 }
1631 else if (!(args(4).is_real_scalar () || args(4).is_real_matrix ()))
1632 {
1633 error ("quadcc: list of singularities (SING) must be a vector of real values");
1634 return retval;
1635 }
1636 else
1637 {
1638 nivals = 1 + args(4).length ();
1639 if (nivals > cquad_heapsize)
1640 {
1641 error ("quadcc: maximum number of singular points is limited to %i",
1642 cquad_heapsize-1);
1643 return retval;
1644 }
1645 sing = args(4).array_value ().fortran_vec ();
1646 iivals[0] = a;
1647 for (i = 0; i < nivals - 2; i++)
1648 iivals[i + 1] = sing[i];
1649 iivals[nivals] = b;
1650 }
1651
1652 /* If a or b are +/-Inf, transform the integral. */
1653 if (xisinf (a) || xisinf (b))
1654 {
1655 wrap = true;
1656 for (i = 0; i <= nivals; i++)
1657 if (xisinf (iivals[i]))
1658 iivals[i] = gnulib::copysign (1.0, iivals[i]);
1659 else
1660 iivals[i] = 2.0 * atan (iivals[i]) / M_PI;
1661 }
1662
1663
1664 /* Initialize the heaps. */
1665 for (i = 0; i < cquad_heapsize; i++)
1666 heap[i] = i;
1667
1668
1669 /* Create the first interval(s). */
1670 igral = 0.0;
1671 err = 0.0;
1672 for (j = 0; j < nivals; j++)
1673 {
1674
1675 /* Initialize the interval. */
1676 iv = &(ivals[heap[j]]);
1677 m = (iivals[j] + iivals[j + 1]) / 2;
1678 h = (iivals[j + 1] - iivals[j]) / 2;
1679 nnans = 0;
1680 ColumnVector ex (33);
1681 if (wrap)
1682 {
1683 for (i = 0; i <= n[3]; i++)
1684 ex (i) = tan (M_PI / 2 * (m + xi[i] * h));
1685 }
1686 else
1687 {
1688 for (i = 0; i <= n[3]; i++)
1689 ex (i) = m + xi[i] * h;
1690 }
1691 fargs(0) = ex;
1692 fvals = feval (fcn, fargs, 1);
1693 if (fvals.length () != 1 || !fvals(0).is_real_matrix ())
1694 {
1695 error ("quadcc: integrand F must return a single, real-valued vector");
1696 return retval;
1697 }
1698 Matrix effex = fvals(0).matrix_value ();
1699 if (effex.length () != ex.length ())
1700 {
1701 error ("quadcc: integrand F must return a single, real-valued vector of the same size as the input");
1702 return retval;
1703 }
1704 for (i = 0; i <= n[3]; i++)
1705 {
1706 iv->fx[i] = effex (i);
1707 if (wrap)
1708 {
1709 xw = ex(i);
1710 iv->fx[i] *= (1.0 + xw * xw) * M_PI / 2;
1711 }
1712 neval++;
1713 if (!xfinite (iv->fx[i]))
1714 {
1715 nans[nnans++] = i;
1716 iv->fx[i] = 0.0;
1717 }
1718 }
1719 Vinvfx (iv->fx, &(iv->c[idx[3]]), 3);
1720 Vinvfx (iv->fx, &(iv->c[idx[2]]), 2);
1721 Vinvfx (iv->fx, &(iv->c[0]), 0);
1722 for (i = 0; i < nnans; i++)
1723 iv->fx[i] = octave_NaN;
1724 iv->a = iivals[j];
1725 iv->b = iivals[j + 1];
1726 iv->depth = 3;
1727 iv->rdepth = 1;
1728 iv->ndiv = 0;
1729 iv->igral = 2 * h * iv->c[idx[3]] * w;
1730 nc = 0.0;
1731 for (i = n[2] + 1; i <= n[3]; i++)
1732 {
1733 temp = iv->c[idx[3] + i];
1734 nc += temp * temp;
1735 }
1736 ncdiff = nc;
1737 for (i = 0; i <= n[2]; i++)
1738 {
1739 temp = iv->c[idx[2] + i] - iv->c[idx[3] + i];
1740 ncdiff += temp * temp;
1741 nc += iv->c[idx[3] + i] * iv->c[idx[3] + i];
1742 }
1743 ncdiff = sqrt (ncdiff);
1744 nc = sqrt (nc);
1745 iv->err = ncdiff * 2 * h;
1746 if (ncdiff / nc > 0.1 && iv->err < 2 * h * nc)
1747 iv->err = 2 * h * nc;
1748
1749 /* Tabulate this interval's data. */
1750 igral += iv->igral;
1751 err += iv->err;
1752
1753 /* Sift it up the heap. */
1754 i = j;
1755 while (i > 0 && ivals[heap[i / 2]].err < ivals[heap[i]].err)
1756 {
1757 temp = heap[i];
1758 heap[i] = heap[i / 2];
1759 heap[i / 2] = temp;
1760 i /= 2;
1761 }
1762
1763 }
1764
1765
1766 /* Initialize some global values. */
1767 igral_final = 0.0;
1768 err_final = 0.0;
1769
1770
1771 /* Main loop. */
1772 while (nivals > 0 && err > 0.0 && err > fabs (igral) * tol
1773 && !(err_final > fabs (igral) * tol
1774 && err - err_final < fabs (igral) * tol))
1775 {
1776
1777 /* Allow the user to interrupt. */
1778 OCTAVE_QUIT;
1779
1780 /* Put our finger on the interval with the largest error. */
1781 iv = &(ivals[heap[0]]);
1782 m = (iv->a + iv->b) / 2;
1783 h = (iv->b - iv->a) / 2;
1784
1785 /* printf
1786 ("quadcc: processing ival %i (of %i) with [%e,%e] int=%e, err=%e, depth=%i\n",
1787 heap[0], nivals, iv->a, iv->b, iv->igral, iv->err, iv->depth);
1788 */
1789 /* Should we try to increase the degree? */
1790 if (iv->depth < 3)
1791 {
1792
1793 /* Keep tabs on some variables. */
1794 d = ++iv->depth;
1795
1796 /* Get the new (missing) function values */
1797 {
1798 ColumnVector ex (n[d] / 2);
1799 if (wrap)
1800 {
1801 for (i = 0; i < n[d] / 2; i++)
1802 ex (i) =
1803 tan (M_PI / 2 * (m + xi[(2 * i + 1) * skip[d]] * h));
1804 }
1805 else
1806 {
1807 for (i = 0; i < n[d] / 2; i++)
1808 ex (i) = m + xi[(2 * i + 1) * skip[d]] * h;
1809 }
1810 fargs(0) = ex;
1811 fvals = feval (fcn, fargs, 1);
1812 if (fvals.length () != 1 || !fvals(0).is_real_matrix ())
1813 {
1814 error ("quadcc: integrand F must return a single, real-valued vector");
1815 return retval;
1816 }
1817 Matrix effex = fvals(0).matrix_value ();
1818 if (effex.length () != ex.length ())
1819 {
1820 error ("quadcc: integrand F must return a single, real-valued vector of the same size as the input");
1821 return retval;
1822 }
1823 neval += effex.length ();
1824 for (i = 0; i < n[d] / 2; i++)
1825 {
1826 j = (2 * i + 1) * skip[d];
1827 iv->fx[j] = effex (i);
1828 if (wrap)
1829 {
1830 xw = ex(i);
1831 iv->fx[j] *= (1.0 + xw * xw) * M_PI / 2;
1832 }
1833 }
1834 }
1835 nnans = 0;
1836 for (i = 0; i <= 32; i += skip[d])
1837 {
1838 if (!xfinite (iv->fx[i]))
1839 {
1840 nans[nnans++] = i;
1841 iv->fx[i] = 0.0;
1842 }
1843 }
1844
1845 /* Compute the new coefficients. */
1846 Vinvfx (iv->fx, &(iv->c[idx[d]]), d);
1847 /* Downdate any NaNs. */
1848 if (nnans > 0)
1849 {
1850 downdate (&(iv->c[idx[d]]), n[d], d, nans, nnans);
1851 for (i = 0; i < nnans; i++)
1852 iv->fx[i] = octave_NaN;
1853 }
1854
1855 /* Compute the error estimate. */
1856 nc = 0.0;
1857 for (i = n[d - 1] + 1; i <= n[d]; i++)
1858 {
1859 temp = iv->c[idx[d] + i];
1860 nc += temp * temp;
1861 }
1862 ncdiff = nc;
1863 for (i = 0; i <= n[d - 1]; i++)
1864 {
1865 temp = iv->c[idx[d - 1] + i] - iv->c[idx[d] + i];
1866 ncdiff += temp * temp;
1867 nc += iv->c[idx[d] + i] * iv->c[idx[d] + i];
1868 }
1869 ncdiff = sqrt (ncdiff);
1870 nc = sqrt (nc);
1871 iv->err = ncdiff * 2 * h;
1872 /* Compute the local integral. */
1873 iv->igral = 2 * h * w * iv->c[idx[d]];
1874 /* Split the interval prematurely? */
1875 split = (nc > 0 && ncdiff / nc > 0.1);
1876 }
1877
1878 /* Maximum degree reached, just split. */
1879 else
1880 {
1881 split = 1;
1882 }
1883
1884
1885 /* Should we drop this interval? */
1886 if ((m + h * xi[0]) >= (m + h * xi[1])
1887 || (m + h * xi[31]) >= (m + h * xi[32])
1888 || iv->err < fabs (iv->igral) * DBL_EPSILON * 10)
1889 {
1890
1891 /* printf
1892 ("quadcc: dropping ival %i (of %i) with [%e,%e] int=%e, err=%e, depth=%i\n",
1893 heap[0], nivals, iv->a, iv->b, iv->igral, iv->err,
1894 iv->depth);
1895 */
1896 /* Keep this interval's contribution */
1897 err_final += iv->err;
1898 igral_final += iv->igral;
1899 /* Swap with the last element on the heap */
1900 t = heap[nivals - 1];
1901 heap[nivals - 1] = heap[0];
1902 heap[0] = t;
1903 nivals--;
1904 /* Fix up the heap */
1905 i = 0;
1906 while (2 * i + 1 < nivals)
1907 {
1908
1909 /* Get the kids */
1910 j = 2 * i + 1;
1911 /* If the j+1st entry exists and is larger than the jth,
1912 use it instead. */
1913 if (j + 1 < nivals
1914 && ivals[heap[j + 1]].err >= ivals[heap[j]].err)
1915 j++;
1916 /* Do we need to move the ith entry up? */
1917 if (ivals[heap[j]].err <= ivals[heap[i]].err)
1918 break;
1919 else
1920 {
1921 t = heap[j];
1922 heap[j] = heap[i];
1923 heap[i] = t;
1924 i = j;
1925 }
1926 }
1927
1928 }
1929
1930 /* Do we need to split this interval? */
1931 else if (split)
1932 {
1933
1934 /* Some values we will need often... */
1935 d = iv->depth;
1936 /* Generate the interval on the left */
1937 ivl = &(ivals[heap[nivals++]]);
1938 ivl->a = iv->a;
1939 ivl->b = m;
1940 ml = (ivl->a + ivl->b) / 2;
1941 hl = h / 2;
1942 ivl->depth = 0;
1943 ivl->rdepth = iv->rdepth + 1;
1944 ivl->fx[0] = iv->fx[0];
1945 ivl->fx[32] = iv->fx[16];
1946 {
1947 ColumnVector ex (n[0] - 1);
1948 if (wrap)
1949 {
1950 for (i = 0; i < n[0] - 1; i++)
1951 ex (i) = tan (M_PI / 2 * (ml + xi[(i + 1) * skip[0]] * hl));
1952 }
1953 else
1954 {
1955 for (i = 0; i < n[0] - 1; i++)
1956 ex (i) = ml + xi[(i + 1) * skip[0]] * hl;
1957 }
1958 fargs(0) = ex;
1959 fvals = feval (fcn, fargs, 1);
1960 if (fvals.length () != 1 || !fvals(0).is_real_matrix ())
1961 {
1962 error ("quadcc: integrand F must return a single, real-valued vector");
1963 return retval;
1964 }
1965 Matrix effex = fvals(0).matrix_value ();
1966 if (effex.length () != ex.length ())
1967 {
1968 error ("quadcc: integrand F must return a single, real-valued vector of the same size as the input");
1969 return retval;
1970 }
1971 neval += effex.length ();
1972 for (i = 0; i < n[0] - 1; i++)
1973 {
1974 j = (i + 1) * skip[0];
1975 ivl->fx[j] = effex (i);
1976 if (wrap)
1977 {
1978 xw = ex(i);
1979 ivl->fx[j] *= (1.0 + xw * xw) * M_PI / 2;
1980 }
1981 }
1982 }
1983 nnans = 0;
1984 for (i = 0; i <= 32; i += skip[0])
1985 {
1986 if (!xfinite (ivl->fx[i]))
1987 {
1988 nans[nnans++] = i;
1989 ivl->fx[i] = 0.0;
1990 }
1991 }
1992 Vinvfx (ivl->fx, ivl->c, 0);
1993 if (nnans > 0)
1994 {
1995 downdate (ivl->c, n[0], 0, nans, nnans);
1996 for (i = 0; i < nnans; i++)
1997 ivl->fx[i] = octave_NaN;
1998 }
1999 for (i = 0; i <= n[d]; i++)
2000 {
2001 ivl->c[idx[d] + i] = 0.0;
2002 for (j = i; j <= n[d]; j++)
2003 ivl->c[idx[d] + i] += Tleft[i * 33 + j] * iv->c[idx[d] + j];
2004 }
2005 ncdiff = 0.0;
2006 for (i = 0; i <= n[0]; i++)
2007 {
2008 temp = ivl->c[i] - ivl->c[idx[d] + i];
2009 ncdiff += temp * temp;
2010 }
2011 for (i = n[0] + 1; i <= n[d]; i++)
2012 {
2013 temp = ivl->c[idx[d] + i];
2014 ncdiff += temp * temp;
2015 }
2016 ncdiff = sqrt (ncdiff);
2017 ivl->err = ncdiff * h;
2018 /* Check for divergence. */
2019 ivl->ndiv = iv->ndiv + (fabs (iv->c[0]) > 0
2020 && ivl->c[0] / iv->c[0] > 2);
2021 if (ivl->ndiv > ndiv_max && 2 * ivl->ndiv > ivl->rdepth)
2022 {
2023 igral = gnulib::copysign (octave_Inf, igral);
2024 warning ("quadcc: divergent integral detected");
2025 break;
2026 }
2027
2028 /* Compute the local integral. */
2029 ivl->igral = h * w * ivl->c[0];
2030
2031
2032 /* Generate the interval on the right */
2033 ivr = &(ivals[heap[nivals++]]);
2034 ivr->a = m;
2035 ivr->b = iv->b;
2036 mr = (ivr->a + ivr->b) / 2;
2037 hr = h / 2;
2038 ivr->depth = 0;
2039 ivr->rdepth = iv->rdepth + 1;
2040 ivr->fx[0] = iv->fx[16];
2041 ivr->fx[32] = iv->fx[32];
2042 {
2043 ColumnVector ex (n[0] - 1);
2044 if (wrap)
2045 {
2046 for (i = 0; i < n[0] - 1; i++)
2047 ex (i) = tan (M_PI / 2 * (mr + xi[(i + 1) * skip[0]] * hr));
2048 }
2049 else
2050 {
2051 for (i = 0; i < n[0] - 1; i++)
2052 ex (i) = mr + xi[(i + 1) * skip[0]] * hr;
2053 }
2054 fargs(0) = ex;
2055 fvals = feval (fcn, fargs, 1);
2056 if (fvals.length () != 1 || !fvals(0).is_real_matrix ())
2057 {
2058 error ("quadcc: integrand F must return a single, real-valued vector");
2059 return retval;
2060 }
2061 Matrix effex = fvals(0).matrix_value ();
2062 if (effex.length () != ex.length ())
2063 {
2064 error ("quadcc: integrand F must return a single, real-valued vector of the same size as the input");
2065 return retval;
2066 }
2067 neval += effex.length ();
2068 for (i = 0; i < n[0] - 1; i++)
2069 {
2070 j = (i + 1) * skip[0];
2071 ivr->fx[j] = effex (i);
2072 if (wrap)
2073 {
2074 xw = ex(i);
2075 ivr->fx[j] *= (1.0 + xw * xw) * M_PI / 2;
2076 }
2077 }
2078 }
2079 nnans = 0;
2080 for (i = 0; i <= 32; i += skip[0])
2081 {
2082 if (!xfinite (ivr->fx[i]))
2083 {
2084 nans[nnans++] = i;
2085 ivr->fx[i] = 0.0;
2086 }
2087 }
2088 Vinvfx (ivr->fx, ivr->c, 0);
2089 if (nnans > 0)
2090 {
2091 downdate (ivr->c, n[0], 0, nans, nnans);
2092 for (i = 0; i < nnans; i++)
2093 ivr->fx[i] = octave_NaN;
2094 }
2095 for (i = 0; i <= n[d]; i++)
2096 {
2097 ivr->c[idx[d] + i] = 0.0;
2098 for (j = i; j <= n[d]; j++)
2099 ivr->c[idx[d] + i] += Tright[i * 33 + j] * iv->c[idx[d] + j];
2100 }
2101 ncdiff = 0.0;
2102 for (i = 0; i <= n[0]; i++)
2103 {
2104 temp = ivr->c[i] - ivr->c[idx[d] + i];
2105 ncdiff += temp * temp;
2106 }
2107 for (i = n[0] + 1; i <= n[d]; i++)
2108 {
2109 temp = ivr->c[idx[d] + i];
2110 ncdiff += temp * temp;
2111 }
2112 ncdiff = sqrt (ncdiff);
2113 ivr->err = ncdiff * h;
2114 /* Check for divergence. */
2115 ivr->ndiv = iv->ndiv + (fabs (iv->c[0]) > 0
2116 && ivr->c[0] / iv->c[0] > 2);
2117 if (ivr->ndiv > ndiv_max && 2 * ivr->ndiv > ivr->rdepth)
2118 {
2119 igral = gnulib::copysign (octave_Inf, igral);
2120 warning ("quadcc: divergent integral detected");
2121 break;
2122 }
2123
2124 /* Compute the local integral. */
2125 ivr->igral = h * w * ivr->c[0];
2126
2127
2128 /* Fix-up the heap: we now have one interval on top
2129 that we don't need any more and two new, unsorted
2130 ones at the bottom. */
2131 /* Flip the last interval to the top of the heap and
2132 sift down. */
2133 t = heap[nivals - 1];
2134 heap[nivals - 1] = heap[0];
2135 heap[0] = t;
2136 nivals--;
2137 /* Sift this interval back down the heap. */
2138 i = 0;
2139 while (2 * i + 1 < nivals - 1)
2140 {
2141 j = 2 * i + 1;
2142 if (j + 1 < nivals - 1
2143 && ivals[heap[j + 1]].err >= ivals[heap[j]].err)
2144 j++;
2145 if (ivals[heap[j]].err <= ivals[heap[i]].err)
2146 break;
2147 else
2148 {
2149 t = heap[j];
2150 heap[j] = heap[i];
2151 heap[i] = t;
2152 i = j;
2153 }
2154 }
2155
2156 /* Now grab the last interval and sift it up the heap. */
2157 i = nivals - 1;
2158 while (i > 0)
2159 {
2160 j = (i - 1) / 2;
2161 if (ivals[heap[j]].err < ivals[heap[i]].err)
2162 {
2163 t = heap[j];
2164 heap[j] = heap[i];
2165 heap[i] = t;
2166 i = j;
2167 }
2168 else
2169 break;
2170 }
2171
2172
2173 }
2174
2175 /* Otherwise, just fix-up the heap. */
2176 else
2177 {
2178 i = 0;
2179 while (2 * i + 1 < nivals)
2180 {
2181 j = 2 * i + 1;
2182 if (j + 1 < nivals
2183 && ivals[heap[j + 1]].err >= ivals[heap[j]].err)
2184 j++;
2185 if (ivals[heap[j]].err <= ivals[heap[i]].err)
2186 break;
2187 else
2188 {
2189 t = heap[j];
2190 heap[j] = heap[i];
2191 heap[i] = t;
2192 i = j;
2193 }
2194 }
2195
2196 }
2197
2198 /* If the heap is about to overflow, remove the last two
2199 intervals. */
2200 while (nivals > cquad_heapsize - 2)
2201 {
2202 iv = &(ivals[heap[nivals - 1]]);
2203 /* printf
2204 ("quadcc: dropping ival %i (of %i) with [%e,%e] int=%e, err=%e, depth=%i\n",
2205 heap[0], nivals, iv->a, iv->b, iv->igral, iv->err,
2206 iv->depth);
2207 */
2208 err_final += iv->err;
2209 igral_final += iv->igral;
2210 nivals--;
2211 }
2212
2213 /* Collect the value of the integral and error. */
2214 igral = igral_final;
2215 err = err_final;
2216 for (i = 0; i < nivals; i++)
2217 {
2218 igral += ivals[heap[i]].igral;
2219 err += ivals[heap[i]].err;
2220 }
2221
2222 }
2223
2224 /* Dump the contents of the heap. */
2225 /* for (i = 0; i < nivals; i++)
2226 {
2227 iv = &(ivals[heap[i]]);
2228 printf
2229 ("quadcc: ival %i (%i) with [%e,%e], int=%e, err=%e, depth=%i, rdepth=%i, ndiv=%i\n",
2230 i, heap[i], iv->a, iv->b, iv->igral, iv->err, iv->depth,
2231 iv->rdepth, iv->ndiv);
2232 }
2233 */
2234 /* Clean up and present the results. */
2235 if (nargout > 2)
2236 retval(2) = neval;
2237 if (nargout > 1)
2238 retval(1) = err;
2239 retval(0) = igral;
2240 /* All is well that ends well. */
2241 return retval;
2242 }
2243
2244
2245 /*
2246 %!assert (quadcc (@sin, -pi, pi), 0, 1e-6)
2247 %!assert (quadcc (inline ("sin"),- pi, pi), 0, 1e-6)
2248 %!assert (quadcc ("sin", -pi, pi), 0, 1e-6)
2249
2250 %!assert (quadcc (@sin, -pi, 0), -2, 1e-6)
2251 %!assert (quadcc (@sin, 0, pi), 2, 1e-6)
2252 %!assert (quadcc (@(x) 1./sqrt (x), 0, 1), 2, 1e-6)
2253 %!assert (quadcc (@(x) 1./(sqrt (x).*(x+1)), 0, Inf), pi, 1e-6)
2254
2255 %!assert (quadcc (@(x) exp (-x .^ 2), -Inf, Inf), sqrt (pi), 1e-6)
2256 %!assert (quadcc (@(x) exp (-x .^ 2), -Inf, 0), sqrt (pi)/2, 1e-6)
2257
2258 %% Test input validation
2259 %!error (quadcc ())
2260 %!error (quadcc (@sin))
2261 %!error (quadcc (@sin, 0))
2262 %!error (quadcc (@sin, ones (2), pi))
2263 %!error (quadcc (@sin, -i, pi))
2264 %!error (quadcc (@sin, 0, ones (2)))
2265 %!error (quadcc (@sin, 0, i))
2266 %!error (quadcc (@sin, 0, pi, 0))
2267 %!error (quadcc (@sin, 0, pi, 1e-6, [ i ]))
2268 */