Class documentation of Concepts

Loading...
Searching...
No Matches
integration.h
1/*************************************************************************
2ALGLIB 3.11.0 (source code generated 2017-05-11)
3Copyright (c) Sergey Bochkanov (ALGLIB project).
4
5>>> SOURCE LICENSE >>>
6This program is free software; you can redistribute it and/or modify
7it under the terms of the GNU General Public License as published by
8the Free Software Foundation (www.fsf.org); either version 2 of the
9License, or (at your option) any later version.
10
11This program is distributed in the hope that it will be useful,
12but WITHOUT ANY WARRANTY; without even the implied warranty of
13MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14GNU General Public License for more details.
15
16A copy of the GNU General Public License is available at
17http://www.fsf.org/licensing/licenses
18>>> END OF LICENSE >>>
19*************************************************************************/
20#ifndef _integration_pkg_h
21#define _integration_pkg_h
22#include "ap.h"
23#include "alglibinternal.h"
24#include "linalg.h"
25#include "alglibmisc.h"
26#include "specialfunctions.h"
27
29//
30// THIS SECTION CONTAINS COMPUTATIONAL CORE DECLARATIONS (DATATYPES)
31//
33namespace alglib_impl
34{
35typedef struct
36{
37 ae_int_t terminationtype;
38 ae_int_t nfev;
39 ae_int_t nintervals;
41typedef struct
42{
43 double a;
44 double b;
45 double eps;
46 double xwidth;
47 double x;
48 double f;
49 ae_int_t info;
50 double r;
51 ae_matrix heap;
52 ae_int_t heapsize;
53 ae_int_t heapwidth;
54 ae_int_t heapused;
55 double sumerr;
56 double sumabs;
57 ae_vector qn;
58 ae_vector wg;
59 ae_vector wk;
60 ae_vector wr;
61 ae_int_t n;
62 rcommstate rstate;
64typedef struct
65{
66 double a;
67 double b;
68 double alpha;
69 double beta;
70 double xwidth;
71 double x;
72 double xminusa;
73 double bminusx;
74 ae_bool needf;
75 double f;
76 ae_int_t wrappermode;
77 autogkinternalstate internalstate;
78 rcommstate rstate;
79 double v;
80 ae_int_t terminationtype;
81 ae_int_t nfev;
82 ae_int_t nintervals;
84
85}
86
88//
89// THIS SECTION CONTAINS C++ INTERFACE
90//
92namespace alglib
93{
94
95
96
97
98
99/*************************************************************************
100Integration report:
101* TerminationType = completetion code:
102 * -5 non-convergence of Gauss-Kronrod nodes
103 calculation subroutine.
104 * -1 incorrect parameters were specified
105 * 1 OK
106* Rep.NFEV countains number of function calculations
107* Rep.NIntervals contains number of intervals [a,b]
108 was partitioned into.
109*************************************************************************/
111{
112public:
115 _autogkreport_owner& operator=(const _autogkreport_owner &rhs);
116 virtual ~_autogkreport_owner();
118 alglib_impl::autogkreport* c_ptr() const;
119protected:
121};
123{
124public:
125 autogkreport();
126 autogkreport(const autogkreport &rhs);
127 autogkreport& operator=(const autogkreport &rhs);
128 virtual ~autogkreport();
129 ae_int_t &terminationtype;
130 ae_int_t &nfev;
131 ae_int_t &nintervals;
132
133};
134
135
136/*************************************************************************
137This structure stores state of the integration algorithm.
138
139Although this class has public fields, they are not intended for external
140use. You should use ALGLIB functions to work with this class:
141* autogksmooth()/AutoGKSmoothW()/... to create objects
142* autogkintegrate() to begin integration
143* autogkresults() to get results
144*************************************************************************/
146{
147public:
150 _autogkstate_owner& operator=(const _autogkstate_owner &rhs);
151 virtual ~_autogkstate_owner();
153 alglib_impl::autogkstate* c_ptr() const;
154protected:
155 alglib_impl::autogkstate *p_struct;
156};
158{
159public:
160 autogkstate();
161 autogkstate(const autogkstate &rhs);
162 autogkstate& operator=(const autogkstate &rhs);
163 virtual ~autogkstate();
164 ae_bool &needf;
165 double &x;
166 double &xminusa;
167 double &bminusx;
168 double &f;
169
170};
171
172/*************************************************************************
173Computation of nodes and weights for a Gauss quadrature formula
174
175The algorithm generates the N-point Gauss quadrature formula with weight
176function given by coefficients alpha and beta of a recurrence relation
177which generates a system of orthogonal polynomials:
178
179P-1(x) = 0
180P0(x) = 1
181Pn+1(x) = (x-alpha(n))*Pn(x) - beta(n)*Pn-1(x)
182
183and zeroth moment Mu0
184
185Mu0 = integral(W(x)dx,a,b)
186
187INPUT PARAMETERS:
188 Alpha - array[0..N-1], alpha coefficients
189 Beta - array[0..N-1], beta coefficients
190 Zero-indexed element is not used and may be arbitrary.
191 Beta[I]>0.
192 Mu0 - zeroth moment of the weight function.
193 N - number of nodes of the quadrature formula, N>=1
194
195OUTPUT PARAMETERS:
196 Info - error code:
197 * -3 internal eigenproblem solver hasn't converged
198 * -2 Beta[i]<=0
199 * -1 incorrect N was passed
200 * 1 OK
201 X - array[0..N-1] - array of quadrature nodes,
202 in ascending order.
203 W - array[0..N-1] - array of quadrature weights.
204
205 -- ALGLIB --
206 Copyright 2005-2009 by Bochkanov Sergey
207*************************************************************************/
208void gqgeneraterec(const real_1d_array &alpha, const real_1d_array &beta, const double mu0, const ae_int_t n, ae_int_t &info, real_1d_array &x, real_1d_array &w);
209
210
211/*************************************************************************
212Computation of nodes and weights for a Gauss-Lobatto quadrature formula
213
214The algorithm generates the N-point Gauss-Lobatto quadrature formula with
215weight function given by coefficients alpha and beta of a recurrence which
216generates a system of orthogonal polynomials.
217
218P-1(x) = 0
219P0(x) = 1
220Pn+1(x) = (x-alpha(n))*Pn(x) - beta(n)*Pn-1(x)
221
222and zeroth moment Mu0
223
224Mu0 = integral(W(x)dx,a,b)
225
226INPUT PARAMETERS:
227 Alpha - array[0..N-2], alpha coefficients
228 Beta - array[0..N-2], beta coefficients.
229 Zero-indexed element is not used, may be arbitrary.
230 Beta[I]>0
231 Mu0 - zeroth moment of the weighting function.
232 A - left boundary of the integration interval.
233 B - right boundary of the integration interval.
234 N - number of nodes of the quadrature formula, N>=3
235 (including the left and right boundary nodes).
236
237OUTPUT PARAMETERS:
238 Info - error code:
239 * -3 internal eigenproblem solver hasn't converged
240 * -2 Beta[i]<=0
241 * -1 incorrect N was passed
242 * 1 OK
243 X - array[0..N-1] - array of quadrature nodes,
244 in ascending order.
245 W - array[0..N-1] - array of quadrature weights.
246
247 -- ALGLIB --
248 Copyright 2005-2009 by Bochkanov Sergey
249*************************************************************************/
250void gqgenerategausslobattorec(const real_1d_array &alpha, const real_1d_array &beta, const double mu0, const double a, const double b, const ae_int_t n, ae_int_t &info, real_1d_array &x, real_1d_array &w);
251
252
253/*************************************************************************
254Computation of nodes and weights for a Gauss-Radau quadrature formula
255
256The algorithm generates the N-point Gauss-Radau quadrature formula with
257weight function given by the coefficients alpha and beta of a recurrence
258which generates a system of orthogonal polynomials.
259
260P-1(x) = 0
261P0(x) = 1
262Pn+1(x) = (x-alpha(n))*Pn(x) - beta(n)*Pn-1(x)
263
264and zeroth moment Mu0
265
266Mu0 = integral(W(x)dx,a,b)
267
268INPUT PARAMETERS:
269 Alpha - array[0..N-2], alpha coefficients.
270 Beta - array[0..N-1], beta coefficients
271 Zero-indexed element is not used.
272 Beta[I]>0
273 Mu0 - zeroth moment of the weighting function.
274 A - left boundary of the integration interval.
275 N - number of nodes of the quadrature formula, N>=2
276 (including the left boundary node).
277
278OUTPUT PARAMETERS:
279 Info - error code:
280 * -3 internal eigenproblem solver hasn't converged
281 * -2 Beta[i]<=0
282 * -1 incorrect N was passed
283 * 1 OK
284 X - array[0..N-1] - array of quadrature nodes,
285 in ascending order.
286 W - array[0..N-1] - array of quadrature weights.
287
288
289 -- ALGLIB --
290 Copyright 2005-2009 by Bochkanov Sergey
291*************************************************************************/
292void gqgenerategaussradaurec(const real_1d_array &alpha, const real_1d_array &beta, const double mu0, const double a, const ae_int_t n, ae_int_t &info, real_1d_array &x, real_1d_array &w);
293
294
295/*************************************************************************
296Returns nodes/weights for Gauss-Legendre quadrature on [-1,1] with N
297nodes.
298
299INPUT PARAMETERS:
300 N - number of nodes, >=1
301
302OUTPUT PARAMETERS:
303 Info - error code:
304 * -4 an error was detected when calculating
305 weights/nodes. N is too large to obtain
306 weights/nodes with high enough accuracy.
307 Try to use multiple precision version.
308 * -3 internal eigenproblem solver hasn't converged
309 * -1 incorrect N was passed
310 * +1 OK
311 X - array[0..N-1] - array of quadrature nodes,
312 in ascending order.
313 W - array[0..N-1] - array of quadrature weights.
314
315
316 -- ALGLIB --
317 Copyright 12.05.2009 by Bochkanov Sergey
318*************************************************************************/
319void gqgenerategausslegendre(const ae_int_t n, ae_int_t &info, real_1d_array &x, real_1d_array &w);
320
321
322/*************************************************************************
323Returns nodes/weights for Gauss-Jacobi quadrature on [-1,1] with weight
324function W(x)=Power(1-x,Alpha)*Power(1+x,Beta).
325
326INPUT PARAMETERS:
327 N - number of nodes, >=1
328 Alpha - power-law coefficient, Alpha>-1
329 Beta - power-law coefficient, Beta>-1
330
331OUTPUT PARAMETERS:
332 Info - error code:
333 * -4 an error was detected when calculating
334 weights/nodes. Alpha or Beta are too close
335 to -1 to obtain weights/nodes with high enough
336 accuracy, or, may be, N is too large. Try to
337 use multiple precision version.
338 * -3 internal eigenproblem solver hasn't converged
339 * -1 incorrect N/Alpha/Beta was passed
340 * +1 OK
341 X - array[0..N-1] - array of quadrature nodes,
342 in ascending order.
343 W - array[0..N-1] - array of quadrature weights.
344
345
346 -- ALGLIB --
347 Copyright 12.05.2009 by Bochkanov Sergey
348*************************************************************************/
349void gqgenerategaussjacobi(const ae_int_t n, const double alpha, const double beta, ae_int_t &info, real_1d_array &x, real_1d_array &w);
350
351
352/*************************************************************************
353Returns nodes/weights for Gauss-Laguerre quadrature on [0,+inf) with
354weight function W(x)=Power(x,Alpha)*Exp(-x)
355
356INPUT PARAMETERS:
357 N - number of nodes, >=1
358 Alpha - power-law coefficient, Alpha>-1
359
360OUTPUT PARAMETERS:
361 Info - error code:
362 * -4 an error was detected when calculating
363 weights/nodes. Alpha is too close to -1 to
364 obtain weights/nodes with high enough accuracy
365 or, may be, N is too large. Try to use
366 multiple precision version.
367 * -3 internal eigenproblem solver hasn't converged
368 * -1 incorrect N/Alpha was passed
369 * +1 OK
370 X - array[0..N-1] - array of quadrature nodes,
371 in ascending order.
372 W - array[0..N-1] - array of quadrature weights.
373
374
375 -- ALGLIB --
376 Copyright 12.05.2009 by Bochkanov Sergey
377*************************************************************************/
378void gqgenerategausslaguerre(const ae_int_t n, const double alpha, ae_int_t &info, real_1d_array &x, real_1d_array &w);
379
380
381/*************************************************************************
382Returns nodes/weights for Gauss-Hermite quadrature on (-inf,+inf) with
383weight function W(x)=Exp(-x*x)
384
385INPUT PARAMETERS:
386 N - number of nodes, >=1
387
388OUTPUT PARAMETERS:
389 Info - error code:
390 * -4 an error was detected when calculating
391 weights/nodes. May be, N is too large. Try to
392 use multiple precision version.
393 * -3 internal eigenproblem solver hasn't converged
394 * -1 incorrect N/Alpha was passed
395 * +1 OK
396 X - array[0..N-1] - array of quadrature nodes,
397 in ascending order.
398 W - array[0..N-1] - array of quadrature weights.
399
400
401 -- ALGLIB --
402 Copyright 12.05.2009 by Bochkanov Sergey
403*************************************************************************/
404void gqgenerategausshermite(const ae_int_t n, ae_int_t &info, real_1d_array &x, real_1d_array &w);
405
406/*************************************************************************
407Computation of nodes and weights of a Gauss-Kronrod quadrature formula
408
409The algorithm generates the N-point Gauss-Kronrod quadrature formula with
410weight function given by coefficients alpha and beta of a recurrence
411relation which generates a system of orthogonal polynomials:
412
413 P-1(x) = 0
414 P0(x) = 1
415 Pn+1(x) = (x-alpha(n))*Pn(x) - beta(n)*Pn-1(x)
416
417and zero moment Mu0
418
419 Mu0 = integral(W(x)dx,a,b)
420
421
422INPUT PARAMETERS:
423 Alpha - alpha coefficients, array[0..floor(3*K/2)].
424 Beta - beta coefficients, array[0..ceil(3*K/2)].
425 Beta[0] is not used and may be arbitrary.
426 Beta[I]>0.
427 Mu0 - zeroth moment of the weight function.
428 N - number of nodes of the Gauss-Kronrod quadrature formula,
429 N >= 3,
430 N = 2*K+1.
431
432OUTPUT PARAMETERS:
433 Info - error code:
434 * -5 no real and positive Gauss-Kronrod formula can
435 be created for such a weight function with a
436 given number of nodes.
437 * -4 N is too large, task may be ill conditioned -
438 x[i]=x[i+1] found.
439 * -3 internal eigenproblem solver hasn't converged
440 * -2 Beta[i]<=0
441 * -1 incorrect N was passed
442 * +1 OK
443 X - array[0..N-1] - array of quadrature nodes,
444 in ascending order.
445 WKronrod - array[0..N-1] - Kronrod weights
446 WGauss - array[0..N-1] - Gauss weights (interleaved with zeros
447 corresponding to extended Kronrod nodes).
448
449 -- ALGLIB --
450 Copyright 08.05.2009 by Bochkanov Sergey
451*************************************************************************/
452void gkqgeneraterec(const real_1d_array &alpha, const real_1d_array &beta, const double mu0, const ae_int_t n, ae_int_t &info, real_1d_array &x, real_1d_array &wkronrod, real_1d_array &wgauss);
453
454
455/*************************************************************************
456Returns Gauss and Gauss-Kronrod nodes/weights for Gauss-Legendre
457quadrature with N points.
458
459GKQLegendreCalc (calculation) or GKQLegendreTbl (precomputed table) is
460used depending on machine precision and number of nodes.
461
462INPUT PARAMETERS:
463 N - number of Kronrod nodes, must be odd number, >=3.
464
465OUTPUT PARAMETERS:
466 Info - error code:
467 * -4 an error was detected when calculating
468 weights/nodes. N is too large to obtain
469 weights/nodes with high enough accuracy.
470 Try to use multiple precision version.
471 * -3 internal eigenproblem solver hasn't converged
472 * -1 incorrect N was passed
473 * +1 OK
474 X - array[0..N-1] - array of quadrature nodes, ordered in
475 ascending order.
476 WKronrod - array[0..N-1] - Kronrod weights
477 WGauss - array[0..N-1] - Gauss weights (interleaved with zeros
478 corresponding to extended Kronrod nodes).
479
480
481 -- ALGLIB --
482 Copyright 12.05.2009 by Bochkanov Sergey
483*************************************************************************/
484void gkqgenerategausslegendre(const ae_int_t n, ae_int_t &info, real_1d_array &x, real_1d_array &wkronrod, real_1d_array &wgauss);
485
486
487/*************************************************************************
488Returns Gauss and Gauss-Kronrod nodes/weights for Gauss-Jacobi
489quadrature on [-1,1] with weight function
490
491 W(x)=Power(1-x,Alpha)*Power(1+x,Beta).
492
493INPUT PARAMETERS:
494 N - number of Kronrod nodes, must be odd number, >=3.
495 Alpha - power-law coefficient, Alpha>-1
496 Beta - power-law coefficient, Beta>-1
497
498OUTPUT PARAMETERS:
499 Info - error code:
500 * -5 no real and positive Gauss-Kronrod formula can
501 be created for such a weight function with a
502 given number of nodes.
503 * -4 an error was detected when calculating
504 weights/nodes. Alpha or Beta are too close
505 to -1 to obtain weights/nodes with high enough
506 accuracy, or, may be, N is too large. Try to
507 use multiple precision version.
508 * -3 internal eigenproblem solver hasn't converged
509 * -1 incorrect N was passed
510 * +1 OK
511 * +2 OK, but quadrature rule have exterior nodes,
512 x[0]<-1 or x[n-1]>+1
513 X - array[0..N-1] - array of quadrature nodes, ordered in
514 ascending order.
515 WKronrod - array[0..N-1] - Kronrod weights
516 WGauss - array[0..N-1] - Gauss weights (interleaved with zeros
517 corresponding to extended Kronrod nodes).
518
519
520 -- ALGLIB --
521 Copyright 12.05.2009 by Bochkanov Sergey
522*************************************************************************/
523void gkqgenerategaussjacobi(const ae_int_t n, const double alpha, const double beta, ae_int_t &info, real_1d_array &x, real_1d_array &wkronrod, real_1d_array &wgauss);
524
525
526/*************************************************************************
527Returns Gauss and Gauss-Kronrod nodes for quadrature with N points.
528
529Reduction to tridiagonal eigenproblem is used.
530
531INPUT PARAMETERS:
532 N - number of Kronrod nodes, must be odd number, >=3.
533
534OUTPUT PARAMETERS:
535 Info - error code:
536 * -4 an error was detected when calculating
537 weights/nodes. N is too large to obtain
538 weights/nodes with high enough accuracy.
539 Try to use multiple precision version.
540 * -3 internal eigenproblem solver hasn't converged
541 * -1 incorrect N was passed
542 * +1 OK
543 X - array[0..N-1] - array of quadrature nodes, ordered in
544 ascending order.
545 WKronrod - array[0..N-1] - Kronrod weights
546 WGauss - array[0..N-1] - Gauss weights (interleaved with zeros
547 corresponding to extended Kronrod nodes).
548
549 -- ALGLIB --
550 Copyright 12.05.2009 by Bochkanov Sergey
551*************************************************************************/
552void gkqlegendrecalc(const ae_int_t n, ae_int_t &info, real_1d_array &x, real_1d_array &wkronrod, real_1d_array &wgauss);
553
554
555/*************************************************************************
556Returns Gauss and Gauss-Kronrod nodes for quadrature with N points using
557pre-calculated table. Nodes/weights were computed with accuracy up to
5581.0E-32 (if MPFR version of ALGLIB is used). In standard double precision
559accuracy reduces to something about 2.0E-16 (depending on your compiler's
560handling of long floating point constants).
561
562INPUT PARAMETERS:
563 N - number of Kronrod nodes.
564 N can be 15, 21, 31, 41, 51, 61.
565
566OUTPUT PARAMETERS:
567 X - array[0..N-1] - array of quadrature nodes, ordered in
568 ascending order.
569 WKronrod - array[0..N-1] - Kronrod weights
570 WGauss - array[0..N-1] - Gauss weights (interleaved with zeros
571 corresponding to extended Kronrod nodes).
572
573
574 -- ALGLIB --
575 Copyright 12.05.2009 by Bochkanov Sergey
576*************************************************************************/
577void gkqlegendretbl(const ae_int_t n, real_1d_array &x, real_1d_array &wkronrod, real_1d_array &wgauss, double &eps);
578
579/*************************************************************************
580Integration of a smooth function F(x) on a finite interval [a,b].
581
582Fast-convergent algorithm based on a Gauss-Kronrod formula is used. Result
583is calculated with accuracy close to the machine precision.
584
585Algorithm works well only with smooth integrands. It may be used with
586continuous non-smooth integrands, but with less performance.
587
588It should never be used with integrands which have integrable singularities
589at lower or upper limits - algorithm may crash. Use AutoGKSingular in such
590cases.
591
592INPUT PARAMETERS:
593 A, B - interval boundaries (A<B, A=B or A>B)
594
595OUTPUT PARAMETERS
596 State - structure which stores algorithm state
597
598SEE ALSO
599 AutoGKSmoothW, AutoGKSingular, AutoGKResults.
600
601
602 -- ALGLIB --
603 Copyright 06.05.2009 by Bochkanov Sergey
604*************************************************************************/
605void autogksmooth(const double a, const double b, autogkstate &state);
606
607
608/*************************************************************************
609Integration of a smooth function F(x) on a finite interval [a,b].
610
611This subroutine is same as AutoGKSmooth(), but it guarantees that interval
612[a,b] is partitioned into subintervals which have width at most XWidth.
613
614Subroutine can be used when integrating nearly-constant function with
615narrow "bumps" (about XWidth wide). If "bumps" are too narrow, AutoGKSmooth
616subroutine can overlook them.
617
618INPUT PARAMETERS:
619 A, B - interval boundaries (A<B, A=B or A>B)
620
621OUTPUT PARAMETERS
622 State - structure which stores algorithm state
623
624SEE ALSO
625 AutoGKSmooth, AutoGKSingular, AutoGKResults.
626
627
628 -- ALGLIB --
629 Copyright 06.05.2009 by Bochkanov Sergey
630*************************************************************************/
631void autogksmoothw(const double a, const double b, const double xwidth, autogkstate &state);
632
633
634/*************************************************************************
635Integration on a finite interval [A,B].
636Integrand have integrable singularities at A/B.
637
638F(X) must diverge as "(x-A)^alpha" at A, as "(B-x)^beta" at B, with known
639alpha/beta (alpha>-1, beta>-1). If alpha/beta are not known, estimates
640from below can be used (but these estimates should be greater than -1 too).
641
642One of alpha/beta variables (or even both alpha/beta) may be equal to 0,
643which means than function F(x) is non-singular at A/B. Anyway (singular at
644bounds or not), function F(x) is supposed to be continuous on (A,B).
645
646Fast-convergent algorithm based on a Gauss-Kronrod formula is used. Result
647is calculated with accuracy close to the machine precision.
648
649INPUT PARAMETERS:
650 A, B - interval boundaries (A<B, A=B or A>B)
651 Alpha - power-law coefficient of the F(x) at A,
652 Alpha>-1
653 Beta - power-law coefficient of the F(x) at B,
654 Beta>-1
655
656OUTPUT PARAMETERS
657 State - structure which stores algorithm state
658
659SEE ALSO
660 AutoGKSmooth, AutoGKSmoothW, AutoGKResults.
661
662
663 -- ALGLIB --
664 Copyright 06.05.2009 by Bochkanov Sergey
665*************************************************************************/
666void autogksingular(const double a, const double b, const double alpha, const double beta, autogkstate &state);
667
668
669/*************************************************************************
670This function provides reverse communication interface
671Reverse communication interface is not documented or recommended to use.
672See below for functions which provide better documented API
673*************************************************************************/
674bool autogkiteration(const autogkstate &state);
675
676
677/*************************************************************************
678This function is used to launcn iterations of the 1-dimensional integrator
679
680It accepts following parameters:
681 func - callback which calculates f(x) for given x
682 ptr - optional pointer which is passed to func; can be NULL
683
684
685 -- ALGLIB --
686 Copyright 07.05.2009 by Bochkanov Sergey
687
688*************************************************************************/
689void autogkintegrate(autogkstate &state,
690 void (*func)(double x, double xminusa, double bminusx, double &y, void *ptr),
691 void *ptr = NULL);
692
693
694/*************************************************************************
695Adaptive integration results
696
697Called after AutoGKIteration returned False.
698
699Input parameters:
700 State - algorithm state (used by AutoGKIteration).
701
702Output parameters:
703 V - integral(f(x)dx,a,b)
704 Rep - optimization report (see AutoGKReport description)
705
706 -- ALGLIB --
707 Copyright 14.11.2007 by Bochkanov Sergey
708*************************************************************************/
709void autogkresults(const autogkstate &state, double &v, autogkreport &rep);
710}
711
713//
714// THIS SECTION CONTAINS COMPUTATIONAL CORE DECLARATIONS (FUNCTIONS)
715//
717namespace alglib_impl
718{
719void gqgeneraterec(/* Real */ ae_vector* alpha,
720 /* Real */ ae_vector* beta,
721 double mu0,
722 ae_int_t n,
723 ae_int_t* info,
724 /* Real */ ae_vector* x,
725 /* Real */ ae_vector* w,
726 ae_state *_state);
727void gqgenerategausslobattorec(/* Real */ ae_vector* alpha,
728 /* Real */ ae_vector* beta,
729 double mu0,
730 double a,
731 double b,
732 ae_int_t n,
733 ae_int_t* info,
734 /* Real */ ae_vector* x,
735 /* Real */ ae_vector* w,
736 ae_state *_state);
737void gqgenerategaussradaurec(/* Real */ ae_vector* alpha,
738 /* Real */ ae_vector* beta,
739 double mu0,
740 double a,
741 ae_int_t n,
742 ae_int_t* info,
743 /* Real */ ae_vector* x,
744 /* Real */ ae_vector* w,
745 ae_state *_state);
746void gqgenerategausslegendre(ae_int_t n,
747 ae_int_t* info,
748 /* Real */ ae_vector* x,
749 /* Real */ ae_vector* w,
750 ae_state *_state);
751void gqgenerategaussjacobi(ae_int_t n,
752 double alpha,
753 double beta,
754 ae_int_t* info,
755 /* Real */ ae_vector* x,
756 /* Real */ ae_vector* w,
757 ae_state *_state);
758void gqgenerategausslaguerre(ae_int_t n,
759 double alpha,
760 ae_int_t* info,
761 /* Real */ ae_vector* x,
762 /* Real */ ae_vector* w,
763 ae_state *_state);
764void gqgenerategausshermite(ae_int_t n,
765 ae_int_t* info,
766 /* Real */ ae_vector* x,
767 /* Real */ ae_vector* w,
768 ae_state *_state);
769void gkqgeneraterec(/* Real */ ae_vector* alpha,
770 /* Real */ ae_vector* beta,
771 double mu0,
772 ae_int_t n,
773 ae_int_t* info,
774 /* Real */ ae_vector* x,
775 /* Real */ ae_vector* wkronrod,
776 /* Real */ ae_vector* wgauss,
777 ae_state *_state);
778void gkqgenerategausslegendre(ae_int_t n,
779 ae_int_t* info,
780 /* Real */ ae_vector* x,
781 /* Real */ ae_vector* wkronrod,
782 /* Real */ ae_vector* wgauss,
783 ae_state *_state);
784void gkqgenerategaussjacobi(ae_int_t n,
785 double alpha,
786 double beta,
787 ae_int_t* info,
788 /* Real */ ae_vector* x,
789 /* Real */ ae_vector* wkronrod,
790 /* Real */ ae_vector* wgauss,
791 ae_state *_state);
792void gkqlegendrecalc(ae_int_t n,
793 ae_int_t* info,
794 /* Real */ ae_vector* x,
795 /* Real */ ae_vector* wkronrod,
796 /* Real */ ae_vector* wgauss,
797 ae_state *_state);
798void gkqlegendretbl(ae_int_t n,
799 /* Real */ ae_vector* x,
800 /* Real */ ae_vector* wkronrod,
801 /* Real */ ae_vector* wgauss,
802 double* eps,
803 ae_state *_state);
804void autogksmooth(double a,
805 double b,
806 autogkstate* state,
807 ae_state *_state);
808void autogksmoothw(double a,
809 double b,
810 double xwidth,
811 autogkstate* state,
812 ae_state *_state);
813void autogksingular(double a,
814 double b,
815 double alpha,
816 double beta,
817 autogkstate* state,
818 ae_state *_state);
819ae_bool autogkiteration(autogkstate* state, ae_state *_state);
820void autogkresults(autogkstate* state,
821 double* v,
822 autogkreport* rep,
823 ae_state *_state);
824void _autogkreport_init(void* _p, ae_state *_state);
825void _autogkreport_init_copy(void* _dst, void* _src, ae_state *_state);
826void _autogkreport_clear(void* _p);
827void _autogkreport_destroy(void* _p);
828void _autogkinternalstate_init(void* _p, ae_state *_state);
829void _autogkinternalstate_init_copy(void* _dst, void* _src, ae_state *_state);
830void _autogkinternalstate_clear(void* _p);
831void _autogkinternalstate_destroy(void* _p);
832void _autogkstate_init(void* _p, ae_state *_state);
833void _autogkstate_init_copy(void* _dst, void* _src, ae_state *_state);
834void _autogkstate_clear(void* _p);
835void _autogkstate_destroy(void* _p);
836
837}
838#endif
839