Free cookie consent management tool by TermsFeed Policy Generator

source: tags/3.3.1/HeuristicLab.ExtLibs/HeuristicLab.ALGLIB/2.5.0/ALGLIB-2.5.0/spline1d.cs @ 4537

Last change on this file since 4537 was 3839, checked in by mkommend, 14 years ago

implemented first version of LR (ticket #1012)

File size: 60.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 spline1d
26    {
27        /*************************************************************************
28        1-dimensional spline inteprolant
29        *************************************************************************/
30        public struct spline1dinterpolant
31        {
32            public int n;
33            public int k;
34            public double[] x;
35            public double[] c;
36        };
37
38
39        /*************************************************************************
40        Spline fitting report:
41            TaskRCond       reciprocal of task's condition number
42            RMSError        RMS error
43            AvgError        average error
44            AvgRelError     average relative error (for non-zero Y[I])
45            MaxError        maximum error
46        *************************************************************************/
47        public struct spline1dfitreport
48        {
49            public double taskrcond;
50            public double rmserror;
51            public double avgerror;
52            public double avgrelerror;
53            public double maxerror;
54        };
55
56
57
58
59        public const int spline1dvnum = 11;
60
61
62        /*************************************************************************
63        This subroutine builds linear spline interpolant
64
65        INPUT PARAMETERS:
66            X   -   spline nodes, array[0..N-1]
67            Y   -   function values, array[0..N-1]
68            N   -   points count, N>=2
69           
70        OUTPUT PARAMETERS:
71            C   -   spline interpolant
72
73          -- ALGLIB PROJECT --
74             Copyright 24.06.2007 by Bochkanov Sergey
75        *************************************************************************/
76        public static void spline1dbuildlinear(double[] x,
77            double[] y,
78            int n,
79            ref spline1dinterpolant c)
80        {
81            int i = 0;
82
83            x = (double[])x.Clone();
84            y = (double[])y.Clone();
85
86            System.Diagnostics.Debug.Assert(n>1, "Spline1DBuildLinear: N<2!");
87           
88            //
89            // Sort points
90            //
91            heapsortpoints(ref x, ref y, n);
92           
93            //
94            // Build
95            //
96            c.n = n;
97            c.k = 3;
98            c.x = new double[n];
99            c.c = new double[4*(n-1)];
100            for(i=0; i<=n-1; i++)
101            {
102                c.x[i] = x[i];
103            }
104            for(i=0; i<=n-2; i++)
105            {
106                c.c[4*i+0] = y[i];
107                c.c[4*i+1] = (y[i+1]-y[i])/(x[i+1]-x[i]);
108                c.c[4*i+2] = 0;
109                c.c[4*i+3] = 0;
110            }
111        }
112
113
114        /*************************************************************************
115        This subroutine builds cubic spline interpolant.
116
117        INPUT PARAMETERS:
118            X           -   spline nodes, array[0..N-1]
119            Y           -   function values, array[0..N-1]
120            N           -   points count, N>=2
121            BoundLType  -   boundary condition type for the left boundary
122            BoundL      -   left boundary condition (first or second derivative,
123                            depending on the BoundLType)
124            BoundRType  -   boundary condition type for the right boundary
125            BoundR      -   right boundary condition (first or second derivative,
126                            depending on the BoundRType)
127
128        OUTPUT PARAMETERS:
129            C           -   spline interpolant
130                           
131        The BoundLType/BoundRType parameters can have the following values:
132            * 0, which  corresponds  to  the  parabolically   terminated  spline
133                 (BoundL/BoundR are ignored).
134            * 1, which corresponds to the first derivative boundary condition
135            * 2, which corresponds to the second derivative boundary condition
136
137          -- ALGLIB PROJECT --
138             Copyright 23.06.2007 by Bochkanov Sergey
139        *************************************************************************/
140        public static void spline1dbuildcubic(double[] x,
141            double[] y,
142            int n,
143            int boundltype,
144            double boundl,
145            int boundrtype,
146            double boundr,
147            ref spline1dinterpolant c)
148        {
149            double[] a1 = new double[0];
150            double[] a2 = new double[0];
151            double[] a3 = new double[0];
152            double[] b = new double[0];
153            double[] d = new double[0];
154            int i = 0;
155
156            x = (double[])x.Clone();
157            y = (double[])y.Clone();
158
159            System.Diagnostics.Debug.Assert(n>=2, "BuildCubicSpline: N<2!");
160            System.Diagnostics.Debug.Assert(boundltype==0 | boundltype==1 | boundltype==2, "BuildCubicSpline: incorrect BoundLType!");
161            System.Diagnostics.Debug.Assert(boundrtype==0 | boundrtype==1 | boundrtype==2, "BuildCubicSpline: incorrect BoundRType!");
162            a1 = new double[n];
163            a2 = new double[n];
164            a3 = new double[n];
165            b = new double[n];
166           
167            //
168            // Special case:
169            // * N=2
170            // * parabolic terminated boundary condition on both ends
171            //
172            if( n==2 & boundltype==0 & boundrtype==0 )
173            {
174               
175                //
176                // Change task type
177                //
178                boundltype = 2;
179                boundl = 0;
180                boundrtype = 2;
181                boundr = 0;
182            }
183           
184            //
185            //
186            // Sort points
187            //
188            heapsortpoints(ref x, ref y, n);
189           
190            //
191            // Left boundary conditions
192            //
193            if( boundltype==0 )
194            {
195                a1[0] = 0;
196                a2[0] = 1;
197                a3[0] = 1;
198                b[0] = 2*(y[1]-y[0])/(x[1]-x[0]);
199            }
200            if( boundltype==1 )
201            {
202                a1[0] = 0;
203                a2[0] = 1;
204                a3[0] = 0;
205                b[0] = boundl;
206            }
207            if( boundltype==2 )
208            {
209                a1[0] = 0;
210                a2[0] = 2;
211                a3[0] = 1;
212                b[0] = 3*(y[1]-y[0])/(x[1]-x[0])-0.5*boundl*(x[1]-x[0]);
213            }
214           
215            //
216            // Central conditions
217            //
218            for(i=1; i<=n-2; i++)
219            {
220                a1[i] = x[i+1]-x[i];
221                a2[i] = 2*(x[i+1]-x[i-1]);
222                a3[i] = x[i]-x[i-1];
223                b[i] = 3*(y[i]-y[i-1])/(x[i]-x[i-1])*(x[i+1]-x[i])+3*(y[i+1]-y[i])/(x[i+1]-x[i])*(x[i]-x[i-1]);
224            }
225           
226            //
227            // Right boundary conditions
228            //
229            if( boundrtype==0 )
230            {
231                a1[n-1] = 1;
232                a2[n-1] = 1;
233                a3[n-1] = 0;
234                b[n-1] = 2*(y[n-1]-y[n-2])/(x[n-1]-x[n-2]);
235            }
236            if( boundrtype==1 )
237            {
238                a1[n-1] = 0;
239                a2[n-1] = 1;
240                a3[n-1] = 0;
241                b[n-1] = boundr;
242            }
243            if( boundrtype==2 )
244            {
245                a1[n-1] = 1;
246                a2[n-1] = 2;
247                a3[n-1] = 0;
248                b[n-1] = 3*(y[n-1]-y[n-2])/(x[n-1]-x[n-2])+0.5*boundr*(x[n-1]-x[n-2]);
249            }
250           
251            //
252            // Solve
253            //
254            solvetridiagonal(a1, a2, a3, b, n, ref d);
255           
256            //
257            // Now problem is reduced to the cubic Hermite spline
258            //
259            spline1dbuildhermite(x, y, d, n, ref c);
260        }
261
262
263        /*************************************************************************
264        This subroutine builds Hermite spline interpolant.
265
266        INPUT PARAMETERS:
267            X           -   spline nodes, array[0..N-1]
268            Y           -   function values, array[0..N-1]
269            D           -   derivatives, array[0..N-1]
270            N           -   points count, N>=2
271
272        OUTPUT PARAMETERS:
273            C           -   spline interpolant.
274
275          -- ALGLIB PROJECT --
276             Copyright 23.06.2007 by Bochkanov Sergey
277        *************************************************************************/
278        public static void spline1dbuildhermite(double[] x,
279            double[] y,
280            double[] d,
281            int n,
282            ref spline1dinterpolant c)
283        {
284            int i = 0;
285            double delta = 0;
286            double delta2 = 0;
287            double delta3 = 0;
288
289            x = (double[])x.Clone();
290            y = (double[])y.Clone();
291            d = (double[])d.Clone();
292
293            System.Diagnostics.Debug.Assert(n>=2, "BuildHermiteSpline: N<2!");
294           
295            //
296            // Sort points
297            //
298            heapsortdpoints(ref x, ref y, ref d, n);
299           
300            //
301            // Build
302            //
303            c.x = new double[n];
304            c.c = new double[4*(n-1)];
305            c.k = 3;
306            c.n = n;
307            for(i=0; i<=n-1; i++)
308            {
309                c.x[i] = x[i];
310            }
311            for(i=0; i<=n-2; i++)
312            {
313                delta = x[i+1]-x[i];
314                delta2 = AP.Math.Sqr(delta);
315                delta3 = delta*delta2;
316                c.c[4*i+0] = y[i];
317                c.c[4*i+1] = d[i];
318                c.c[4*i+2] = (3*(y[i+1]-y[i])-2*d[i]*delta-d[i+1]*delta)/delta2;
319                c.c[4*i+3] = (2*(y[i]-y[i+1])+d[i]*delta+d[i+1]*delta)/delta3;
320            }
321        }
322
323
324        /*************************************************************************
325        This subroutine builds Akima spline interpolant
326
327        INPUT PARAMETERS:
328            X           -   spline nodes, array[0..N-1]
329            Y           -   function values, array[0..N-1]
330            N           -   points count, N>=5
331
332        OUTPUT PARAMETERS:
333            C           -   spline interpolant
334
335          -- ALGLIB PROJECT --
336             Copyright 24.06.2007 by Bochkanov Sergey
337        *************************************************************************/
338        public static void spline1dbuildakima(double[] x,
339            double[] y,
340            int n,
341            ref spline1dinterpolant c)
342        {
343            int i = 0;
344            double[] d = new double[0];
345            double[] w = new double[0];
346            double[] diff = new double[0];
347
348            x = (double[])x.Clone();
349            y = (double[])y.Clone();
350
351            System.Diagnostics.Debug.Assert(n>=5, "BuildAkimaSpline: N<5!");
352           
353            //
354            // Sort points
355            //
356            heapsortpoints(ref x, ref y, n);
357           
358            //
359            // Prepare W (weights), Diff (divided differences)
360            //
361            w = new double[n-1];
362            diff = new double[n-1];
363            for(i=0; i<=n-2; i++)
364            {
365                diff[i] = (y[i+1]-y[i])/(x[i+1]-x[i]);
366            }
367            for(i=1; i<=n-2; i++)
368            {
369                w[i] = Math.Abs(diff[i]-diff[i-1]);
370            }
371           
372            //
373            // Prepare Hermite interpolation scheme
374            //
375            d = new double[n];
376            for(i=2; i<=n-3; i++)
377            {
378                if( (double)(Math.Abs(w[i-1])+Math.Abs(w[i+1]))!=(double)(0) )
379                {
380                    d[i] = (w[i+1]*diff[i-1]+w[i-1]*diff[i])/(w[i+1]+w[i-1]);
381                }
382                else
383                {
384                    d[i] = ((x[i+1]-x[i])*diff[i-1]+(x[i]-x[i-1])*diff[i])/(x[i+1]-x[i-1]);
385                }
386            }
387            d[0] = diffthreepoint(x[0], x[0], y[0], x[1], y[1], x[2], y[2]);
388            d[1] = diffthreepoint(x[1], x[0], y[0], x[1], y[1], x[2], y[2]);
389            d[n-2] = diffthreepoint(x[n-2], x[n-3], y[n-3], x[n-2], y[n-2], x[n-1], y[n-1]);
390            d[n-1] = diffthreepoint(x[n-1], x[n-3], y[n-3], x[n-2], y[n-2], x[n-1], y[n-1]);
391           
392            //
393            // Build Akima spline using Hermite interpolation scheme
394            //
395            spline1dbuildhermite(x, y, d, n, ref c);
396        }
397
398
399        /*************************************************************************
400        Weighted fitting by cubic  spline,  with constraints on function values or
401        derivatives.
402
403        Equidistant grid with M-2 nodes on [min(x,xc),max(x,xc)] is  used to build
404        basis functions. Basis functions are cubic splines with continuous  second
405        derivatives  and  non-fixed first  derivatives  at  interval  ends.  Small
406        regularizing term is used  when  solving  constrained  tasks  (to  improve
407        stability).
408
409        Task is linear, so linear least squares solver is used. Complexity of this
410        computational scheme is O(N*M^2), mostly dominated by least squares solver
411
412        SEE ALSO
413            Spline1DFitHermiteWC()  -   fitting by Hermite splines (more flexible,
414                                        less smooth)
415            Spline1DFitCubic()      -   "lightweight" fitting  by  cubic  splines,
416                                        without invididual weights and constraints
417
418        INPUT PARAMETERS:
419            X   -   points, array[0..N-1].
420            Y   -   function values, array[0..N-1].
421            W   -   weights, array[0..N-1]
422                    Each summand in square  sum  of  approximation deviations from
423                    given  values  is  multiplied  by  the square of corresponding
424                    weight. Fill it by 1's if you don't  want  to  solve  weighted
425                    task.
426            N   -   number of points, N>0.
427            XC  -   points where spline values/derivatives are constrained,
428                    array[0..K-1].
429            YC  -   values of constraints, array[0..K-1]
430            DC  -   array[0..K-1], types of constraints:
431                    * DC[i]=0   means that S(XC[i])=YC[i]
432                    * DC[i]=1   means that S'(XC[i])=YC[i]
433                    SEE BELOW FOR IMPORTANT INFORMATION ON CONSTRAINTS
434            K   -   number of constraints, 0<=K<M.
435                    K=0 means no constraints (XC/YC/DC are not used in such cases)
436            M   -   number of basis functions ( = number_of_nodes+2), M>=4.
437
438        OUTPUT PARAMETERS:
439            Info-   same format as in LSFitLinearWC() subroutine.
440                    * Info>0    task is solved
441                    * Info<=0   an error occured:
442                                -4 means inconvergence of internal SVD
443                                -3 means inconsistent constraints
444                                -1 means another errors in parameters passed
445                                   (N<=0, for example)
446            S   -   spline interpolant.
447            Rep -   report, same format as in LSFitLinearWC() subroutine.
448                    Following fields are set:
449                    * RMSError      rms error on the (X,Y).
450                    * AvgError      average error on the (X,Y).
451                    * AvgRelError   average relative error on the non-zero Y
452                    * MaxError      maximum error
453                                    NON-WEIGHTED ERRORS ARE CALCULATED
454
455        IMPORTANT:
456            this subroitine doesn't calculate task's condition number for K<>0.
457
458        SETTING CONSTRAINTS - DANGERS AND OPPORTUNITIES:
459
460        Setting constraints can lead  to undesired  results,  like ill-conditioned
461        behavior, or inconsistency being detected. From the other side,  it allows
462        us to improve quality of the fit. Here we summarize  our  experience  with
463        constrained regression splines:
464        * excessive constraints can be inconsistent. Splines are  piecewise  cubic
465          functions, and it is easy to create an example, where  large  number  of
466          constraints  concentrated  in  small  area will result in inconsistency.
467          Just because spline is not flexible enough to satisfy all of  them.  And
468          same constraints spread across the  [min(x),max(x)]  will  be  perfectly
469          consistent.
470        * the more evenly constraints are spread across [min(x),max(x)],  the more
471          chances that they will be consistent
472        * the  greater  is  M (given  fixed  constraints),  the  more chances that
473          constraints will be consistent
474        * in the general case, consistency of constraints IS NOT GUARANTEED.
475        * in the several special cases, however, we CAN guarantee consistency.
476        * one of this cases is constraints  on  the  function  values  AND/OR  its
477          derivatives at the interval boundaries.
478        * another  special  case  is ONE constraint on the function value (OR, but
479          not AND, derivative) anywhere in the interval
480
481        Our final recommendation is to use constraints  WHEN  AND  ONLY  WHEN  you
482        can't solve your task without them. Anything beyond  special  cases  given
483        above is not guaranteed and may result in inconsistency.
484
485
486          -- ALGLIB PROJECT --
487             Copyright 18.08.2009 by Bochkanov Sergey
488        *************************************************************************/
489        public static void spline1dfitcubicwc(ref double[] x,
490            ref double[] y,
491            ref double[] w,
492            int n,
493            ref double[] xc,
494            ref double[] yc,
495            ref int[] dc,
496            int k,
497            int m,
498            ref int info,
499            ref spline1dinterpolant s,
500            ref spline1dfitreport rep)
501        {
502            spline1dfitinternal(0, x, y, ref w, n, xc, yc, ref dc, k, m, ref info, ref s, ref rep);
503        }
504
505
506        /*************************************************************************
507        Weighted  fitting  by Hermite spline,  with constraints on function values
508        or first derivatives.
509
510        Equidistant grid with M nodes on [min(x,xc),max(x,xc)] is  used  to  build
511        basis functions. Basis functions are Hermite splines.  Small  regularizing
512        term is used when solving constrained tasks (to improve stability).
513
514        Task is linear, so linear least squares solver is used. Complexity of this
515        computational scheme is O(N*M^2), mostly dominated by least squares solver
516
517        SEE ALSO
518            Spline1DFitCubicWC()    -   fitting by Cubic splines (less flexible,
519                                        more smooth)
520            Spline1DFitHermite()    -   "lightweight" Hermite fitting, without
521                                        invididual weights and constraints
522
523        INPUT PARAMETERS:
524            X   -   points, array[0..N-1].
525            Y   -   function values, array[0..N-1].
526            W   -   weights, array[0..N-1]
527                    Each summand in square  sum  of  approximation deviations from
528                    given  values  is  multiplied  by  the square of corresponding
529                    weight. Fill it by 1's if you don't  want  to  solve  weighted
530                    task.
531            N   -   number of points, N>0.
532            XC  -   points where spline values/derivatives are constrained,
533                    array[0..K-1].
534            YC  -   values of constraints, array[0..K-1]
535            DC  -   array[0..K-1], types of constraints:
536                    * DC[i]=0   means that S(XC[i])=YC[i]
537                    * DC[i]=1   means that S'(XC[i])=YC[i]
538                    SEE BELOW FOR IMPORTANT INFORMATION ON CONSTRAINTS
539            K   -   number of constraints, 0<=K<M.
540                    K=0 means no constraints (XC/YC/DC are not used in such cases)
541            M   -   number of basis functions (= 2 * number of nodes),
542                    M>=4,
543                    M IS EVEN!
544
545        OUTPUT PARAMETERS:
546            Info-   same format as in LSFitLinearW() subroutine:
547                    * Info>0    task is solved
548                    * Info<=0   an error occured:
549                                -4 means inconvergence of internal SVD
550                                -3 means inconsistent constraints
551                                -2 means odd M was passed (which is not supported)
552                                -1 means another errors in parameters passed
553                                   (N<=0, for example)
554            S   -   spline interpolant.
555            Rep -   report, same format as in LSFitLinearW() subroutine.
556                    Following fields are set:
557                    * RMSError      rms error on the (X,Y).
558                    * AvgError      average error on the (X,Y).
559                    * AvgRelError   average relative error on the non-zero Y
560                    * MaxError      maximum error
561                                    NON-WEIGHTED ERRORS ARE CALCULATED
562
563        IMPORTANT:
564            this subroitine doesn't calculate task's condition number for K<>0.
565
566        IMPORTANT:
567            this subroitine supports only even M's
568
569        SETTING CONSTRAINTS - DANGERS AND OPPORTUNITIES:
570
571        Setting constraints can lead  to undesired  results,  like ill-conditioned
572        behavior, or inconsistency being detected. From the other side,  it allows
573        us to improve quality of the fit. Here we summarize  our  experience  with
574        constrained regression splines:
575        * excessive constraints can be inconsistent. Splines are  piecewise  cubic
576          functions, and it is easy to create an example, where  large  number  of
577          constraints  concentrated  in  small  area will result in inconsistency.
578          Just because spline is not flexible enough to satisfy all of  them.  And
579          same constraints spread across the  [min(x),max(x)]  will  be  perfectly
580          consistent.
581        * the more evenly constraints are spread across [min(x),max(x)],  the more
582          chances that they will be consistent
583        * the  greater  is  M (given  fixed  constraints),  the  more chances that
584          constraints will be consistent
585        * in the general case, consistency of constraints is NOT GUARANTEED.
586        * in the several special cases, however, we can guarantee consistency.
587        * one of this cases is  M>=4  and   constraints  on   the  function  value
588          (AND/OR its derivative) at the interval boundaries.
589        * another special case is M>=4  and  ONE  constraint on the function value
590          (OR, BUT NOT AND, derivative) anywhere in [min(x),max(x)]
591
592        Our final recommendation is to use constraints  WHEN  AND  ONLY  when  you
593        can't solve your task without them. Anything beyond  special  cases  given
594        above is not guaranteed and may result in inconsistency.
595
596          -- ALGLIB PROJECT --
597             Copyright 18.08.2009 by Bochkanov Sergey
598        *************************************************************************/
599        public static void spline1dfithermitewc(ref double[] x,
600            ref double[] y,
601            ref double[] w,
602            int n,
603            ref double[] xc,
604            ref double[] yc,
605            ref int[] dc,
606            int k,
607            int m,
608            ref int info,
609            ref spline1dinterpolant s,
610            ref spline1dfitreport rep)
611        {
612            spline1dfitinternal(1, x, y, ref w, n, xc, yc, ref dc, k, m, ref info, ref s, ref rep);
613        }
614
615
616        /*************************************************************************
617        Least squares fitting by cubic spline.
618
619        This subroutine is "lightweight" alternative for more complex and feature-
620        rich Spline1DFitCubicWC().  See  Spline1DFitCubicWC() for more information
621        about subroutine parameters (we don't duplicate it here because of length)
622
623          -- ALGLIB PROJECT --
624             Copyright 18.08.2009 by Bochkanov Sergey
625        *************************************************************************/
626        public static void spline1dfitcubic(ref double[] x,
627            ref double[] y,
628            int n,
629            int m,
630            ref int info,
631            ref spline1dinterpolant s,
632            ref spline1dfitreport rep)
633        {
634            int i = 0;
635            double[] w = new double[0];
636            double[] xc = new double[0];
637            double[] yc = new double[0];
638            int[] dc = new int[0];
639
640            if( n>0 )
641            {
642                w = new double[n];
643                for(i=0; i<=n-1; i++)
644                {
645                    w[i] = 1;
646                }
647            }
648            spline1dfitcubicwc(ref x, ref y, ref w, n, ref xc, ref yc, ref dc, 0, m, ref info, ref s, ref rep);
649        }
650
651
652        /*************************************************************************
653        Least squares fitting by Hermite spline.
654
655        This subroutine is "lightweight" alternative for more complex and feature-
656        rich Spline1DFitHermiteWC().  See Spline1DFitHermiteWC()  description  for
657        more information about subroutine parameters (we don't duplicate  it  here
658        because of length).
659
660          -- ALGLIB PROJECT --
661             Copyright 18.08.2009 by Bochkanov Sergey
662        *************************************************************************/
663        public static void spline1dfithermite(ref double[] x,
664            ref double[] y,
665            int n,
666            int m,
667            ref int info,
668            ref spline1dinterpolant s,
669            ref spline1dfitreport rep)
670        {
671            int i = 0;
672            double[] w = new double[0];
673            double[] xc = new double[0];
674            double[] yc = new double[0];
675            int[] dc = new int[0];
676
677            if( n>0 )
678            {
679                w = new double[n];
680                for(i=0; i<=n-1; i++)
681                {
682                    w[i] = 1;
683                }
684            }
685            spline1dfithermitewc(ref x, ref y, ref w, n, ref xc, ref yc, ref dc, 0, m, ref info, ref s, ref rep);
686        }
687
688
689        /*************************************************************************
690        This subroutine calculates the value of the spline at the given point X.
691
692        INPUT PARAMETERS:
693            C   -   spline interpolant
694            X   -   point
695
696        Result:
697            S(x)
698
699          -- ALGLIB PROJECT --
700             Copyright 23.06.2007 by Bochkanov Sergey
701        *************************************************************************/
702        public static double spline1dcalc(ref spline1dinterpolant c,
703            double x)
704        {
705            double result = 0;
706            int l = 0;
707            int r = 0;
708            int m = 0;
709
710            System.Diagnostics.Debug.Assert(c.k==3, "Spline1DCalc: internal error");
711           
712            //
713            // Binary search in the [ x[0], ..., x[n-2] ] (x[n-1] is not included)
714            //
715            l = 0;
716            r = c.n-2+1;
717            while( l!=r-1 )
718            {
719                m = (l+r)/2;
720                if( (double)(c.x[m])>=(double)(x) )
721                {
722                    r = m;
723                }
724                else
725                {
726                    l = m;
727                }
728            }
729           
730            //
731            // Interpolation
732            //
733            x = x-c.x[l];
734            m = 4*l;
735            result = c.c[m]+x*(c.c[m+1]+x*(c.c[m+2]+x*c.c[m+3]));
736            return result;
737        }
738
739
740        /*************************************************************************
741        This subroutine differentiates the spline.
742
743        INPUT PARAMETERS:
744            C   -   spline interpolant.
745            X   -   point
746
747        Result:
748            S   -   S(x)
749            DS  -   S'(x)
750            D2S -   S''(x)
751
752          -- ALGLIB PROJECT --
753             Copyright 24.06.2007 by Bochkanov Sergey
754        *************************************************************************/
755        public static void spline1ddiff(ref spline1dinterpolant c,
756            double x,
757            ref double s,
758            ref double ds,
759            ref double d2s)
760        {
761            int l = 0;
762            int r = 0;
763            int m = 0;
764
765            System.Diagnostics.Debug.Assert(c.k==3, "Spline1DCalc: internal error");
766           
767            //
768            // Binary search
769            //
770            l = 0;
771            r = c.n-2+1;
772            while( l!=r-1 )
773            {
774                m = (l+r)/2;
775                if( (double)(c.x[m])>=(double)(x) )
776                {
777                    r = m;
778                }
779                else
780                {
781                    l = m;
782                }
783            }
784           
785            //
786            // Differentiation
787            //
788            x = x-c.x[l];
789            m = 4*l;
790            s = c.c[m]+x*(c.c[m+1]+x*(c.c[m+2]+x*c.c[m+3]));
791            ds = c.c[m+1]+2*x*c.c[m+2]+3*AP.Math.Sqr(x)*c.c[m+3];
792            d2s = 2*c.c[m+2]+6*x*c.c[m+3];
793        }
794
795
796        /*************************************************************************
797        This subroutine makes the copy of the spline.
798
799        INPUT PARAMETERS:
800            C   -   spline interpolant.
801
802        Result:
803            CC  -   spline copy
804
805          -- ALGLIB PROJECT --
806             Copyright 29.06.2007 by Bochkanov Sergey
807        *************************************************************************/
808        public static void spline1dcopy(ref spline1dinterpolant c,
809            ref spline1dinterpolant cc)
810        {
811            int i_ = 0;
812
813            cc.n = c.n;
814            cc.k = c.k;
815            cc.x = new double[cc.n];
816            for(i_=0; i_<=cc.n-1;i_++)
817            {
818                cc.x[i_] = c.x[i_];
819            }
820            cc.c = new double[(cc.k+1)*(cc.n-1)];
821            for(i_=0; i_<=(cc.k+1)*(cc.n-1)-1;i_++)
822            {
823                cc.c[i_] = c.c[i_];
824            }
825        }
826
827
828        /*************************************************************************
829        Serialization of the spline interpolant
830
831        INPUT PARAMETERS:
832            B   -   spline interpolant
833
834        OUTPUT PARAMETERS:
835            RA      -   array of real numbers which contains interpolant,
836                        array[0..RLen-1]
837            RLen    -   RA lenght
838
839          -- ALGLIB --
840             Copyright 17.08.2009 by Bochkanov Sergey
841        *************************************************************************/
842        public static void spline1dserialize(ref spline1dinterpolant c,
843            ref double[] ra,
844            ref int ralen)
845        {
846            int i_ = 0;
847            int i1_ = 0;
848
849            ralen = 2+2+c.n+(c.k+1)*(c.n-1);
850            ra = new double[ralen];
851            ra[0] = ralen;
852            ra[1] = spline1dvnum;
853            ra[2] = c.n;
854            ra[3] = c.k;
855            i1_ = (0) - (4);
856            for(i_=4; i_<=4+c.n-1;i_++)
857            {
858                ra[i_] = c.x[i_+i1_];
859            }
860            i1_ = (0) - (4+c.n);
861            for(i_=4+c.n; i_<=4+c.n+(c.k+1)*(c.n-1)-1;i_++)
862            {
863                ra[i_] = c.c[i_+i1_];
864            }
865        }
866
867
868        /*************************************************************************
869        Unserialization of the spline interpolant
870
871        INPUT PARAMETERS:
872            RA  -   array of real numbers which contains interpolant,
873
874        OUTPUT PARAMETERS:
875            B   -   spline interpolant
876
877          -- ALGLIB --
878             Copyright 17.08.2009 by Bochkanov Sergey
879        *************************************************************************/
880        public static void spline1dunserialize(ref double[] ra,
881            ref spline1dinterpolant c)
882        {
883            int i_ = 0;
884            int i1_ = 0;
885
886            System.Diagnostics.Debug.Assert((int)Math.Round(ra[1])==spline1dvnum, "Spline1DUnserialize: corrupted array!");
887            c.n = (int)Math.Round(ra[2]);
888            c.k = (int)Math.Round(ra[3]);
889            c.x = new double[c.n];
890            c.c = new double[(c.k+1)*(c.n-1)];
891            i1_ = (4) - (0);
892            for(i_=0; i_<=c.n-1;i_++)
893            {
894                c.x[i_] = ra[i_+i1_];
895            }
896            i1_ = (4+c.n) - (0);
897            for(i_=0; i_<=(c.k+1)*(c.n-1)-1;i_++)
898            {
899                c.c[i_] = ra[i_+i1_];
900            }
901        }
902
903
904        /*************************************************************************
905        This subroutine unpacks the spline into the coefficients table.
906
907        INPUT PARAMETERS:
908            C   -   spline interpolant.
909            X   -   point
910
911        Result:
912            Tbl -   coefficients table, unpacked format, array[0..N-2, 0..5].
913                    For I = 0...N-2:
914                        Tbl[I,0] = X[i]
915                        Tbl[I,1] = X[i+1]
916                        Tbl[I,2] = C0
917                        Tbl[I,3] = C1
918                        Tbl[I,4] = C2
919                        Tbl[I,5] = C3
920                    On [x[i], x[i+1]] spline is equals to:
921                        S(x) = C0 + C1*t + C2*t^2 + C3*t^3
922                        t = x-x[i]
923
924          -- ALGLIB PROJECT --
925             Copyright 29.06.2007 by Bochkanov Sergey
926        *************************************************************************/
927        public static void spline1dunpack(ref spline1dinterpolant c,
928            ref int n,
929            ref double[,] tbl)
930        {
931            int i = 0;
932            int j = 0;
933
934            tbl = new double[c.n-2+1, 2+c.k+1];
935            n = c.n;
936           
937            //
938            // Fill
939            //
940            for(i=0; i<=n-2; i++)
941            {
942                tbl[i,0] = c.x[i];
943                tbl[i,1] = c.x[i+1];
944                for(j=0; j<=c.k; j++)
945                {
946                    tbl[i,2+j] = c.c[(c.k+1)*i+j];
947                }
948            }
949        }
950
951
952        /*************************************************************************
953        This subroutine performs linear transformation of the spline argument.
954
955        INPUT PARAMETERS:
956            C   -   spline interpolant.
957            A, B-   transformation coefficients: x = A*t + B
958        Result:
959            C   -   transformed spline
960
961          -- ALGLIB PROJECT --
962             Copyright 30.06.2007 by Bochkanov Sergey
963        *************************************************************************/
964        public static void spline1dlintransx(ref spline1dinterpolant c,
965            double a,
966            double b)
967        {
968            int i = 0;
969            int j = 0;
970            int n = 0;
971            double v = 0;
972            double dv = 0;
973            double d2v = 0;
974            double[] x = new double[0];
975            double[] y = new double[0];
976            double[] d = new double[0];
977
978            n = c.n;
979           
980            //
981            // Special case: A=0
982            //
983            if( (double)(a)==(double)(0) )
984            {
985                v = spline1dcalc(ref c, b);
986                for(i=0; i<=n-2; i++)
987                {
988                    c.c[(c.k+1)*i] = v;
989                    for(j=1; j<=c.k; j++)
990                    {
991                        c.c[(c.k+1)*i+j] = 0;
992                    }
993                }
994                return;
995            }
996           
997            //
998            // General case: A<>0.
999            // Unpack, X, Y, dY/dX.
1000            // Scale and pack again.
1001            //
1002            System.Diagnostics.Debug.Assert(c.k==3, "Spline1DLinTransX: internal error");
1003            x = new double[n-1+1];
1004            y = new double[n-1+1];
1005            d = new double[n-1+1];
1006            for(i=0; i<=n-1; i++)
1007            {
1008                x[i] = c.x[i];
1009                spline1ddiff(ref c, x[i], ref v, ref dv, ref d2v);
1010                x[i] = (x[i]-b)/a;
1011                y[i] = v;
1012                d[i] = a*dv;
1013            }
1014            spline1dbuildhermite(x, y, d, n, ref c);
1015        }
1016
1017
1018        /*************************************************************************
1019        This subroutine performs linear transformation of the spline.
1020
1021        INPUT PARAMETERS:
1022            C   -   spline interpolant.
1023            A, B-   transformation coefficients: S2(x) = A*S(x) + B
1024        Result:
1025            C   -   transformed spline
1026
1027          -- ALGLIB PROJECT --
1028             Copyright 30.06.2007 by Bochkanov Sergey
1029        *************************************************************************/
1030        public static void spline1dlintransy(ref spline1dinterpolant c,
1031            double a,
1032            double b)
1033        {
1034            int i = 0;
1035            int j = 0;
1036            int n = 0;
1037
1038            n = c.n;
1039            for(i=0; i<=n-2; i++)
1040            {
1041                c.c[(c.k+1)*i] = a*c.c[(c.k+1)*i]+b;
1042                for(j=1; j<=c.k; j++)
1043                {
1044                    c.c[(c.k+1)*i+j] = a*c.c[(c.k+1)*i+j];
1045                }
1046            }
1047        }
1048
1049
1050        /*************************************************************************
1051        This subroutine integrates the spline.
1052
1053        INPUT PARAMETERS:
1054            C   -   spline interpolant.
1055            X   -   right bound of the integration interval [a, x]
1056        Result:
1057            integral(S(t)dt,a,x)
1058
1059          -- ALGLIB PROJECT --
1060             Copyright 23.06.2007 by Bochkanov Sergey
1061        *************************************************************************/
1062        public static double spline1dintegrate(ref spline1dinterpolant c,
1063            double x)
1064        {
1065            double result = 0;
1066            int n = 0;
1067            int i = 0;
1068            int j = 0;
1069            int l = 0;
1070            int r = 0;
1071            int m = 0;
1072            double w = 0;
1073            double v = 0;
1074
1075            n = c.n;
1076           
1077            //
1078            // Binary search in the [ x[0], ..., x[n-2] ] (x[n-1] is not included)
1079            //
1080            l = 0;
1081            r = n-2+1;
1082            while( l!=r-1 )
1083            {
1084                m = (l+r)/2;
1085                if( (double)(c.x[m])>=(double)(x) )
1086                {
1087                    r = m;
1088                }
1089                else
1090                {
1091                    l = m;
1092                }
1093            }
1094           
1095            //
1096            // Integration
1097            //
1098            result = 0;
1099            for(i=0; i<=l-1; i++)
1100            {
1101                w = c.x[i+1]-c.x[i];
1102                m = (c.k+1)*i;
1103                result = result+c.c[m]*w;
1104                v = w;
1105                for(j=1; j<=c.k; j++)
1106                {
1107                    v = v*w;
1108                    result = result+c.c[m+j]*v/(j+1);
1109                }
1110            }
1111            w = x-c.x[l];
1112            m = (c.k+1)*l;
1113            v = w;
1114            result = result+c.c[m]*w;
1115            for(j=1; j<=c.k; j++)
1116            {
1117                v = v*w;
1118                result = result+c.c[m+j]*v/(j+1);
1119            }
1120            return result;
1121        }
1122
1123
1124        /*************************************************************************
1125        Internal spline fitting subroutine
1126
1127          -- ALGLIB PROJECT --
1128             Copyright 08.09.2009 by Bochkanov Sergey
1129        *************************************************************************/
1130        private static void spline1dfitinternal(int st,
1131            double[] x,
1132            double[] y,
1133            ref double[] w,
1134            int n,
1135            double[] xc,
1136            double[] yc,
1137            ref int[] dc,
1138            int k,
1139            int m,
1140            ref int info,
1141            ref spline1dinterpolant s,
1142            ref spline1dfitreport rep)
1143        {
1144            double[,] fmatrix = new double[0,0];
1145            double[,] cmatrix = new double[0,0];
1146            double[] y2 = new double[0];
1147            double[] w2 = new double[0];
1148            double[] sx = new double[0];
1149            double[] sy = new double[0];
1150            double[] sd = new double[0];
1151            double[] tmp = new double[0];
1152            double[] xoriginal = new double[0];
1153            double[] yoriginal = new double[0];
1154            lsfit.lsfitreport lrep = new lsfit.lsfitreport();
1155            double v0 = 0;
1156            double v1 = 0;
1157            double v2 = 0;
1158            double mx = 0;
1159            spline1dinterpolant s2 = new spline1dinterpolant();
1160            int i = 0;
1161            int j = 0;
1162            int relcnt = 0;
1163            double xa = 0;
1164            double xb = 0;
1165            double sa = 0;
1166            double sb = 0;
1167            double bl = 0;
1168            double br = 0;
1169            double decay = 0;
1170            int i_ = 0;
1171
1172            x = (double[])x.Clone();
1173            y = (double[])y.Clone();
1174            xc = (double[])xc.Clone();
1175            yc = (double[])yc.Clone();
1176
1177            System.Diagnostics.Debug.Assert(st==0 | st==1, "Spline1DFit: internal error!");
1178            if( st==0 & m<4 )
1179            {
1180                info = -1;
1181                return;
1182            }
1183            if( st==1 & m<4 )
1184            {
1185                info = -1;
1186                return;
1187            }
1188            if( n<1 | k<0 | k>=m )
1189            {
1190                info = -1;
1191                return;
1192            }
1193            for(i=0; i<=k-1; i++)
1194            {
1195                info = 0;
1196                if( dc[i]<0 )
1197                {
1198                    info = -1;
1199                }
1200                if( dc[i]>1 )
1201                {
1202                    info = -1;
1203                }
1204                if( info<0 )
1205                {
1206                    return;
1207                }
1208            }
1209            if( st==1 & m%2!=0 )
1210            {
1211               
1212                //
1213                // Hermite fitter must have even number of basis functions
1214                //
1215                info = -2;
1216                return;
1217            }
1218           
1219            //
1220            // weight decay for correct handling of task which becomes
1221            // degenerate after constraints are applied
1222            //
1223            decay = 10000*AP.Math.MachineEpsilon;
1224           
1225            //
1226            // Scale X, Y, XC, YC
1227            //
1228            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);
1229           
1230            //
1231            // allocate space, initialize:
1232            // * SX     -   grid for basis functions
1233            // * SY     -   values of basis functions at grid points
1234            // * FMatrix-   values of basis functions at X[]
1235            // * CMatrix-   values (derivatives) of basis functions at XC[]
1236            //
1237            y2 = new double[n+m];
1238            w2 = new double[n+m];
1239            fmatrix = new double[n+m, m];
1240            if( k>0 )
1241            {
1242                cmatrix = new double[k, m+1];
1243            }
1244            if( st==0 )
1245            {
1246               
1247                //
1248                // allocate space for cubic spline
1249                //
1250                sx = new double[m-2];
1251                sy = new double[m-2];
1252                for(j=0; j<=m-2-1; j++)
1253                {
1254                    sx[j] = (double)(2*j)/((double)(m-2-1))-1;
1255                }
1256            }
1257            if( st==1 )
1258            {
1259               
1260                //
1261                // allocate space for Hermite spline
1262                //
1263                sx = new double[m/2];
1264                sy = new double[m/2];
1265                sd = new double[m/2];
1266                for(j=0; j<=m/2-1; j++)
1267                {
1268                    sx[j] = (double)(2*j)/((double)(m/2-1))-1;
1269                }
1270            }
1271           
1272            //
1273            // Prepare design and constraints matrices:
1274            // * fill constraints matrix
1275            // * fill first N rows of design matrix with values
1276            // * fill next M rows of design matrix with regularizing term
1277            // * append M zeros to Y
1278            // * append M elements, mean(abs(W)) each, to W
1279            //
1280            for(j=0; j<=m-1; j++)
1281            {
1282               
1283                //
1284                // prepare Jth basis function
1285                //
1286                if( st==0 )
1287                {
1288                   
1289                    //
1290                    // cubic spline basis
1291                    //
1292                    for(i=0; i<=m-2-1; i++)
1293                    {
1294                        sy[i] = 0;
1295                    }
1296                    bl = 0;
1297                    br = 0;
1298                    if( j<m-2 )
1299                    {
1300                        sy[j] = 1;
1301                    }
1302                    if( j==m-2 )
1303                    {
1304                        bl = 1;
1305                    }
1306                    if( j==m-1 )
1307                    {
1308                        br = 1;
1309                    }
1310                    spline1dbuildcubic(sx, sy, m-2, 1, bl, 1, br, ref s2);
1311                }
1312                if( st==1 )
1313                {
1314                   
1315                    //
1316                    // Hermite basis
1317                    //
1318                    for(i=0; i<=m/2-1; i++)
1319                    {
1320                        sy[i] = 0;
1321                        sd[i] = 0;
1322                    }
1323                    if( j%2==0 )
1324                    {
1325                        sy[j/2] = 1;
1326                    }
1327                    else
1328                    {
1329                        sd[j/2] = 1;
1330                    }
1331                    spline1dbuildhermite(sx, sy, sd, m/2, ref s2);
1332                }
1333               
1334                //
1335                // values at X[], XC[]
1336                //
1337                for(i=0; i<=n-1; i++)
1338                {
1339                    fmatrix[i,j] = spline1dcalc(ref s2, x[i]);
1340                }
1341                for(i=0; i<=k-1; i++)
1342                {
1343                    System.Diagnostics.Debug.Assert(dc[i]>=0 & dc[i]<=2, "Spline1DFit: internal error!");
1344                    spline1ddiff(ref s2, xc[i], ref v0, ref v1, ref v2);
1345                    if( dc[i]==0 )
1346                    {
1347                        cmatrix[i,j] = v0;
1348                    }
1349                    if( dc[i]==1 )
1350                    {
1351                        cmatrix[i,j] = v1;
1352                    }
1353                    if( dc[i]==2 )
1354                    {
1355                        cmatrix[i,j] = v2;
1356                    }
1357                }
1358            }
1359            for(i=0; i<=k-1; i++)
1360            {
1361                cmatrix[i,m] = yc[i];
1362            }
1363            for(i=0; i<=m-1; i++)
1364            {
1365                for(j=0; j<=m-1; j++)
1366                {
1367                    if( i==j )
1368                    {
1369                        fmatrix[n+i,j] = decay;
1370                    }
1371                    else
1372                    {
1373                        fmatrix[n+i,j] = 0;
1374                    }
1375                }
1376            }
1377            y2 = new double[n+m];
1378            w2 = new double[n+m];
1379            for(i_=0; i_<=n-1;i_++)
1380            {
1381                y2[i_] = y[i_];
1382            }
1383            for(i_=0; i_<=n-1;i_++)
1384            {
1385                w2[i_] = w[i_];
1386            }
1387            mx = 0;
1388            for(i=0; i<=n-1; i++)
1389            {
1390                mx = mx+Math.Abs(w[i]);
1391            }
1392            mx = mx/n;
1393            for(i=0; i<=m-1; i++)
1394            {
1395                y2[n+i] = 0;
1396                w2[n+i] = mx;
1397            }
1398           
1399            //
1400            // Solve constrained task
1401            //
1402            if( k>0 )
1403            {
1404               
1405                //
1406                // solve using regularization
1407                //
1408                lsfit.lsfitlinearwc(y2, ref w2, ref fmatrix, cmatrix, n+m, m, k, ref info, ref tmp, ref lrep);
1409            }
1410            else
1411            {
1412               
1413                //
1414                // no constraints, no regularization needed
1415                //
1416                lsfit.lsfitlinearwc(y, ref w, ref fmatrix, cmatrix, n, m, k, ref info, ref tmp, ref lrep);
1417            }
1418            if( info<0 )
1419            {
1420                return;
1421            }
1422           
1423            //
1424            // Generate spline and scale it
1425            //
1426            if( st==0 )
1427            {
1428               
1429                //
1430                // cubic spline basis
1431                //
1432                for(i_=0; i_<=m-2-1;i_++)
1433                {
1434                    sy[i_] = tmp[i_];
1435                }
1436                spline1dbuildcubic(sx, sy, m-2, 1, tmp[m-2], 1, tmp[m-1], ref s);
1437            }
1438            if( st==1 )
1439            {
1440               
1441                //
1442                // Hermite basis
1443                //
1444                for(i=0; i<=m/2-1; i++)
1445                {
1446                    sy[i] = tmp[2*i];
1447                    sd[i] = tmp[2*i+1];
1448                }
1449                spline1dbuildhermite(sx, sy, sd, m/2, ref s);
1450            }
1451            spline1dlintransx(ref s, 2/(xb-xa), -((xa+xb)/(xb-xa)));
1452            spline1dlintransy(ref s, sb-sa, sa);
1453           
1454            //
1455            // Scale absolute errors obtained from LSFitLinearW.
1456            // Relative error should be calculated separately
1457            // (because of shifting/scaling of the task)
1458            //
1459            rep.taskrcond = lrep.taskrcond;
1460            rep.rmserror = lrep.rmserror*(sb-sa);
1461            rep.avgerror = lrep.avgerror*(sb-sa);
1462            rep.maxerror = lrep.maxerror*(sb-sa);
1463            rep.avgrelerror = 0;
1464            relcnt = 0;
1465            for(i=0; i<=n-1; i++)
1466            {
1467                if( (double)(yoriginal[i])!=(double)(0) )
1468                {
1469                    rep.avgrelerror = rep.avgrelerror+Math.Abs(spline1dcalc(ref s, xoriginal[i])-yoriginal[i])/Math.Abs(yoriginal[i]);
1470                    relcnt = relcnt+1;
1471                }
1472            }
1473            if( relcnt!=0 )
1474            {
1475                rep.avgrelerror = rep.avgrelerror/relcnt;
1476            }
1477        }
1478
1479
1480        /*************************************************************************
1481        Internal subroutine. Heap sort.
1482        *************************************************************************/
1483        private static void heapsortpoints(ref double[] x,
1484            ref double[] y,
1485            int n)
1486        {
1487            int i = 0;
1488            int j = 0;
1489            int k = 0;
1490            int t = 0;
1491            double tmp = 0;
1492            bool isascending = new bool();
1493            bool isdescending = new bool();
1494
1495           
1496            //
1497            // Test for already sorted set
1498            //
1499            isascending = true;
1500            isdescending = true;
1501            for(i=1; i<=n-1; i++)
1502            {
1503                isascending = isascending & (double)(x[i])>(double)(x[i-1]);
1504                isdescending = isdescending & (double)(x[i])<(double)(x[i-1]);
1505            }
1506            if( isascending )
1507            {
1508                return;
1509            }
1510            if( isdescending )
1511            {
1512                for(i=0; i<=n-1; i++)
1513                {
1514                    j = n-1-i;
1515                    if( j<=i )
1516                    {
1517                        break;
1518                    }
1519                    tmp = x[i];
1520                    x[i] = x[j];
1521                    x[j] = tmp;
1522                    tmp = y[i];
1523                    y[i] = y[j];
1524                    y[j] = tmp;
1525                }
1526                return;
1527            }
1528           
1529            //
1530            // Special case: N=1
1531            //
1532            if( n==1 )
1533            {
1534                return;
1535            }
1536           
1537            //
1538            // General case
1539            //
1540            i = 2;
1541            do
1542            {
1543                t = i;
1544                while( t!=1 )
1545                {
1546                    k = t/2;
1547                    if( (double)(x[k-1])>=(double)(x[t-1]) )
1548                    {
1549                        t = 1;
1550                    }
1551                    else
1552                    {
1553                        tmp = x[k-1];
1554                        x[k-1] = x[t-1];
1555                        x[t-1] = tmp;
1556                        tmp = y[k-1];
1557                        y[k-1] = y[t-1];
1558                        y[t-1] = tmp;
1559                        t = k;
1560                    }
1561                }
1562                i = i+1;
1563            }
1564            while( i<=n );
1565            i = n-1;
1566            do
1567            {
1568                tmp = x[i];
1569                x[i] = x[0];
1570                x[0] = tmp;
1571                tmp = y[i];
1572                y[i] = y[0];
1573                y[0] = tmp;
1574                t = 1;
1575                while( t!=0 )
1576                {
1577                    k = 2*t;
1578                    if( k>i )
1579                    {
1580                        t = 0;
1581                    }
1582                    else
1583                    {
1584                        if( k<i )
1585                        {
1586                            if( (double)(x[k])>(double)(x[k-1]) )
1587                            {
1588                                k = k+1;
1589                            }
1590                        }
1591                        if( (double)(x[t-1])>=(double)(x[k-1]) )
1592                        {
1593                            t = 0;
1594                        }
1595                        else
1596                        {
1597                            tmp = x[k-1];
1598                            x[k-1] = x[t-1];
1599                            x[t-1] = tmp;
1600                            tmp = y[k-1];
1601                            y[k-1] = y[t-1];
1602                            y[t-1] = tmp;
1603                            t = k;
1604                        }
1605                    }
1606                }
1607                i = i-1;
1608            }
1609            while( i>=1 );
1610        }
1611
1612
1613        /*************************************************************************
1614        Internal subroutine. Heap sort.
1615        *************************************************************************/
1616        private static void heapsortdpoints(ref double[] x,
1617            ref double[] y,
1618            ref double[] d,
1619            int n)
1620        {
1621            int i = 0;
1622            int j = 0;
1623            int k = 0;
1624            int t = 0;
1625            double tmp = 0;
1626            bool isascending = new bool();
1627            bool isdescending = new bool();
1628
1629           
1630            //
1631            // Test for already sorted set
1632            //
1633            isascending = true;
1634            isdescending = true;
1635            for(i=1; i<=n-1; i++)
1636            {
1637                isascending = isascending & (double)(x[i])>(double)(x[i-1]);
1638                isdescending = isdescending & (double)(x[i])<(double)(x[i-1]);
1639            }
1640            if( isascending )
1641            {
1642                return;
1643            }
1644            if( isdescending )
1645            {
1646                for(i=0; i<=n-1; i++)
1647                {
1648                    j = n-1-i;
1649                    if( j<=i )
1650                    {
1651                        break;
1652                    }
1653                    tmp = x[i];
1654                    x[i] = x[j];
1655                    x[j] = tmp;
1656                    tmp = y[i];
1657                    y[i] = y[j];
1658                    y[j] = tmp;
1659                    tmp = d[i];
1660                    d[i] = d[j];
1661                    d[j] = tmp;
1662                }
1663                return;
1664            }
1665           
1666            //
1667            // Special case: N=1
1668            //
1669            if( n==1 )
1670            {
1671                return;
1672            }
1673           
1674            //
1675            // General case
1676            //
1677            i = 2;
1678            do
1679            {
1680                t = i;
1681                while( t!=1 )
1682                {
1683                    k = t/2;
1684                    if( (double)(x[k-1])>=(double)(x[t-1]) )
1685                    {
1686                        t = 1;
1687                    }
1688                    else
1689                    {
1690                        tmp = x[k-1];
1691                        x[k-1] = x[t-1];
1692                        x[t-1] = tmp;
1693                        tmp = y[k-1];
1694                        y[k-1] = y[t-1];
1695                        y[t-1] = tmp;
1696                        tmp = d[k-1];
1697                        d[k-1] = d[t-1];
1698                        d[t-1] = tmp;
1699                        t = k;
1700                    }
1701                }
1702                i = i+1;
1703            }
1704            while( i<=n );
1705            i = n-1;
1706            do
1707            {
1708                tmp = x[i];
1709                x[i] = x[0];
1710                x[0] = tmp;
1711                tmp = y[i];
1712                y[i] = y[0];
1713                y[0] = tmp;
1714                tmp = d[i];
1715                d[i] = d[0];
1716                d[0] = tmp;
1717                t = 1;
1718                while( t!=0 )
1719                {
1720                    k = 2*t;
1721                    if( k>i )
1722                    {
1723                        t = 0;
1724                    }
1725                    else
1726                    {
1727                        if( k<i )
1728                        {
1729                            if( (double)(x[k])>(double)(x[k-1]) )
1730                            {
1731                                k = k+1;
1732                            }
1733                        }
1734                        if( (double)(x[t-1])>=(double)(x[k-1]) )
1735                        {
1736                            t = 0;
1737                        }
1738                        else
1739                        {
1740                            tmp = x[k-1];
1741                            x[k-1] = x[t-1];
1742                            x[t-1] = tmp;
1743                            tmp = y[k-1];
1744                            y[k-1] = y[t-1];
1745                            y[t-1] = tmp;
1746                            tmp = d[k-1];
1747                            d[k-1] = d[t-1];
1748                            d[t-1] = tmp;
1749                            t = k;
1750                        }
1751                    }
1752                }
1753                i = i-1;
1754            }
1755            while( i>=1 );
1756        }
1757
1758
1759        /*************************************************************************
1760        Internal subroutine. Tridiagonal solver.
1761        *************************************************************************/
1762        private static void solvetridiagonal(double[] a,
1763            double[] b,
1764            double[] c,
1765            double[] d,
1766            int n,
1767            ref double[] x)
1768        {
1769            int k = 0;
1770            double t = 0;
1771
1772            a = (double[])a.Clone();
1773            b = (double[])b.Clone();
1774            c = (double[])c.Clone();
1775            d = (double[])d.Clone();
1776
1777            x = new double[n-1+1];
1778            a[0] = 0;
1779            c[n-1] = 0;
1780            for(k=1; k<=n-1; k++)
1781            {
1782                t = a[k]/b[k-1];
1783                b[k] = b[k]-t*c[k-1];
1784                d[k] = d[k]-t*d[k-1];
1785            }
1786            x[n-1] = d[n-1]/b[n-1];
1787            for(k=n-2; k>=0; k--)
1788            {
1789                x[k] = (d[k]-c[k]*x[k+1])/b[k];
1790            }
1791        }
1792
1793
1794        /*************************************************************************
1795        Internal subroutine. Three-point differentiation
1796        *************************************************************************/
1797        private static double diffthreepoint(double t,
1798            double x0,
1799            double f0,
1800            double x1,
1801            double f1,
1802            double x2,
1803            double f2)
1804        {
1805            double result = 0;
1806            double a = 0;
1807            double b = 0;
1808
1809            t = t-x0;
1810            x1 = x1-x0;
1811            x2 = x2-x0;
1812            a = (f2-f0-x2/x1*(f1-f0))/(AP.Math.Sqr(x2)-x1*x2);
1813            b = (f1-f0-a*AP.Math.Sqr(x1))/x1;
1814            result = 2*a*t+b;
1815            return result;
1816        }
1817    }
1818}
Note: See TracBrowser for help on using the repository browser.