Free cookie consent management tool by TermsFeed Policy Generator

source: trunk/sources/ALGLIB/polint.cs @ 2575

Last change on this file since 2575 was 2563, checked in by gkronber, 15 years ago

Updated ALGLIB to latest version. #751 (Plugin for for data-modeling with ANN (integrated into CEDMA))

File size: 40.8 KB
Line 
1/*************************************************************************
2Copyright (c) 2006-2009, 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 polint
26    {
27        /*************************************************************************
28        Polynomial fitting report:
29            TaskRCond       reciprocal of task's condition number
30            RMSError        RMS error
31            AvgError        average error
32            AvgRelError     average relative error (for non-zero Y[I])
33            MaxError        maximum error
34        *************************************************************************/
35        public struct polynomialfitreport
36        {
37            public double taskrcond;
38            public double rmserror;
39            public double avgerror;
40            public double avgrelerror;
41            public double maxerror;
42        };
43
44
45
46
47        /*************************************************************************
48        Lagrange intepolant: generation of the model on the general grid.
49        This function has O(N^2) complexity.
50
51        INPUT PARAMETERS:
52            X   -   abscissas, array[0..N-1]
53            Y   -   function values, array[0..N-1]
54            N   -   number of points, N>=1
55
56        OIYTPUT PARAMETERS
57            P   -   barycentric model which represents Lagrange interpolant
58                    (see ratint unit info and BarycentricCalc() description for
59                    more information).
60
61          -- ALGLIB --
62             Copyright 02.12.2009 by Bochkanov Sergey
63        *************************************************************************/
64        public static void polynomialbuild(ref double[] x,
65            ref double[] y,
66            int n,
67            ref ratint.barycentricinterpolant p)
68        {
69            int j = 0;
70            int k = 0;
71            double[] w = new double[0];
72            double b = 0;
73            double a = 0;
74            double v = 0;
75            double mx = 0;
76            int i_ = 0;
77
78            System.Diagnostics.Debug.Assert(n>0, "PolIntBuild: N<=0!");
79           
80            //
81            // calculate W[j]
82            // multi-pass algorithm is used to avoid overflow
83            //
84            w = new double[n];
85            a = x[0];
86            b = x[0];
87            for(j=0; j<=n-1; j++)
88            {
89                w[j] = 1;
90                a = Math.Min(a, x[j]);
91                b = Math.Max(b, x[j]);
92            }
93            for(k=0; k<=n-1; k++)
94            {
95               
96                //
97                // W[K] is used instead of 0.0 because
98                // cycle on J does not touch K-th element
99                // and we MUST get maximum from ALL elements
100                //
101                mx = Math.Abs(w[k]);
102                for(j=0; j<=n-1; j++)
103                {
104                    if( j!=k )
105                    {
106                        v = (b-a)/(x[j]-x[k]);
107                        w[j] = w[j]*v;
108                        mx = Math.Max(mx, Math.Abs(w[j]));
109                    }
110                }
111                if( k%5==0 )
112                {
113                   
114                    //
115                    // every 5-th run we renormalize W[]
116                    //
117                    v = 1/mx;
118                    for(i_=0; i_<=n-1;i_++)
119                    {
120                        w[i_] = v*w[i_];
121                    }
122                }
123            }
124            ratint.barycentricbuildxyw(ref x, ref y, ref w, n, ref p);
125        }
126
127
128        /*************************************************************************
129        Lagrange intepolant: generation of the model on equidistant grid.
130        This function has O(N) complexity.
131
132        INPUT PARAMETERS:
133            A   -   left boundary of [A,B]
134            B   -   right boundary of [A,B]
135            Y   -   function values at the nodes, array[0..N-1]
136            N   -   number of points, N>=1
137                    for N=1 a constant model is constructed.
138
139        OIYTPUT PARAMETERS
140            P   -   barycentric model which represents Lagrange interpolant
141                    (see ratint unit info and BarycentricCalc() description for
142                    more information).
143
144          -- ALGLIB --
145             Copyright 03.12.2009 by Bochkanov Sergey
146        *************************************************************************/
147        public static void polynomialbuildeqdist(double a,
148            double b,
149            ref double[] y,
150            int n,
151            ref ratint.barycentricinterpolant p)
152        {
153            int i = 0;
154            double[] w = new double[0];
155            double[] x = new double[0];
156            double v = 0;
157
158            System.Diagnostics.Debug.Assert(n>0, "PolIntBuildEqDist: N<=0!");
159           
160            //
161            // Special case: N=1
162            //
163            if( n==1 )
164            {
165                x = new double[1];
166                w = new double[1];
167                x[0] = 0.5*(b+a);
168                w[0] = 1;
169                ratint.barycentricbuildxyw(ref x, ref y, ref w, 1, ref p);
170                return;
171            }
172           
173            //
174            // general case
175            //
176            x = new double[n];
177            w = new double[n];
178            v = 1;
179            for(i=0; i<=n-1; i++)
180            {
181                w[i] = v;
182                x[i] = a+(b-a)*i/(n-1);
183                v = -(v*(n-1-i));
184                v = v/(i+1);
185            }
186            ratint.barycentricbuildxyw(ref x, ref y, ref w, n, ref p);
187        }
188
189
190        /*************************************************************************
191        Lagrange intepolant on Chebyshev grid (first kind).
192        This function has O(N) complexity.
193
194        INPUT PARAMETERS:
195            A   -   left boundary of [A,B]
196            B   -   right boundary of [A,B]
197            Y   -   function values at the nodes, array[0..N-1],
198                    Y[I] = Y(0.5*(B+A) + 0.5*(B-A)*Cos(PI*(2*i+1)/(2*n)))
199            N   -   number of points, N>=1
200                    for N=1 a constant model is constructed.
201
202        OIYTPUT PARAMETERS
203            P   -   barycentric model which represents Lagrange interpolant
204                    (see ratint unit info and BarycentricCalc() description for
205                    more information).
206
207          -- ALGLIB --
208             Copyright 03.12.2009 by Bochkanov Sergey
209        *************************************************************************/
210        public static void polynomialbuildcheb1(double a,
211            double b,
212            ref double[] y,
213            int n,
214            ref ratint.barycentricinterpolant p)
215        {
216            int i = 0;
217            double[] w = new double[0];
218            double[] x = new double[0];
219            double v = 0;
220            double t = 0;
221
222            System.Diagnostics.Debug.Assert(n>0, "PolIntBuildCheb1: N<=0!");
223           
224            //
225            // Special case: N=1
226            //
227            if( n==1 )
228            {
229                x = new double[1];
230                w = new double[1];
231                x[0] = 0.5*(b+a);
232                w[0] = 1;
233                ratint.barycentricbuildxyw(ref x, ref y, ref w, 1, ref p);
234                return;
235            }
236           
237            //
238            // general case
239            //
240            x = new double[n];
241            w = new double[n];
242            v = 1;
243            for(i=0; i<=n-1; i++)
244            {
245                t = Math.Tan(0.5*Math.PI*(2*i+1)/(2*n));
246                w[i] = 2*v*t/(1+AP.Math.Sqr(t));
247                x[i] = 0.5*(b+a)+0.5*(b-a)*(1-AP.Math.Sqr(t))/(1+AP.Math.Sqr(t));
248                v = -v;
249            }
250            ratint.barycentricbuildxyw(ref x, ref y, ref w, n, ref p);
251        }
252
253
254        /*************************************************************************
255        Lagrange intepolant on Chebyshev grid (second kind).
256        This function has O(N) complexity.
257
258        INPUT PARAMETERS:
259            A   -   left boundary of [A,B]
260            B   -   right boundary of [A,B]
261            Y   -   function values at the nodes, array[0..N-1],
262                    Y[I] = Y(0.5*(B+A) + 0.5*(B-A)*Cos(PI*i/(n-1)))
263            N   -   number of points, N>=1
264                    for N=1 a constant model is constructed.
265
266        OIYTPUT PARAMETERS
267            P   -   barycentric model which represents Lagrange interpolant
268                    (see ratint unit info and BarycentricCalc() description for
269                    more information).
270
271          -- ALGLIB --
272             Copyright 03.12.2009 by Bochkanov Sergey
273        *************************************************************************/
274        public static void polynomialbuildcheb2(double a,
275            double b,
276            ref double[] y,
277            int n,
278            ref ratint.barycentricinterpolant p)
279        {
280            int i = 0;
281            double[] w = new double[0];
282            double[] x = new double[0];
283            double v = 0;
284            double t = 0;
285
286            System.Diagnostics.Debug.Assert(n>0, "PolIntBuildCheb2: N<=0!");
287           
288            //
289            // Special case: N=1
290            //
291            if( n==1 )
292            {
293                x = new double[1];
294                w = new double[1];
295                x[0] = 0.5*(b+a);
296                w[0] = 1;
297                ratint.barycentricbuildxyw(ref x, ref y, ref w, 1, ref p);
298                return;
299            }
300           
301            //
302            // general case
303            //
304            x = new double[n];
305            w = new double[n];
306            v = 1;
307            for(i=0; i<=n-1; i++)
308            {
309                if( i==0 | i==n-1 )
310                {
311                    w[i] = v*0.5;
312                }
313                else
314                {
315                    w[i] = v;
316                }
317                x[i] = 0.5*(b+a)+0.5*(b-a)*Math.Cos(Math.PI*i/(n-1));
318                v = -v;
319            }
320            ratint.barycentricbuildxyw(ref x, ref y, ref w, n, ref p);
321        }
322
323
324        /*************************************************************************
325        Fast equidistant polynomial interpolation function with O(N) complexity
326
327        INPUT PARAMETERS:
328            A   -   left boundary of [A,B]
329            B   -   right boundary of [A,B]
330            F   -   function values, array[0..N-1]
331            N   -   number of points on equidistant grid, N>=1
332                    for N=1 a constant model is constructed.
333            T   -   position where P(x) is calculated
334
335        RESULT
336            value of the Lagrange interpolant at T
337           
338        IMPORTANT
339            this function provides fast interface which is not overflow-safe
340            nor it is very precise.
341            the best option is to use  PolynomialBuildEqDist()/BarycentricCalc()
342            subroutines unless you are pretty sure that your data will not result
343            in overflow.
344
345          -- ALGLIB --
346             Copyright 02.12.2009 by Bochkanov Sergey
347        *************************************************************************/
348        public static double polynomialcalceqdist(double a,
349            double b,
350            ref double[] f,
351            int n,
352            double t)
353        {
354            double result = 0;
355            double s1 = 0;
356            double s2 = 0;
357            double v = 0;
358            double threshold = 0;
359            double s = 0;
360            double h = 0;
361            int i = 0;
362            int j = 0;
363            double w = 0;
364            double x = 0;
365
366            System.Diagnostics.Debug.Assert(n>0, "PolIntEqDist: N<=0!");
367            threshold = Math.Sqrt(AP.Math.MinRealNumber);
368           
369            //
370            // Special case: N=1
371            //
372            if( n==1 )
373            {
374                result = f[0];
375                return result;
376            }
377           
378            //
379            // First, decide: should we use "safe" formula (guarded
380            // against overflow) or fast one?
381            //
382            j = 0;
383            s = t-a;
384            for(i=1; i<=n-1; i++)
385            {
386                x = a+(double)(i)/((double)(n-1))*(b-a);
387                if( (double)(Math.Abs(t-x))<(double)(Math.Abs(s)) )
388                {
389                    s = t-x;
390                    j = i;
391                }
392            }
393            if( (double)(s)==(double)(0) )
394            {
395                result = f[j];
396                return result;
397            }
398            if( (double)(Math.Abs(s))>(double)(threshold) )
399            {
400               
401                //
402                // use fast formula
403                //
404                j = -1;
405                s = 1.0;
406            }
407           
408            //
409            // Calculate using safe or fast barycentric formula
410            //
411            s1 = 0;
412            s2 = 0;
413            w = 1.0;
414            h = (b-a)/(n-1);
415            for(i=0; i<=n-1; i++)
416            {
417                if( i!=j )
418                {
419                    v = s*w/(t-(a+i*h));
420                    s1 = s1+v*f[i];
421                    s2 = s2+v;
422                }
423                else
424                {
425                    v = w;
426                    s1 = s1+v*f[i];
427                    s2 = s2+v;
428                }
429                w = -(w*(n-1-i));
430                w = w/(i+1);
431            }
432            result = s1/s2;
433            return result;
434        }
435
436
437        /*************************************************************************
438        Fast polynomial interpolation function on Chebyshev points (first kind)
439        with O(N) complexity.
440
441        INPUT PARAMETERS:
442            A   -   left boundary of [A,B]
443            B   -   right boundary of [A,B]
444            F   -   function values, array[0..N-1]
445            N   -   number of points on Chebyshev grid (first kind),
446                    X[i] = 0.5*(B+A) + 0.5*(B-A)*Cos(PI*(2*i+1)/(2*n))
447                    for N=1 a constant model is constructed.
448            T   -   position where P(x) is calculated
449
450        RESULT
451            value of the Lagrange interpolant at T
452
453        IMPORTANT
454            this function provides fast interface which is not overflow-safe
455            nor it is very precise.
456            the best option is to use  PolIntBuildCheb1()/BarycentricCalc()
457            subroutines unless you are pretty sure that your data will not result
458            in overflow.
459
460          -- ALGLIB --
461             Copyright 02.12.2009 by Bochkanov Sergey
462        *************************************************************************/
463        public static double polynomialcalccheb1(double a,
464            double b,
465            ref double[] f,
466            int n,
467            double t)
468        {
469            double result = 0;
470            double s1 = 0;
471            double s2 = 0;
472            double v = 0;
473            double threshold = 0;
474            double s = 0;
475            int i = 0;
476            int j = 0;
477            double a0 = 0;
478            double delta = 0;
479            double alpha = 0;
480            double beta = 0;
481            double ca = 0;
482            double sa = 0;
483            double tempc = 0;
484            double temps = 0;
485            double x = 0;
486            double w = 0;
487            double p1 = 0;
488
489            System.Diagnostics.Debug.Assert(n>0, "PolIntCheb1: N<=0!");
490            threshold = Math.Sqrt(AP.Math.MinRealNumber);
491            t = (t-0.5*(a+b))/(0.5*(b-a));
492           
493            //
494            // Fast exit
495            //
496            if( n==1 )
497            {
498                result = f[0];
499                return result;
500            }
501           
502            //
503            // Prepare information for the recurrence formula
504            // used to calculate sin(pi*(2j+1)/(2n+2)) and
505            // cos(pi*(2j+1)/(2n+2)):
506            //
507            // A0    = pi/(2n+2)
508            // Delta = pi/(n+1)
509            // Alpha = 2 sin^2 (Delta/2)
510            // Beta  = sin(Delta)
511            //
512            // so that sin(..) = sin(A0+j*delta) and cos(..) = cos(A0+j*delta).
513            // Then we use
514            //
515            // sin(x+delta) = sin(x) - (alpha*sin(x) - beta*cos(x))
516            // cos(x+delta) = cos(x) - (alpha*cos(x) - beta*sin(x))
517            //
518            // to repeatedly calculate sin(..) and cos(..).
519            //
520            a0 = Math.PI/(2*(n-1)+2);
521            delta = 2*Math.PI/(2*(n-1)+2);
522            alpha = 2*AP.Math.Sqr(Math.Sin(delta/2));
523            beta = Math.Sin(delta);
524           
525            //
526            // First, decide: should we use "safe" formula (guarded
527            // against overflow) or fast one?
528            //
529            ca = Math.Cos(a0);
530            sa = Math.Sin(a0);
531            j = 0;
532            x = ca;
533            s = t-x;
534            for(i=1; i<=n-1; i++)
535            {
536               
537                //
538                // Next X[i]
539                //
540                temps = sa-(alpha*sa-beta*ca);
541                tempc = ca-(alpha*ca+beta*sa);
542                sa = temps;
543                ca = tempc;
544                x = ca;
545               
546                //
547                // Use X[i]
548                //
549                if( (double)(Math.Abs(t-x))<(double)(Math.Abs(s)) )
550                {
551                    s = t-x;
552                    j = i;
553                }
554            }
555            if( (double)(s)==(double)(0) )
556            {
557                result = f[j];
558                return result;
559            }
560            if( (double)(Math.Abs(s))>(double)(threshold) )
561            {
562               
563                //
564                // use fast formula
565                //
566                j = -1;
567                s = 1.0;
568            }
569           
570            //
571            // Calculate using safe or fast barycentric formula
572            //
573            s1 = 0;
574            s2 = 0;
575            ca = Math.Cos(a0);
576            sa = Math.Sin(a0);
577            p1 = 1.0;
578            for(i=0; i<=n-1; i++)
579            {
580               
581                //
582                // Calculate X[i], W[i]
583                //
584                x = ca;
585                w = p1*sa;
586               
587                //
588                // Proceed
589                //
590                if( i!=j )
591                {
592                    v = s*w/(t-x);
593                    s1 = s1+v*f[i];
594                    s2 = s2+v;
595                }
596                else
597                {
598                    v = w;
599                    s1 = s1+v*f[i];
600                    s2 = s2+v;
601                }
602               
603                //
604                // Next CA, SA, P1
605                //
606                temps = sa-(alpha*sa-beta*ca);
607                tempc = ca-(alpha*ca+beta*sa);
608                sa = temps;
609                ca = tempc;
610                p1 = -p1;
611            }
612            result = s1/s2;
613            return result;
614        }
615
616
617        /*************************************************************************
618        Fast polynomial interpolation function on Chebyshev points (second kind)
619        with O(N) complexity.
620
621        INPUT PARAMETERS:
622            A   -   left boundary of [A,B]
623            B   -   right boundary of [A,B]
624            F   -   function values, array[0..N-1]
625            N   -   number of points on Chebyshev grid (second kind),
626                    X[i] = 0.5*(B+A) + 0.5*(B-A)*Cos(PI*i/(n-1))
627                    for N=1 a constant model is constructed.
628            T   -   position where P(x) is calculated
629
630        RESULT
631            value of the Lagrange interpolant at T
632
633        IMPORTANT
634            this function provides fast interface which is not overflow-safe
635            nor it is very precise.
636            the best option is to use PolIntBuildCheb2()/BarycentricCalc()
637            subroutines unless you are pretty sure that your data will not result
638            in overflow.
639
640          -- ALGLIB --
641             Copyright 02.12.2009 by Bochkanov Sergey
642        *************************************************************************/
643        public static double polynomialcalccheb2(double a,
644            double b,
645            ref double[] f,
646            int n,
647            double t)
648        {
649            double result = 0;
650            double s1 = 0;
651            double s2 = 0;
652            double v = 0;
653            double threshold = 0;
654            double s = 0;
655            int i = 0;
656            int j = 0;
657            double a0 = 0;
658            double delta = 0;
659            double alpha = 0;
660            double beta = 0;
661            double ca = 0;
662            double sa = 0;
663            double tempc = 0;
664            double temps = 0;
665            double x = 0;
666            double w = 0;
667            double p1 = 0;
668
669            System.Diagnostics.Debug.Assert(n>0, "PolIntCheb2: N<=0!");
670            threshold = Math.Sqrt(AP.Math.MinRealNumber);
671            t = (t-0.5*(a+b))/(0.5*(b-a));
672           
673            //
674            // Fast exit
675            //
676            if( n==1 )
677            {
678                result = f[0];
679                return result;
680            }
681           
682            //
683            // Prepare information for the recurrence formula
684            // used to calculate sin(pi*i/n) and
685            // cos(pi*i/n):
686            //
687            // A0    = 0
688            // Delta = pi/n
689            // Alpha = 2 sin^2 (Delta/2)
690            // Beta  = sin(Delta)
691            //
692            // so that sin(..) = sin(A0+j*delta) and cos(..) = cos(A0+j*delta).
693            // Then we use
694            //
695            // sin(x+delta) = sin(x) - (alpha*sin(x) - beta*cos(x))
696            // cos(x+delta) = cos(x) - (alpha*cos(x) - beta*sin(x))
697            //
698            // to repeatedly calculate sin(..) and cos(..).
699            //
700            a0 = 0.0;
701            delta = Math.PI/(n-1);
702            alpha = 2*AP.Math.Sqr(Math.Sin(delta/2));
703            beta = Math.Sin(delta);
704           
705            //
706            // First, decide: should we use "safe" formula (guarded
707            // against overflow) or fast one?
708            //
709            ca = Math.Cos(a0);
710            sa = Math.Sin(a0);
711            j = 0;
712            x = ca;
713            s = t-x;
714            for(i=1; i<=n-1; i++)
715            {
716               
717                //
718                // Next X[i]
719                //
720                temps = sa-(alpha*sa-beta*ca);
721                tempc = ca-(alpha*ca+beta*sa);
722                sa = temps;
723                ca = tempc;
724                x = ca;
725               
726                //
727                // Use X[i]
728                //
729                if( (double)(Math.Abs(t-x))<(double)(Math.Abs(s)) )
730                {
731                    s = t-x;
732                    j = i;
733                }
734            }
735            if( (double)(s)==(double)(0) )
736            {
737                result = f[j];
738                return result;
739            }
740            if( (double)(Math.Abs(s))>(double)(threshold) )
741            {
742               
743                //
744                // use fast formula
745                //
746                j = -1;
747                s = 1.0;
748            }
749           
750            //
751            // Calculate using safe or fast barycentric formula
752            //
753            s1 = 0;
754            s2 = 0;
755            ca = Math.Cos(a0);
756            sa = Math.Sin(a0);
757            p1 = 1.0;
758            for(i=0; i<=n-1; i++)
759            {
760               
761                //
762                // Calculate X[i], W[i]
763                //
764                x = ca;
765                if( i==0 | i==n-1 )
766                {
767                    w = 0.5*p1;
768                }
769                else
770                {
771                    w = 1.0*p1;
772                }
773               
774                //
775                // Proceed
776                //
777                if( i!=j )
778                {
779                    v = s*w/(t-x);
780                    s1 = s1+v*f[i];
781                    s2 = s2+v;
782                }
783                else
784                {
785                    v = w;
786                    s1 = s1+v*f[i];
787                    s2 = s2+v;
788                }
789               
790                //
791                // Next CA, SA, P1
792                //
793                temps = sa-(alpha*sa-beta*ca);
794                tempc = ca-(alpha*ca+beta*sa);
795                sa = temps;
796                ca = tempc;
797                p1 = -p1;
798            }
799            result = s1/s2;
800            return result;
801        }
802
803
804        /*************************************************************************
805        Least squares fitting by polynomial.
806
807        This subroutine is "lightweight" alternative for more complex and feature-
808        rich PolynomialFitWC().  See  PolynomialFitWC() for more information about
809        subroutine parameters (we don't duplicate it here because of length)
810
811          -- ALGLIB PROJECT --
812             Copyright 12.10.2009 by Bochkanov Sergey
813        *************************************************************************/
814        public static void polynomialfit(ref double[] x,
815            ref double[] y,
816            int n,
817            int m,
818            ref int info,
819            ref ratint.barycentricinterpolant p,
820            ref polynomialfitreport rep)
821        {
822            int i = 0;
823            double[] w = new double[0];
824            double[] xc = new double[0];
825            double[] yc = new double[0];
826            int[] dc = new int[0];
827
828            if( n>0 )
829            {
830                w = new double[n];
831                for(i=0; i<=n-1; i++)
832                {
833                    w[i] = 1;
834                }
835            }
836            polynomialfitwc(x, y, ref w, n, xc, yc, ref dc, 0, m, ref info, ref p, ref rep);
837        }
838
839
840        /*************************************************************************
841        Weighted  fitting  by  Chebyshev  polynomial  in  barycentric  form,  with
842        constraints on function values or first derivatives.
843
844        Small regularizing term is used when solving constrained tasks (to improve
845        stability).
846
847        Task is linear, so linear least squares solver is used. Complexity of this
848        computational scheme is O(N*M^2), mostly dominated by least squares solver
849
850        SEE ALSO:
851            PolynomialFit()
852
853        INPUT PARAMETERS:
854            X   -   points, array[0..N-1].
855            Y   -   function values, array[0..N-1].
856            W   -   weights, array[0..N-1]
857                    Each summand in square  sum  of  approximation deviations from
858                    given  values  is  multiplied  by  the square of corresponding
859                    weight. Fill it by 1's if you don't  want  to  solve  weighted
860                    task.
861            N   -   number of points, N>0.
862            XC  -   points where polynomial values/derivatives are constrained,
863                    array[0..K-1].
864            YC  -   values of constraints, array[0..K-1]
865            DC  -   array[0..K-1], types of constraints:
866                    * DC[i]=0   means that P(XC[i])=YC[i]
867                    * DC[i]=1   means that P'(XC[i])=YC[i]
868                    SEE BELOW FOR IMPORTANT INFORMATION ON CONSTRAINTS
869            K   -   number of constraints, 0<=K<M.
870                    K=0 means no constraints (XC/YC/DC are not used in such cases)
871            M   -   number of basis functions (= polynomial_degree + 1), M>=1
872
873        OUTPUT PARAMETERS:
874            Info-   same format as in LSFitLinearW() subroutine:
875                    * Info>0    task is solved
876                    * Info<=0   an error occured:
877                                -4 means inconvergence of internal SVD
878                                -3 means inconsistent constraints
879                                -1 means another errors in parameters passed
880                                   (N<=0, for example)
881            P   -   interpolant in barycentric form.
882            Rep -   report, same format as in LSFitLinearW() subroutine.
883                    Following fields are set:
884                    * RMSError      rms error on the (X,Y).
885                    * AvgError      average error on the (X,Y).
886                    * AvgRelError   average relative error on the non-zero Y
887                    * MaxError      maximum error
888                                    NON-WEIGHTED ERRORS ARE CALCULATED
889
890        IMPORTANT:
891            this subroitine doesn't calculate task's condition number for K<>0.
892
893        SETTING CONSTRAINTS - DANGERS AND OPPORTUNITIES:
894
895        Setting constraints can lead  to undesired  results,  like ill-conditioned
896        behavior, or inconsistency being detected. From the other side,  it allows
897        us to improve quality of the fit. Here we summarize  our  experience  with
898        constrained regression splines:
899        * even simple constraints can be inconsistent, see  Wikipedia  article  on
900          this subject: http://en.wikipedia.org/wiki/Birkhoff_interpolation
901        * the  greater  is  M (given  fixed  constraints),  the  more chances that
902          constraints will be consistent
903        * in the general case, consistency of constraints is NOT GUARANTEED.
904        * in the one special cases, however, we can  guarantee  consistency.  This
905          case  is:  M>1  and constraints on the function values (NOT DERIVATIVES)
906
907        Our final recommendation is to use constraints  WHEN  AND  ONLY  when  you
908        can't solve your task without them. Anything beyond  special  cases  given
909        above is not guaranteed and may result in inconsistency.
910
911          -- ALGLIB PROJECT --
912             Copyright 10.12.2009 by Bochkanov Sergey
913        *************************************************************************/
914        public static void polynomialfitwc(double[] x,
915            double[] y,
916            ref double[] w,
917            int n,
918            double[] xc,
919            double[] yc,
920            ref int[] dc,
921            int k,
922            int m,
923            ref int info,
924            ref ratint.barycentricinterpolant p,
925            ref polynomialfitreport rep)
926        {
927            double xa = 0;
928            double xb = 0;
929            double sa = 0;
930            double sb = 0;
931            double[] xoriginal = new double[0];
932            double[] yoriginal = new double[0];
933            double[] y2 = new double[0];
934            double[] w2 = new double[0];
935            double[] tmp = new double[0];
936            double[] tmp2 = new double[0];
937            double[] tmpdiff = new double[0];
938            double[] bx = new double[0];
939            double[] by = new double[0];
940            double[] bw = new double[0];
941            double[,] fmatrix = new double[0,0];
942            double[,] cmatrix = new double[0,0];
943            int i = 0;
944            int j = 0;
945            double mx = 0;
946            double decay = 0;
947            double u = 0;
948            double v = 0;
949            double s = 0;
950            int relcnt = 0;
951            lsfit.lsfitreport lrep = new lsfit.lsfitreport();
952            int i_ = 0;
953
954            x = (double[])x.Clone();
955            y = (double[])y.Clone();
956            xc = (double[])xc.Clone();
957            yc = (double[])yc.Clone();
958
959            if( m<1 | n<1 | k<0 | k>=m )
960            {
961                info = -1;
962                return;
963            }
964            for(i=0; i<=k-1; i++)
965            {
966                info = 0;
967                if( dc[i]<0 )
968                {
969                    info = -1;
970                }
971                if( dc[i]>1 )
972                {
973                    info = -1;
974                }
975                if( info<0 )
976                {
977                    return;
978                }
979            }
980           
981            //
982            // weight decay for correct handling of task which becomes
983            // degenerate after constraints are applied
984            //
985            decay = 10000*AP.Math.MachineEpsilon;
986           
987            //
988            // Scale X, Y, XC, YC
989            //
990            lsfit.lsfitscalexy(ref x, ref y, n, ref xc, ref yc, ref dc, k, ref xa, ref xb, ref sa, ref sb, ref xoriginal, ref yoriginal);
991           
992            //
993            // allocate space, initialize/fill:
994            // * FMatrix-   values of basis functions at X[]
995            // * CMatrix-   values (derivatives) of basis functions at XC[]
996            // * fill constraints matrix
997            // * fill first N rows of design matrix with values
998            // * fill next M rows of design matrix with regularizing term
999            // * append M zeros to Y
1000            // * append M elements, mean(abs(W)) each, to W
1001            //
1002            y2 = new double[n+m];
1003            w2 = new double[n+m];
1004            tmp = new double[m];
1005            tmpdiff = new double[m];
1006            fmatrix = new double[n+m, m];
1007            if( k>0 )
1008            {
1009                cmatrix = new double[k, m+1];
1010            }
1011           
1012            //
1013            // Fill design matrix, Y2, W2:
1014            // * first N rows with basis functions for original points
1015            // * next M rows with decay terms
1016            //
1017            for(i=0; i<=n-1; i++)
1018            {
1019               
1020                //
1021                // prepare Ith row
1022                // use Tmp for calculations to avoid multidimensional arrays overhead
1023                //
1024                for(j=0; j<=m-1; j++)
1025                {
1026                    if( j==0 )
1027                    {
1028                        tmp[j] = 1;
1029                    }
1030                    else
1031                    {
1032                        if( j==1 )
1033                        {
1034                            tmp[j] = x[i];
1035                        }
1036                        else
1037                        {
1038                            tmp[j] = 2*x[i]*tmp[j-1]-tmp[j-2];
1039                        }
1040                    }
1041                }
1042                for(i_=0; i_<=m-1;i_++)
1043                {
1044                    fmatrix[i,i_] = tmp[i_];
1045                }
1046            }
1047            for(i=0; i<=m-1; i++)
1048            {
1049                for(j=0; j<=m-1; j++)
1050                {
1051                    if( i==j )
1052                    {
1053                        fmatrix[n+i,j] = decay;
1054                    }
1055                    else
1056                    {
1057                        fmatrix[n+i,j] = 0;
1058                    }
1059                }
1060            }
1061            for(i_=0; i_<=n-1;i_++)
1062            {
1063                y2[i_] = y[i_];
1064            }
1065            for(i_=0; i_<=n-1;i_++)
1066            {
1067                w2[i_] = w[i_];
1068            }
1069            mx = 0;
1070            for(i=0; i<=n-1; i++)
1071            {
1072                mx = mx+Math.Abs(w[i]);
1073            }
1074            mx = mx/n;
1075            for(i=0; i<=m-1; i++)
1076            {
1077                y2[n+i] = 0;
1078                w2[n+i] = mx;
1079            }
1080           
1081            //
1082            // fill constraints matrix
1083            //
1084            for(i=0; i<=k-1; i++)
1085            {
1086               
1087                //
1088                // prepare Ith row
1089                // use Tmp for basis function values,
1090                // TmpDiff for basos function derivatives
1091                //
1092                for(j=0; j<=m-1; j++)
1093                {
1094                    if( j==0 )
1095                    {
1096                        tmp[j] = 1;
1097                        tmpdiff[j] = 0;
1098                    }
1099                    else
1100                    {
1101                        if( j==1 )
1102                        {
1103                            tmp[j] = xc[i];
1104                            tmpdiff[j] = 1;
1105                        }
1106                        else
1107                        {
1108                            tmp[j] = 2*xc[i]*tmp[j-1]-tmp[j-2];
1109                            tmpdiff[j] = 2*(tmp[j-1]+xc[i]*tmpdiff[j-1])-tmpdiff[j-2];
1110                        }
1111                    }
1112                }
1113                if( dc[i]==0 )
1114                {
1115                    for(i_=0; i_<=m-1;i_++)
1116                    {
1117                        cmatrix[i,i_] = tmp[i_];
1118                    }
1119                }
1120                if( dc[i]==1 )
1121                {
1122                    for(i_=0; i_<=m-1;i_++)
1123                    {
1124                        cmatrix[i,i_] = tmpdiff[i_];
1125                    }
1126                }
1127                cmatrix[i,m] = yc[i];
1128            }
1129           
1130            //
1131            // Solve constrained task
1132            //
1133            if( k>0 )
1134            {
1135               
1136                //
1137                // solve using regularization
1138                //
1139                lsfit.lsfitlinearwc(y2, ref w2, ref fmatrix, cmatrix, n+m, m, k, ref info, ref tmp, ref lrep);
1140            }
1141            else
1142            {
1143               
1144                //
1145                // no constraints, no regularization needed
1146                //
1147                lsfit.lsfitlinearwc(y, ref w, ref fmatrix, cmatrix, n, m, 0, ref info, ref tmp, ref lrep);
1148            }
1149            if( info<0 )
1150            {
1151                return;
1152            }
1153           
1154            //
1155            // Generate barycentric model and scale it
1156            // * BX, BY store barycentric model nodes
1157            // * FMatrix is reused (remember - it is at least MxM, what we need)
1158            //
1159            // Model intialization is done in O(M^2). In principle, it can be
1160            // done in O(M*log(M)), but before it we solved task with O(N*M^2)
1161            // complexity, so it is only a small amount of total time spent.
1162            //
1163            bx = new double[m];
1164            by = new double[m];
1165            bw = new double[m];
1166            tmp2 = new double[m];
1167            s = 1;
1168            for(i=0; i<=m-1; i++)
1169            {
1170                if( m!=1 )
1171                {
1172                    u = Math.Cos(Math.PI*i/(m-1));
1173                }
1174                else
1175                {
1176                    u = 0;
1177                }
1178                v = 0;
1179                for(j=0; j<=m-1; j++)
1180                {
1181                    if( j==0 )
1182                    {
1183                        tmp2[j] = 1;
1184                    }
1185                    else
1186                    {
1187                        if( j==1 )
1188                        {
1189                            tmp2[j] = u;
1190                        }
1191                        else
1192                        {
1193                            tmp2[j] = 2*u*tmp2[j-1]-tmp2[j-2];
1194                        }
1195                    }
1196                    v = v+tmp[j]*tmp2[j];
1197                }
1198                bx[i] = u;
1199                by[i] = v;
1200                bw[i] = s;
1201                if( i==0 | i==m-1 )
1202                {
1203                    bw[i] = 0.5*bw[i];
1204                }
1205                s = -s;
1206            }
1207            ratint.barycentricbuildxyw(ref bx, ref by, ref bw, m, ref p);
1208            ratint.barycentriclintransx(ref p, 2/(xb-xa), -((xa+xb)/(xb-xa)));
1209            ratint.barycentriclintransy(ref p, sb-sa, sa);
1210           
1211            //
1212            // Scale absolute errors obtained from LSFitLinearW.
1213            // Relative error should be calculated separately
1214            // (because of shifting/scaling of the task)
1215            //
1216            rep.taskrcond = lrep.taskrcond;
1217            rep.rmserror = lrep.rmserror*(sb-sa);
1218            rep.avgerror = lrep.avgerror*(sb-sa);
1219            rep.maxerror = lrep.maxerror*(sb-sa);
1220            rep.avgrelerror = 0;
1221            relcnt = 0;
1222            for(i=0; i<=n-1; i++)
1223            {
1224                if( (double)(yoriginal[i])!=(double)(0) )
1225                {
1226                    rep.avgrelerror = rep.avgrelerror+Math.Abs(ratint.barycentriccalc(ref p, xoriginal[i])-yoriginal[i])/Math.Abs(yoriginal[i]);
1227                    relcnt = relcnt+1;
1228                }
1229            }
1230            if( relcnt!=0 )
1231            {
1232                rep.avgrelerror = rep.avgrelerror/relcnt;
1233            }
1234        }
1235    }
1236}
Note: See TracBrowser for help on using the repository browser.