Free cookie consent management tool by TermsFeed Policy Generator

source: branches/3.2/sources/HeuristicLab.ExtLibs/HeuristicLab.ALGLIB/2.1.2.2591/ALGLIB-2.1.2.2591/gq.cs @ 4539

Last change on this file since 4539 was 2645, checked in by mkommend, 15 years ago

extracted external libraries and adapted dependent plugins (ticket #837)

File size: 24.4 KB
Line 
1/*************************************************************************
2Copyright (c) 2005-2007, Sergey Bochkanov (ALGLIB project).
3
4>>> SOURCE LICENSE >>>
5This program is free software; you can redistribute it and/or modify
6it under the terms of the GNU General Public License as published by
7the Free Software Foundation (www.fsf.org); either version 2 of the
8License, or (at your option) any later version.
9
10This program is distributed in the hope that it will be useful,
11but WITHOUT ANY WARRANTY; without even the implied warranty of
12MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
13GNU General Public License for more details.
14
15A copy of the GNU General Public License is available at
16http://www.fsf.org/licensing/licenses
17
18>>> END OF LICENSE >>>
19*************************************************************************/
20
21using System;
22
23namespace alglib
24{
25    public class gq
26    {
27        /*************************************************************************
28        Computation of nodes and weights for a Gauss quadrature formula
29
30        The algorithm generates the N-point Gauss quadrature formula  with  weight
31        function given by coefficients alpha and beta  of  a  recurrence  relation
32        which generates a system of orthogonal polynomials:
33
34        P-1(x)   =  0
35        P0(x)    =  1
36        Pn+1(x)  =  (x-alpha(n))*Pn(x)  -  beta(n)*Pn-1(x)
37
38        and zeroth moment Mu0
39
40        Mu0 = integral(W(x)dx,a,b)
41
42        INPUT PARAMETERS:
43            Alpha   –   array[0..N-1], alpha coefficients
44            Beta    –   array[0..N-1], beta coefficients
45                        Zero-indexed element is not used and may be arbitrary.
46                        Beta[I]>0.
47            Mu0     –   zeroth moment of the weight function.
48            N       –   number of nodes of the quadrature formula, N>=1
49
50        OUTPUT PARAMETERS:
51            Info    -   error code:
52                        * -3    internal eigenproblem solver hasn't converged
53                        * -2    Beta[i]<=0
54                        * -1    incorrect N was passed
55                        *  1    OK
56            X       -   array[0..N-1] - array of quadrature nodes,
57                        in ascending order.
58            W       -   array[0..N-1] - array of quadrature weights.
59
60          -- ALGLIB --
61             Copyright 2005-2009 by Bochkanov Sergey
62        *************************************************************************/
63        public static void gqgeneraterec(ref double[] alpha,
64            ref double[] beta,
65            double mu0,
66            int n,
67            ref int info,
68            ref double[] x,
69            ref double[] w)
70        {
71            int i = 0;
72            double[] d = new double[0];
73            double[] e = new double[0];
74            double[,] z = new double[0,0];
75
76            if( n<1 )
77            {
78                info = -1;
79                return;
80            }
81            info = 1;
82           
83            //
84            // Initialize
85            //
86            d = new double[n+1];
87            e = new double[n+1];
88            for(i=1; i<=n-1; i++)
89            {
90                d[i] = alpha[i-1];
91                if( (double)(beta[i])<=(double)(0) )
92                {
93                    info = -2;
94                    return;
95                }
96                e[i] = Math.Sqrt(beta[i]);
97            }
98            d[n] = alpha[n-1];
99           
100            //
101            // EVD
102            //
103            if( !tdevd.tridiagonalevd(ref d, e, n, 3, ref z) )
104            {
105                info = -3;
106                return;
107            }
108           
109            //
110            // Generate
111            //
112            x = new double[n];
113            w = new double[n];
114            for(i=1; i<=n; i++)
115            {
116                x[i-1] = d[i];
117                w[i-1] = mu0*AP.Math.Sqr(z[1,i]);
118            }
119        }
120
121
122        /*************************************************************************
123        Computation of nodes and weights for a Gauss-Lobatto quadrature formula
124
125        The algorithm generates the N-point Gauss-Lobatto quadrature formula  with
126        weight function given by coefficients alpha and beta of a recurrence which
127        generates a system of orthogonal polynomials.
128
129        P-1(x)   =  0
130        P0(x)    =  1
131        Pn+1(x)  =  (x-alpha(n))*Pn(x)  -  beta(n)*Pn-1(x)
132
133        and zeroth moment Mu0
134
135        Mu0 = integral(W(x)dx,a,b)
136
137        INPUT PARAMETERS:
138            Alpha   –   array[0..N-2], alpha coefficients
139            Beta    –   array[0..N-2], beta coefficients.
140                        Zero-indexed element is not used, may be arbitrary.
141                        Beta[I]>0
142            Mu0     –   zeroth moment of the weighting function.
143            A       –   left boundary of the integration interval.
144            B       –   right boundary of the integration interval.
145            N       –   number of nodes of the quadrature formula, N>=3
146                        (including the left and right boundary nodes).
147
148        OUTPUT PARAMETERS:
149            Info    -   error code:
150                        * -3    internal eigenproblem solver hasn't converged
151                        * -2    Beta[i]<=0
152                        * -1    incorrect N was passed
153                        *  1    OK
154            X       -   array[0..N-1] - array of quadrature nodes,
155                        in ascending order.
156            W       -   array[0..N-1] - array of quadrature weights.
157
158          -- ALGLIB --
159             Copyright 2005-2009 by Bochkanov Sergey
160        *************************************************************************/
161        public static void gqgenerategausslobattorec(double[] alpha,
162            double[] beta,
163            double mu0,
164            double a,
165            double b,
166            int n,
167            ref int info,
168            ref double[] x,
169            ref double[] w)
170        {
171            int i = 0;
172            double[] d = new double[0];
173            double[] e = new double[0];
174            double[,] z = new double[0,0];
175            double pim1a = 0;
176            double pia = 0;
177            double pim1b = 0;
178            double pib = 0;
179            double t = 0;
180            double a11 = 0;
181            double a12 = 0;
182            double a21 = 0;
183            double a22 = 0;
184            double b1 = 0;
185            double b2 = 0;
186            double alph = 0;
187            double bet = 0;
188
189            alpha = (double[])alpha.Clone();
190            beta = (double[])beta.Clone();
191
192            if( n<=2 )
193            {
194                info = -1;
195                return;
196            }
197            info = 1;
198           
199            //
200            // Initialize, D[1:N+1], E[1:N]
201            //
202            n = n-2;
203            d = new double[n+3];
204            e = new double[n+2];
205            for(i=1; i<=n+1; i++)
206            {
207                d[i] = alpha[i-1];
208            }
209            for(i=1; i<=n; i++)
210            {
211                if( (double)(beta[i])<=(double)(0) )
212                {
213                    info = -2;
214                    return;
215                }
216                e[i] = Math.Sqrt(beta[i]);
217            }
218           
219            //
220            // Caclulate Pn(a), Pn+1(a), Pn(b), Pn+1(b)
221            //
222            beta[0] = 0;
223            pim1a = 0;
224            pia = 1;
225            pim1b = 0;
226            pib = 1;
227            for(i=1; i<=n+1; i++)
228            {
229               
230                //
231                // Pi(a)
232                //
233                t = (a-alpha[i-1])*pia-beta[i-1]*pim1a;
234                pim1a = pia;
235                pia = t;
236               
237                //
238                // Pi(b)
239                //
240                t = (b-alpha[i-1])*pib-beta[i-1]*pim1b;
241                pim1b = pib;
242                pib = t;
243            }
244           
245            //
246            // Calculate alpha'(n+1), beta'(n+1)
247            //
248            a11 = pia;
249            a12 = pim1a;
250            a21 = pib;
251            a22 = pim1b;
252            b1 = a*pia;
253            b2 = b*pib;
254            if( (double)(Math.Abs(a11))>(double)(Math.Abs(a21)) )
255            {
256                a22 = a22-a12*a21/a11;
257                b2 = b2-b1*a21/a11;
258                bet = b2/a22;
259                alph = (b1-bet*a12)/a11;
260            }
261            else
262            {
263                a12 = a12-a22*a11/a21;
264                b1 = b1-b2*a11/a21;
265                bet = b1/a12;
266                alph = (b2-bet*a22)/a21;
267            }
268            if( (double)(bet)<(double)(0) )
269            {
270                info = -3;
271                return;
272            }
273            d[n+2] = alph;
274            e[n+1] = Math.Sqrt(bet);
275           
276            //
277            // EVD
278            //
279            if( !tdevd.tridiagonalevd(ref d, e, n+2, 3, ref z) )
280            {
281                info = -3;
282                return;
283            }
284           
285            //
286            // Generate
287            //
288            x = new double[n+2];
289            w = new double[n+2];
290            for(i=1; i<=n+2; i++)
291            {
292                x[i-1] = d[i];
293                w[i-1] = mu0*AP.Math.Sqr(z[1,i]);
294            }
295        }
296
297
298        /*************************************************************************
299        Computation of nodes and weights for a Gauss-Radau quadrature formula
300
301        The algorithm generates the N-point Gauss-Radau  quadrature  formula  with
302        weight function given by the coefficients alpha and  beta  of a recurrence
303        which generates a system of orthogonal polynomials.
304
305        P-1(x)   =  0
306        P0(x)    =  1
307        Pn+1(x)  =  (x-alpha(n))*Pn(x)  -  beta(n)*Pn-1(x)
308
309        and zeroth moment Mu0
310
311        Mu0 = integral(W(x)dx,a,b)
312
313        INPUT PARAMETERS:
314            Alpha   –   array[0..N-2], alpha coefficients.
315            Beta    –   array[0..N-1], beta coefficients
316                        Zero-indexed element is not used.
317                        Beta[I]>0
318            Mu0     –   zeroth moment of the weighting function.
319            A       –   left boundary of the integration interval.
320            N       –   number of nodes of the quadrature formula, N>=2
321                        (including the left boundary node).
322
323        OUTPUT PARAMETERS:
324            Info    -   error code:
325                        * -3    internal eigenproblem solver hasn't converged
326                        * -2    Beta[i]<=0
327                        * -1    incorrect N was passed
328                        *  1    OK
329            X       -   array[0..N-1] - array of quadrature nodes,
330                        in ascending order.
331            W       -   array[0..N-1] - array of quadrature weights.
332
333
334          -- ALGLIB --
335             Copyright 2005-2009 by Bochkanov Sergey
336        *************************************************************************/
337        public static void gqgenerategaussradaurec(double[] alpha,
338            double[] beta,
339            double mu0,
340            double a,
341            int n,
342            ref int info,
343            ref double[] x,
344            ref double[] w)
345        {
346            int i = 0;
347            double[] d = new double[0];
348            double[] e = new double[0];
349            double[,] z = new double[0,0];
350            double polim1 = 0;
351            double poli = 0;
352            double t = 0;
353
354            alpha = (double[])alpha.Clone();
355            beta = (double[])beta.Clone();
356
357            if( n<2 )
358            {
359                info = -1;
360                return;
361            }
362            info = 1;
363           
364            //
365            // Initialize, D[1:N], E[1:N]
366            //
367            n = n-1;
368            d = new double[n+1+1];
369            e = new double[n+1];
370            for(i=1; i<=n; i++)
371            {
372                d[i] = alpha[i-1];
373                if( (double)(beta[i])<=(double)(0) )
374                {
375                    info = -2;
376                    return;
377                }
378                e[i] = Math.Sqrt(beta[i]);
379            }
380           
381            //
382            // Caclulate Pn(a), Pn-1(a), and D[N+1]
383            //
384            beta[0] = 0;
385            polim1 = 0;
386            poli = 1;
387            for(i=1; i<=n; i++)
388            {
389                t = (a-alpha[i-1])*poli-beta[i-1]*polim1;
390                polim1 = poli;
391                poli = t;
392            }
393            d[n+1] = a-beta[n]*polim1/poli;
394           
395            //
396            // EVD
397            //
398            if( !tdevd.tridiagonalevd(ref d, e, n+1, 3, ref z) )
399            {
400                info = -3;
401                return;
402            }
403           
404            //
405            // Generate
406            //
407            x = new double[n+1];
408            w = new double[n+1];
409            for(i=1; i<=n+1; i++)
410            {
411                x[i-1] = d[i];
412                w[i-1] = mu0*AP.Math.Sqr(z[1,i]);
413            }
414        }
415
416
417        /*************************************************************************
418        Returns nodes/weights for Gauss-Legendre quadrature on [-1,1] with N
419        nodes.
420
421        INPUT PARAMETERS:
422            N           -   number of nodes, >=1
423
424        OUTPUT PARAMETERS:
425            Info        -   error code:
426                            * -4    an  error   was   detected   when  calculating
427                                    weights/nodes.  N  is  too  large   to  obtain
428                                    weights/nodes  with  high   enough   accuracy.
429                                    Try  to   use   multiple   precision  version.
430                            * -3    internal eigenproblem solver hasn't  converged
431                            * -1    incorrect N was passed
432                            * +1    OK
433            X           -   array[0..N-1] - array of quadrature nodes,
434                            in ascending order.
435            W           -   array[0..N-1] - array of quadrature weights.
436
437
438          -- ALGLIB --
439             Copyright 12.05.2009 by Bochkanov Sergey
440        *************************************************************************/
441        public static void gqgenerategausslegendre(int n,
442            ref int info,
443            ref double[] x,
444            ref double[] w)
445        {
446            double[] alpha = new double[0];
447            double[] beta = new double[0];
448            int i = 0;
449
450            if( n<1 )
451            {
452                info = -1;
453                return;
454            }
455            alpha = new double[n];
456            beta = new double[n];
457            for(i=0; i<=n-1; i++)
458            {
459                alpha[i] = 0;
460            }
461            beta[0] = 2;
462            for(i=1; i<=n-1; i++)
463            {
464                beta[i] = 1/(4-1/AP.Math.Sqr(i));
465            }
466            gqgeneraterec(ref alpha, ref beta, beta[0], n, ref info, ref x, ref w);
467           
468            //
469            // test basic properties to detect errors
470            //
471            if( info>0 )
472            {
473                if( (double)(x[0])<(double)(-1) | (double)(x[n-1])>(double)(+1) )
474                {
475                    info = -4;
476                }
477                for(i=0; i<=n-2; i++)
478                {
479                    if( (double)(x[i])>=(double)(x[i+1]) )
480                    {
481                        info = -4;
482                    }
483                }
484            }
485        }
486
487
488        /*************************************************************************
489        Returns  nodes/weights  for  Gauss-Jacobi quadrature on [-1,1] with weight
490        function W(x)=Power(1-x,Alpha)*Power(1+x,Beta).
491
492        INPUT PARAMETERS:
493            N           -   number of nodes, >=1
494            Alpha       -   power-law coefficient, Alpha>-1
495            Beta        -   power-law coefficient, Beta>-1
496
497        OUTPUT PARAMETERS:
498            Info        -   error code:
499                            * -4    an  error  was   detected   when   calculating
500                                    weights/nodes. Alpha or  Beta  are  too  close
501                                    to -1 to obtain weights/nodes with high enough
502                                    accuracy, or, may be, N is too large.  Try  to
503                                    use multiple precision version.
504                            * -3    internal eigenproblem solver hasn't converged
505                            * -1    incorrect N/Alpha/Beta was passed
506                            * +1    OK
507            X           -   array[0..N-1] - array of quadrature nodes,
508                            in ascending order.
509            W           -   array[0..N-1] - array of quadrature weights.
510
511
512          -- ALGLIB --
513             Copyright 12.05.2009 by Bochkanov Sergey
514        *************************************************************************/
515        public static void gqgenerategaussjacobi(int n,
516            double alpha,
517            double beta,
518            ref int info,
519            ref double[] x,
520            ref double[] w)
521        {
522            double[] a = new double[0];
523            double[] b = new double[0];
524            double alpha2 = 0;
525            double beta2 = 0;
526            double apb = 0;
527            double t = 0;
528            int i = 0;
529            double s = 0;
530
531            if( n<1 | (double)(alpha)<=(double)(-1) | (double)(beta)<=(double)(-1) )
532            {
533                info = -1;
534                return;
535            }
536            a = new double[n];
537            b = new double[n];
538            apb = alpha+beta;
539            a[0] = (beta-alpha)/(apb+2);
540            t = (apb+1)*Math.Log(2)+gammafunc.lngamma(alpha+1, ref s)+gammafunc.lngamma(beta+1, ref s)-gammafunc.lngamma(apb+2, ref s);
541            if( (double)(t)>(double)(Math.Log(AP.Math.MaxRealNumber)) )
542            {
543                info = -4;
544                return;
545            }
546            b[0] = Math.Exp(t);
547            if( n>1 )
548            {
549                alpha2 = AP.Math.Sqr(alpha);
550                beta2 = AP.Math.Sqr(beta);
551                a[1] = (beta2-alpha2)/((apb+2)*(apb+4));
552                b[1] = 4*(alpha+1)*(beta+1)/((apb+3)*AP.Math.Sqr(apb+2));
553                for(i=2; i<=n-1; i++)
554                {
555                    a[i] = 0.25*(beta2-alpha2)/(i*i*(1+0.5*apb/i)*(1+0.5*(apb+2)/i));
556                    b[i] = 0.25*(1+alpha/i)*(1+beta/i)*(1+apb/i)/((1+0.5*(apb+1)/i)*(1+0.5*(apb-1)/i)*AP.Math.Sqr(1+0.5*apb/i));
557                }
558            }
559            gqgeneraterec(ref a, ref b, b[0], n, ref info, ref x, ref w);
560           
561            //
562            // test basic properties to detect errors
563            //
564            if( info>0 )
565            {
566                if( (double)(x[0])<(double)(-1) | (double)(x[n-1])>(double)(+1) )
567                {
568                    info = -4;
569                }
570                for(i=0; i<=n-2; i++)
571                {
572                    if( (double)(x[i])>=(double)(x[i+1]) )
573                    {
574                        info = -4;
575                    }
576                }
577            }
578        }
579
580
581        /*************************************************************************
582        Returns  nodes/weights  for  Gauss-Laguerre  quadrature  on  [0,+inf) with
583        weight function W(x)=Power(x,Alpha)*Exp(-x)
584
585        INPUT PARAMETERS:
586            N           -   number of nodes, >=1
587            Alpha       -   power-law coefficient, Alpha>-1
588
589        OUTPUT PARAMETERS:
590            Info        -   error code:
591                            * -4    an  error  was   detected   when   calculating
592                                    weights/nodes. Alpha is too  close  to  -1  to
593                                    obtain weights/nodes with high enough accuracy
594                                    or, may  be,  N  is  too  large.  Try  to  use
595                                    multiple precision version.
596                            * -3    internal eigenproblem solver hasn't converged
597                            * -1    incorrect N/Alpha was passed
598                            * +1    OK
599            X           -   array[0..N-1] - array of quadrature nodes,
600                            in ascending order.
601            W           -   array[0..N-1] - array of quadrature weights.
602
603
604          -- ALGLIB --
605             Copyright 12.05.2009 by Bochkanov Sergey
606        *************************************************************************/
607        public static void gqgenerategausslaguerre(int n,
608            double alpha,
609            ref int info,
610            ref double[] x,
611            ref double[] w)
612        {
613            double[] a = new double[0];
614            double[] b = new double[0];
615            double t = 0;
616            int i = 0;
617            double s = 0;
618
619            if( n<1 | (double)(alpha)<=(double)(-1) )
620            {
621                info = -1;
622                return;
623            }
624            a = new double[n];
625            b = new double[n];
626            a[0] = alpha+1;
627            t = gammafunc.lngamma(alpha+1, ref s);
628            if( (double)(t)>=(double)(Math.Log(AP.Math.MaxRealNumber)) )
629            {
630                info = -4;
631                return;
632            }
633            b[0] = Math.Exp(t);
634            if( n>1 )
635            {
636                for(i=1; i<=n-1; i++)
637                {
638                    a[i] = 2*i+alpha+1;
639                    b[i] = i*(i+alpha);
640                }
641            }
642            gqgeneraterec(ref a, ref b, b[0], n, ref info, ref x, ref w);
643           
644            //
645            // test basic properties to detect errors
646            //
647            if( info>0 )
648            {
649                if( (double)(x[0])<(double)(0) )
650                {
651                    info = -4;
652                }
653                for(i=0; i<=n-2; i++)
654                {
655                    if( (double)(x[i])>=(double)(x[i+1]) )
656                    {
657                        info = -4;
658                    }
659                }
660            }
661        }
662
663
664        /*************************************************************************
665        Returns  nodes/weights  for  Gauss-Hermite  quadrature on (-inf,+inf) with
666        weight function W(x)=Exp(-x*x)
667
668        INPUT PARAMETERS:
669            N           -   number of nodes, >=1
670
671        OUTPUT PARAMETERS:
672            Info        -   error code:
673                            * -4    an  error  was   detected   when   calculating
674                                    weights/nodes.  May be, N is too large. Try to
675                                    use multiple precision version.
676                            * -3    internal eigenproblem solver hasn't converged
677                            * -1    incorrect N/Alpha was passed
678                            * +1    OK
679            X           -   array[0..N-1] - array of quadrature nodes,
680                            in ascending order.
681            W           -   array[0..N-1] - array of quadrature weights.
682
683
684          -- ALGLIB --
685             Copyright 12.05.2009 by Bochkanov Sergey
686        *************************************************************************/
687        public static void gqgenerategausshermite(int n,
688            ref int info,
689            ref double[] x,
690            ref double[] w)
691        {
692            double[] a = new double[0];
693            double[] b = new double[0];
694            double t = 0;
695            int i = 0;
696            double s = 0;
697
698            if( n<1 )
699            {
700                info = -1;
701                return;
702            }
703            a = new double[n];
704            b = new double[n];
705            for(i=0; i<=n-1; i++)
706            {
707                a[i] = 0;
708            }
709            b[0] = Math.Sqrt(4*Math.Atan(1));
710            if( n>1 )
711            {
712                for(i=1; i<=n-1; i++)
713                {
714                    b[i] = 0.5*i;
715                }
716            }
717            gqgeneraterec(ref a, ref b, b[0], n, ref info, ref x, ref w);
718           
719            //
720            // test basic properties to detect errors
721            //
722            if( info>0 )
723            {
724                for(i=0; i<=n-2; i++)
725                {
726                    if( (double)(x[i])>=(double)(x[i+1]) )
727                    {
728                        info = -4;
729                    }
730                }
731            }
732        }
733    }
734}
Note: See TracBrowser for help on using the repository browser.