Free cookie consent management tool by TermsFeed Policy Generator

source: tags/3.3.2/HeuristicLab.ExtLibs/HeuristicLab.ALGLIB/2.5.0/ALGLIB-2.5.0/lsfit.cs @ 13398

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

implemented first version of LR (ticket #1012)

File size: 54.9 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 lsfit
26    {
27        /*************************************************************************
28        Least squares 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 lsfitreport
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        public struct lsfitstate
46        {
47            public int n;
48            public int m;
49            public int k;
50            public double epsf;
51            public double epsx;
52            public int maxits;
53            public double stpmax;
54            public double[,] taskx;
55            public double[] tasky;
56            public double[] w;
57            public bool cheapfg;
58            public bool havehess;
59            public bool needf;
60            public bool needfg;
61            public bool needfgh;
62            public int pointindex;
63            public double[] x;
64            public double[] c;
65            public double f;
66            public double[] g;
67            public double[,] h;
68            public int repterminationtype;
69            public double reprmserror;
70            public double repavgerror;
71            public double repavgrelerror;
72            public double repmaxerror;
73            public minlm.minlmstate optstate;
74            public minlm.minlmreport optrep;
75            public AP.rcommstate rstate;
76        };
77
78
79
80
81        /*************************************************************************
82        Weighted linear least squares fitting.
83
84        QR decomposition is used to reduce task to MxM, then triangular solver  or
85        SVD-based solver is used depending on condition number of the  system.  It
86        allows to maximize speed and retain decent accuracy.
87
88        INPUT PARAMETERS:
89            Y       -   array[0..N-1] Function values in  N  points.
90            W       -   array[0..N-1]  Weights  corresponding to function  values.
91                        Each summand in square  sum  of  approximation  deviations
92                        from  given  values  is  multiplied  by  the   square   of
93                        corresponding weight.
94            FMatrix -   a table of basis functions values, array[0..N-1, 0..M-1].
95                        FMatrix[I, J] - value of J-th basis function in I-th point.
96            N       -   number of points used. N>=1.
97            M       -   number of basis functions, M>=1.
98
99        OUTPUT PARAMETERS:
100            Info    -   error code:
101                        * -4    internal SVD decomposition subroutine failed (very
102                                rare and for degenerate systems only)
103                        * -1    incorrect N/M were specified
104                        *  1    task is solved
105            C       -   decomposition coefficients, array[0..M-1]
106            Rep     -   fitting report. Following fields are set:
107                        * Rep.TaskRCond     reciprocal of condition number
108                        * RMSError          rms error on the (X,Y).
109                        * AvgError          average error on the (X,Y).
110                        * AvgRelError       average relative error on the non-zero Y
111                        * MaxError          maximum error
112                                            NON-WEIGHTED ERRORS ARE CALCULATED
113
114        SEE ALSO
115            LSFitLinear
116            LSFitLinearC
117            LSFitLinearWC
118
119          -- ALGLIB --
120             Copyright 17.08.2009 by Bochkanov Sergey
121        *************************************************************************/
122        public static void lsfitlinearw(ref double[] y,
123            ref double[] w,
124            ref double[,] fmatrix,
125            int n,
126            int m,
127            ref int info,
128            ref double[] c,
129            ref lsfitreport rep)
130        {
131            lsfitlinearinternal(ref y, ref w, ref fmatrix, n, m, ref info, ref c, ref rep);
132        }
133
134
135        /*************************************************************************
136        Weighted constained linear least squares fitting.
137
138        This  is  variation  of LSFitLinearW(), which searchs for min|A*x=b| given
139        that  K  additional  constaints  C*x=bc are satisfied. It reduces original
140        task to modified one: min|B*y-d| WITHOUT constraints,  then LSFitLinearW()
141        is called.
142
143        INPUT PARAMETERS:
144            Y       -   array[0..N-1] Function values in  N  points.
145            W       -   array[0..N-1]  Weights  corresponding to function  values.
146                        Each summand in square  sum  of  approximation  deviations
147                        from  given  values  is  multiplied  by  the   square   of
148                        corresponding weight.
149            FMatrix -   a table of basis functions values, array[0..N-1, 0..M-1].
150                        FMatrix[I,J] - value of J-th basis function in I-th point.
151            CMatrix -   a table of constaints, array[0..K-1,0..M].
152                        I-th row of CMatrix corresponds to I-th linear constraint:
153                        CMatrix[I,0]*C[0] + ... + CMatrix[I,M-1]*C[M-1] = CMatrix[I,M]
154            N       -   number of points used. N>=1.
155            M       -   number of basis functions, M>=1.
156            K       -   number of constraints, 0 <= K < M
157                        K=0 corresponds to absence of constraints.
158
159        OUTPUT PARAMETERS:
160            Info    -   error code:
161                        * -4    internal SVD decomposition subroutine failed (very
162                                rare and for degenerate systems only)
163                        * -3    either   too   many  constraints  (M   or   more),
164                                degenerate  constraints   (some   constraints  are
165                                repetead twice) or inconsistent  constraints  were
166                                specified.
167                        * -1    incorrect N/M/K were specified
168                        *  1    task is solved
169            C       -   decomposition coefficients, array[0..M-1]
170            Rep     -   fitting report. Following fields are set:
171                        * RMSError          rms error on the (X,Y).
172                        * AvgError          average error on the (X,Y).
173                        * AvgRelError       average relative error on the non-zero Y
174                        * MaxError          maximum error
175                                            NON-WEIGHTED ERRORS ARE CALCULATED
176
177        IMPORTANT:
178            this subroitine doesn't calculate task's condition number for K<>0.
179
180        SEE ALSO
181            LSFitLinear
182            LSFitLinearC
183            LSFitLinearWC
184
185          -- ALGLIB --
186             Copyright 07.09.2009 by Bochkanov Sergey
187        *************************************************************************/
188        public static void lsfitlinearwc(double[] y,
189            ref double[] w,
190            ref double[,] fmatrix,
191            double[,] cmatrix,
192            int n,
193            int m,
194            int k,
195            ref int info,
196            ref double[] c,
197            ref lsfitreport rep)
198        {
199            int i = 0;
200            int j = 0;
201            double[] tau = new double[0];
202            double[,] q = new double[0,0];
203            double[,] f2 = new double[0,0];
204            double[] tmp = new double[0];
205            double[] c0 = new double[0];
206            double v = 0;
207            int i_ = 0;
208
209            y = (double[])y.Clone();
210            cmatrix = (double[,])cmatrix.Clone();
211
212            if( n<1 | m<1 | k<0 )
213            {
214                info = -1;
215                return;
216            }
217            if( k>=m )
218            {
219                info = -3;
220                return;
221            }
222           
223            //
224            // Solve
225            //
226            if( k==0 )
227            {
228               
229                //
230                // no constraints
231                //
232                lsfitlinearinternal(ref y, ref w, ref fmatrix, n, m, ref info, ref c, ref rep);
233            }
234            else
235            {
236               
237                //
238                // First, find general form solution of constraints system:
239                // * factorize C = L*Q
240                // * unpack Q
241                // * fill upper part of C with zeros (for RCond)
242                //
243                // We got C=C0+Q2'*y where Q2 is lower M-K rows of Q.
244                //
245                ortfac.rmatrixlq(ref cmatrix, k, m, ref tau);
246                ortfac.rmatrixlqunpackq(ref cmatrix, k, m, ref tau, m, ref q);
247                for(i=0; i<=k-1; i++)
248                {
249                    for(j=i+1; j<=m-1; j++)
250                    {
251                        cmatrix[i,j] = 0.0;
252                    }
253                }
254                if( (double)(rcond.rmatrixlurcondinf(ref cmatrix, k))<(double)(1000*AP.Math.MachineEpsilon) )
255                {
256                    info = -3;
257                    return;
258                }
259                tmp = new double[k];
260                for(i=0; i<=k-1; i++)
261                {
262                    if( i>0 )
263                    {
264                        v = 0.0;
265                        for(i_=0; i_<=i-1;i_++)
266                        {
267                            v += cmatrix[i,i_]*tmp[i_];
268                        }
269                    }
270                    else
271                    {
272                        v = 0;
273                    }
274                    tmp[i] = (cmatrix[i,m]-v)/cmatrix[i,i];
275                }
276                c0 = new double[m];
277                for(i=0; i<=m-1; i++)
278                {
279                    c0[i] = 0;
280                }
281                for(i=0; i<=k-1; i++)
282                {
283                    v = tmp[i];
284                    for(i_=0; i_<=m-1;i_++)
285                    {
286                        c0[i_] = c0[i_] + v*q[i,i_];
287                    }
288                }
289               
290                //
291                // Second, prepare modified matrix F2 = F*Q2' and solve modified task
292                //
293                tmp = new double[Math.Max(n, m)+1];
294                f2 = new double[n, m-k];
295                blas.matrixvectormultiply(ref fmatrix, 0, n-1, 0, m-1, false, ref c0, 0, m-1, -1.0, ref y, 0, n-1, 1.0);
296                blas.matrixmatrixmultiply(ref fmatrix, 0, n-1, 0, m-1, false, ref q, k, m-1, 0, m-1, true, 1.0, ref f2, 0, n-1, 0, m-k-1, 0.0, ref tmp);
297                lsfitlinearinternal(ref y, ref w, ref f2, n, m-k, ref info, ref tmp, ref rep);
298                rep.taskrcond = -1;
299                if( info<=0 )
300                {
301                    return;
302                }
303               
304                //
305                // then, convert back to original answer: C = C0 + Q2'*Y0
306                //
307                c = new double[m];
308                for(i_=0; i_<=m-1;i_++)
309                {
310                    c[i_] = c0[i_];
311                }
312                blas.matrixvectormultiply(ref q, k, m-1, 0, m-1, true, ref tmp, 0, m-k-1, 1.0, ref c, 0, m-1, 1.0);
313            }
314        }
315
316
317        /*************************************************************************
318        Linear least squares fitting, without weights.
319
320        See LSFitLinearW for more information.
321
322          -- ALGLIB --
323             Copyright 17.08.2009 by Bochkanov Sergey
324        *************************************************************************/
325        public static void lsfitlinear(ref double[] y,
326            ref double[,] fmatrix,
327            int n,
328            int m,
329            ref int info,
330            ref double[] c,
331            ref lsfitreport rep)
332        {
333            double[] w = new double[0];
334            int i = 0;
335
336            if( n<1 )
337            {
338                info = -1;
339                return;
340            }
341            w = new double[n];
342            for(i=0; i<=n-1; i++)
343            {
344                w[i] = 1;
345            }
346            lsfitlinearinternal(ref y, ref w, ref fmatrix, n, m, ref info, ref c, ref rep);
347        }
348
349
350        /*************************************************************************
351        Constained linear least squares fitting, without weights.
352
353        See LSFitLinearWC() for more information.
354
355          -- ALGLIB --
356             Copyright 07.09.2009 by Bochkanov Sergey
357        *************************************************************************/
358        public static void lsfitlinearc(double[] y,
359            ref double[,] fmatrix,
360            ref double[,] cmatrix,
361            int n,
362            int m,
363            int k,
364            ref int info,
365            ref double[] c,
366            ref lsfitreport rep)
367        {
368            double[] w = new double[0];
369            int i = 0;
370
371            y = (double[])y.Clone();
372
373            if( n<1 )
374            {
375                info = -1;
376                return;
377            }
378            w = new double[n];
379            for(i=0; i<=n-1; i++)
380            {
381                w[i] = 1;
382            }
383            lsfitlinearwc(y, ref w, ref fmatrix, cmatrix, n, m, k, ref info, ref c, ref rep);
384        }
385
386
387        /*************************************************************************
388        Weighted nonlinear least squares fitting using gradient and Hessian.
389
390        Nonlinear task min(F(c)) is solved, where
391
392            F(c) = (w[0]*(f(x[0],c)-y[0]))^2 + ... + (w[n-1]*(f(x[n-1],c)-y[n-1]))^2,
393           
394            * N is a number of points,
395            * M is a dimension of a space points belong to,
396            * K is a dimension of a space of parameters being fitted,
397            * w is an N-dimensional vector of weight coefficients,
398            * x is a set of N points, each of them is an M-dimensional vector,
399            * c is a K-dimensional vector of parameters being fitted
400           
401        This subroutine uses only f(x[i],c) and its gradient.
402           
403        INPUT PARAMETERS:
404            X       -   array[0..N-1,0..M-1], points (one row = one point)
405            Y       -   array[0..N-1], function values.
406            W       -   weights, array[0..N-1]
407            C       -   array[0..K-1], initial approximation to the solution,
408            N       -   number of points, N>1
409            M       -   dimension of space
410            K       -   number of parameters being fitted
411            CheapFG -   boolean flag, which is:
412                        * True  if both function and gradient calculation complexity
413                                are less than O(M^2).  An improved  algorithm  can
414                                be  used  which corresponds  to  FGJ  scheme  from
415                                MINLM unit.
416                        * False otherwise.
417                                Standard Jacibian-bases  Levenberg-Marquardt  algo
418                                will be used (FJ scheme).
419
420        OUTPUT PARAMETERS:
421            State   -   structure which stores algorithm state between subsequent
422                        calls  of   LSFitNonlinearIteration.   Used  for  reverse
423                        communication.  This  structure   should   be  passed  to
424                        LSFitNonlinearIteration subroutine.
425
426        See also:
427            LSFitNonlinearIteration
428            LSFitNonlinearResults
429            LSFitNonlinearFG (fitting without weights)
430            LSFitNonlinearWFGH (fitting using Hessian)
431            LSFitNonlinearFGH (fitting using Hessian, without weights)
432
433
434          -- ALGLIB --
435             Copyright 17.08.2009 by Bochkanov Sergey
436        *************************************************************************/
437        public static void lsfitnonlinearwfg(ref double[,] x,
438            ref double[] y,
439            ref double[] w,
440            ref double[] c,
441            int n,
442            int m,
443            int k,
444            bool cheapfg,
445            ref lsfitstate state)
446        {
447            int i = 0;
448            int i_ = 0;
449
450            state.n = n;
451            state.m = m;
452            state.k = k;
453            lsfitnonlinearsetcond(ref state, 0.0, 0.0, 0);
454            lsfitnonlinearsetstpmax(ref state, 0.0);
455            state.cheapfg = cheapfg;
456            state.havehess = false;
457            if( n>=1 & m>=1 & k>=1 )
458            {
459                state.taskx = new double[n, m];
460                state.tasky = new double[n];
461                state.w = new double[n];
462                state.c = new double[k];
463                for(i_=0; i_<=k-1;i_++)
464                {
465                    state.c[i_] = c[i_];
466                }
467                for(i_=0; i_<=n-1;i_++)
468                {
469                    state.w[i_] = w[i_];
470                }
471                for(i=0; i<=n-1; i++)
472                {
473                    for(i_=0; i_<=m-1;i_++)
474                    {
475                        state.taskx[i,i_] = x[i,i_];
476                    }
477                    state.tasky[i] = y[i];
478                }
479            }
480            state.rstate.ia = new int[4+1];
481            state.rstate.ra = new double[1+1];
482            state.rstate.stage = -1;
483        }
484
485
486        /*************************************************************************
487        Nonlinear least squares fitting, no individual weights.
488        See LSFitNonlinearWFG for more information.
489
490          -- ALGLIB --
491             Copyright 17.08.2009 by Bochkanov Sergey
492        *************************************************************************/
493        public static void lsfitnonlinearfg(ref double[,] x,
494            ref double[] y,
495            ref double[] c,
496            int n,
497            int m,
498            int k,
499            bool cheapfg,
500            ref lsfitstate state)
501        {
502            int i = 0;
503            int i_ = 0;
504
505            state.n = n;
506            state.m = m;
507            state.k = k;
508            lsfitnonlinearsetcond(ref state, 0.0, 0.0, 0);
509            lsfitnonlinearsetstpmax(ref state, 0.0);
510            state.cheapfg = cheapfg;
511            state.havehess = false;
512            if( n>=1 & m>=1 & k>=1 )
513            {
514                state.taskx = new double[n, m];
515                state.tasky = new double[n];
516                state.w = new double[n];
517                state.c = new double[k];
518                for(i_=0; i_<=k-1;i_++)
519                {
520                    state.c[i_] = c[i_];
521                }
522                for(i=0; i<=n-1; i++)
523                {
524                    for(i_=0; i_<=m-1;i_++)
525                    {
526                        state.taskx[i,i_] = x[i,i_];
527                    }
528                    state.tasky[i] = y[i];
529                    state.w[i] = 1;
530                }
531            }
532            state.rstate.ia = new int[4+1];
533            state.rstate.ra = new double[1+1];
534            state.rstate.stage = -1;
535        }
536
537
538        /*************************************************************************
539        Weighted nonlinear least squares fitting using gradient/Hessian.
540
541        Nonlinear task min(F(c)) is solved, where
542
543            F(c) = (w[0]*(f(x[0],c)-y[0]))^2 + ... + (w[n-1]*(f(x[n-1],c)-y[n-1]))^2,
544
545            * N is a number of points,
546            * M is a dimension of a space points belong to,
547            * K is a dimension of a space of parameters being fitted,
548            * w is an N-dimensional vector of weight coefficients,
549            * x is a set of N points, each of them is an M-dimensional vector,
550            * c is a K-dimensional vector of parameters being fitted
551
552        This subroutine uses f(x[i],c), its gradient and its Hessian.
553
554        See LSFitNonlinearWFG() subroutine for information about function
555        parameters.
556
557          -- ALGLIB --
558             Copyright 17.08.2009 by Bochkanov Sergey
559        *************************************************************************/
560        public static void lsfitnonlinearwfgh(ref double[,] x,
561            ref double[] y,
562            ref double[] w,
563            ref double[] c,
564            int n,
565            int m,
566            int k,
567            ref lsfitstate state)
568        {
569            int i = 0;
570            int i_ = 0;
571
572            state.n = n;
573            state.m = m;
574            state.k = k;
575            lsfitnonlinearsetcond(ref state, 0.0, 0.0, 0);
576            lsfitnonlinearsetstpmax(ref state, 0.0);
577            state.cheapfg = true;
578            state.havehess = true;
579            if( n>=1 & m>=1 & k>=1 )
580            {
581                state.taskx = new double[n, m];
582                state.tasky = new double[n];
583                state.w = new double[n];
584                state.c = new double[k];
585                for(i_=0; i_<=k-1;i_++)
586                {
587                    state.c[i_] = c[i_];
588                }
589                for(i_=0; i_<=n-1;i_++)
590                {
591                    state.w[i_] = w[i_];
592                }
593                for(i=0; i<=n-1; i++)
594                {
595                    for(i_=0; i_<=m-1;i_++)
596                    {
597                        state.taskx[i,i_] = x[i,i_];
598                    }
599                    state.tasky[i] = y[i];
600                }
601            }
602            state.rstate.ia = new int[4+1];
603            state.rstate.ra = new double[1+1];
604            state.rstate.stage = -1;
605        }
606
607
608        /*************************************************************************
609        Nonlinear least squares fitting using gradient/Hessian without  individual
610        weights. See LSFitNonlinearWFGH() for more information.
611
612
613          -- ALGLIB --
614             Copyright 17.08.2009 by Bochkanov Sergey
615        *************************************************************************/
616        public static void lsfitnonlinearfgh(ref double[,] x,
617            ref double[] y,
618            ref double[] c,
619            int n,
620            int m,
621            int k,
622            ref lsfitstate state)
623        {
624            int i = 0;
625            int i_ = 0;
626
627            state.n = n;
628            state.m = m;
629            state.k = k;
630            lsfitnonlinearsetcond(ref state, 0.0, 0.0, 0);
631            lsfitnonlinearsetstpmax(ref state, 0.0);
632            state.cheapfg = true;
633            state.havehess = true;
634            if( n>=1 & m>=1 & k>=1 )
635            {
636                state.taskx = new double[n, m];
637                state.tasky = new double[n];
638                state.w = new double[n];
639                state.c = new double[k];
640                for(i_=0; i_<=k-1;i_++)
641                {
642                    state.c[i_] = c[i_];
643                }
644                for(i=0; i<=n-1; i++)
645                {
646                    for(i_=0; i_<=m-1;i_++)
647                    {
648                        state.taskx[i,i_] = x[i,i_];
649                    }
650                    state.tasky[i] = y[i];
651                    state.w[i] = 1;
652                }
653            }
654            state.rstate.ia = new int[4+1];
655            state.rstate.ra = new double[1+1];
656            state.rstate.stage = -1;
657        }
658
659
660        /*************************************************************************
661        Stopping conditions for nonlinear least squares fitting.
662
663        INPUT PARAMETERS:
664            State   -   structure which stores algorithm state between calls and
665                        which is used for reverse communication. Must be initialized
666                        with LSFitNonLinearCreate???()
667            EpsF    -   stopping criterion. Algorithm stops if
668                        |F(k+1)-F(k)| <= EpsF*max{|F(k)|, |F(k+1)|, 1}
669            EpsX    -   stopping criterion. Algorithm stops if
670                        |X(k+1)-X(k)| <= EpsX*(1+|X(k)|)
671            MaxIts  -   stopping criterion. Algorithm stops after MaxIts iterations.
672                        MaxIts=0 means no stopping criterion.
673
674        NOTE
675
676        Passing EpsF=0, EpsX=0 and MaxIts=0 (simultaneously) will lead to automatic
677        stopping criterion selection (according to the scheme used by MINLM unit).
678
679
680          -- ALGLIB --
681             Copyright 17.08.2009 by Bochkanov Sergey
682        *************************************************************************/
683        public static void lsfitnonlinearsetcond(ref lsfitstate state,
684            double epsf,
685            double epsx,
686            int maxits)
687        {
688            System.Diagnostics.Debug.Assert((double)(epsf)>=(double)(0), "LSFitNonlinearSetCond: negative EpsF!");
689            System.Diagnostics.Debug.Assert((double)(epsx)>=(double)(0), "LSFitNonlinearSetCond: negative EpsX!");
690            System.Diagnostics.Debug.Assert(maxits>=0, "LSFitNonlinearSetCond: negative MaxIts!");
691            state.epsf = epsf;
692            state.epsx = epsx;
693            state.maxits = maxits;
694        }
695
696
697        /*************************************************************************
698        This function sets maximum step length
699
700        INPUT PARAMETERS:
701            State   -   structure which stores algorithm state between calls and
702                        which is used for reverse communication. Must be
703                        initialized with LSFitNonLinearCreate???()
704            StpMax  -   maximum step length, >=0. Set StpMax to 0.0,  if you don't
705                        want to limit step length.
706
707        Use this subroutine when you optimize target function which contains exp()
708        or  other  fast  growing  functions,  and optimization algorithm makes too
709        large  steps  which  leads  to overflow. This function allows us to reject
710        steps  that  are  too  large  (and  therefore  expose  us  to the possible
711        overflow) without actually calculating function value at the x+stp*d.
712
713        NOTE: non-zero StpMax leads to moderate  performance  degradation  because
714        intermediate  step  of  preconditioned L-BFGS optimization is incompatible
715        with limits on step size.
716
717          -- ALGLIB --
718             Copyright 02.04.2010 by Bochkanov Sergey
719        *************************************************************************/
720        public static void lsfitnonlinearsetstpmax(ref lsfitstate state,
721            double stpmax)
722        {
723            System.Diagnostics.Debug.Assert((double)(stpmax)>=(double)(0), "LSFitNonlinearSetStpMax: StpMax<0!");
724            state.stpmax = stpmax;
725        }
726
727
728        /*************************************************************************
729        Nonlinear least squares fitting. Algorithm iteration.
730
731        Called after inialization of the State structure with  LSFitNonlinearXXX()
732        subroutine. See HTML docs for examples.
733
734        INPUT PARAMETERS:
735            State   -   structure which stores algorithm state between  subsequent
736                        calls and which is used for reverse communication. Must be
737                        initialized with LSFitNonlinearXXX() call first.
738
739        RESULT
740        1. If subroutine returned False, iterative algorithm has converged.
741        2. If subroutine returned True, then if:
742        * if State.NeedF=True,      function value F(X,C) is required
743        * if State.NeedFG=True,     function value F(X,C) and gradient  dF/dC(X,C)
744                                    are required
745        * if State.NeedFGH=True     function value F(X,C), gradient dF/dC(X,C) and
746                                    Hessian are required
747
748        One and only one of this fields can be set at time.
749
750        Function, its gradient and Hessian are calculated at  (X,C),  where  X  is
751        stored in State.X[0..M-1] and C is stored in State.C[0..K-1].
752
753        Results are stored:
754        * function value            -   in State.F
755        * gradient                  -   in State.G[0..K-1]
756        * Hessian                   -   in State.H[0..K-1,0..K-1]
757
758          -- ALGLIB --
759             Copyright 17.08.2009 by Bochkanov Sergey
760        *************************************************************************/
761        public static bool lsfitnonlineariteration(ref lsfitstate state)
762        {
763            bool result = new bool();
764            int n = 0;
765            int m = 0;
766            int k = 0;
767            int i = 0;
768            int j = 0;
769            double v = 0;
770            double relcnt = 0;
771            int i_ = 0;
772
773           
774            //
775            // Reverse communication preparations
776            // I know it looks ugly, but it works the same way
777            // anywhere from C++ to Python.
778            //
779            // This code initializes locals by:
780            // * random values determined during code
781            //   generation - on first subroutine call
782            // * values from previous call - on subsequent calls
783            //
784            if( state.rstate.stage>=0 )
785            {
786                n = state.rstate.ia[0];
787                m = state.rstate.ia[1];
788                k = state.rstate.ia[2];
789                i = state.rstate.ia[3];
790                j = state.rstate.ia[4];
791                v = state.rstate.ra[0];
792                relcnt = state.rstate.ra[1];
793            }
794            else
795            {
796                n = -983;
797                m = -989;
798                k = -834;
799                i = 900;
800                j = -287;
801                v = 364;
802                relcnt = 214;
803            }
804            if( state.rstate.stage==0 )
805            {
806                goto lbl_0;
807            }
808            if( state.rstate.stage==1 )
809            {
810                goto lbl_1;
811            }
812            if( state.rstate.stage==2 )
813            {
814                goto lbl_2;
815            }
816            if( state.rstate.stage==3 )
817            {
818                goto lbl_3;
819            }
820            if( state.rstate.stage==4 )
821            {
822                goto lbl_4;
823            }
824           
825            //
826            // Routine body
827            //
828           
829            //
830            // check params
831            //
832            if( state.n<1 | state.m<1 | state.k<1 | (double)(state.epsf)<(double)(0) | (double)(state.epsx)<(double)(0) | state.maxits<0 )
833            {
834                state.repterminationtype = -1;
835                result = false;
836                return result;
837            }
838           
839            //
840            // init
841            //
842            n = state.n;
843            m = state.m;
844            k = state.k;
845            state.x = new double[m];
846            state.g = new double[k];
847            if( state.havehess )
848            {
849                state.h = new double[k, k];
850            }
851           
852            //
853            // initialize LM optimizer
854            //
855            if( state.havehess )
856            {
857               
858                //
859                // use Hessian.
860                // transform stopping conditions.
861                //
862                minlm.minlmcreatefgh(k, ref state.c, ref state.optstate);
863            }
864            else
865            {
866               
867                //
868                // use one of gradient-based schemes (depending on gradient cost).
869                // transform stopping conditions.
870                //
871                if( state.cheapfg )
872                {
873                    minlm.minlmcreatefgj(k, n, ref state.c, ref state.optstate);
874                }
875                else
876                {
877                    minlm.minlmcreatefj(k, n, ref state.c, ref state.optstate);
878                }
879            }
880            minlm.minlmsetcond(ref state.optstate, 0.0, state.epsf, state.epsx, state.maxits);
881            minlm.minlmsetstpmax(ref state.optstate, state.stpmax);
882           
883            //
884            // Optimize
885            //
886        lbl_5:
887            if( ! minlm.minlmiteration(ref state.optstate) )
888            {
889                goto lbl_6;
890            }
891            if( ! state.optstate.needf )
892            {
893                goto lbl_7;
894            }
895           
896            //
897            // calculate F = sum (wi*(f(xi,c)-yi))^2
898            //
899            state.optstate.f = 0;
900            i = 0;
901        lbl_9:
902            if( i>n-1 )
903            {
904                goto lbl_11;
905            }
906            for(i_=0; i_<=k-1;i_++)
907            {
908                state.c[i_] = state.optstate.x[i_];
909            }
910            for(i_=0; i_<=m-1;i_++)
911            {
912                state.x[i_] = state.taskx[i,i_];
913            }
914            state.pointindex = i;
915            lsfitclearrequestfields(ref state);
916            state.needf = true;
917            state.rstate.stage = 0;
918            goto lbl_rcomm;
919        lbl_0:
920            state.optstate.f = state.optstate.f+AP.Math.Sqr(state.w[i]*(state.f-state.tasky[i]));
921            i = i+1;
922            goto lbl_9;
923        lbl_11:
924            goto lbl_5;
925        lbl_7:
926            if( ! state.optstate.needfg )
927            {
928                goto lbl_12;
929            }
930           
931            //
932            // calculate F/gradF
933            //
934            state.optstate.f = 0;
935            for(i=0; i<=k-1; i++)
936            {
937                state.optstate.g[i] = 0;
938            }
939            i = 0;
940        lbl_14:
941            if( i>n-1 )
942            {
943                goto lbl_16;
944            }
945            for(i_=0; i_<=k-1;i_++)
946            {
947                state.c[i_] = state.optstate.x[i_];
948            }
949            for(i_=0; i_<=m-1;i_++)
950            {
951                state.x[i_] = state.taskx[i,i_];
952            }
953            state.pointindex = i;
954            lsfitclearrequestfields(ref state);
955            state.needfg = true;
956            state.rstate.stage = 1;
957            goto lbl_rcomm;
958        lbl_1:
959            state.optstate.f = state.optstate.f+AP.Math.Sqr(state.w[i]*(state.f-state.tasky[i]));
960            v = AP.Math.Sqr(state.w[i])*2*(state.f-state.tasky[i]);
961            for(i_=0; i_<=k-1;i_++)
962            {
963                state.optstate.g[i_] = state.optstate.g[i_] + v*state.g[i_];
964            }
965            i = i+1;
966            goto lbl_14;
967        lbl_16:
968            goto lbl_5;
969        lbl_12:
970            if( ! state.optstate.needfij )
971            {
972                goto lbl_17;
973            }
974           
975            //
976            // calculate Fi/jac(Fi)
977            //
978            i = 0;
979        lbl_19:
980            if( i>n-1 )
981            {
982                goto lbl_21;
983            }
984            for(i_=0; i_<=k-1;i_++)
985            {
986                state.c[i_] = state.optstate.x[i_];
987            }
988            for(i_=0; i_<=m-1;i_++)
989            {
990                state.x[i_] = state.taskx[i,i_];
991            }
992            state.pointindex = i;
993            lsfitclearrequestfields(ref state);
994            state.needfg = true;
995            state.rstate.stage = 2;
996            goto lbl_rcomm;
997        lbl_2:
998            state.optstate.fi[i] = state.w[i]*(state.f-state.tasky[i]);
999            v = state.w[i];
1000            for(i_=0; i_<=k-1;i_++)
1001            {
1002                state.optstate.j[i,i_] = v*state.g[i_];
1003            }
1004            i = i+1;
1005            goto lbl_19;
1006        lbl_21:
1007            goto lbl_5;
1008        lbl_17:
1009            if( ! state.optstate.needfgh )
1010            {
1011                goto lbl_22;
1012            }
1013           
1014            //
1015            // calculate F/grad(F)/hess(F)
1016            //
1017            state.optstate.f = 0;
1018            for(i=0; i<=k-1; i++)
1019            {
1020                state.optstate.g[i] = 0;
1021            }
1022            for(i=0; i<=k-1; i++)
1023            {
1024                for(j=0; j<=k-1; j++)
1025                {
1026                    state.optstate.h[i,j] = 0;
1027                }
1028            }
1029            i = 0;
1030        lbl_24:
1031            if( i>n-1 )
1032            {
1033                goto lbl_26;
1034            }
1035            for(i_=0; i_<=k-1;i_++)
1036            {
1037                state.c[i_] = state.optstate.x[i_];
1038            }
1039            for(i_=0; i_<=m-1;i_++)
1040            {
1041                state.x[i_] = state.taskx[i,i_];
1042            }
1043            state.pointindex = i;
1044            lsfitclearrequestfields(ref state);
1045            state.needfgh = true;
1046            state.rstate.stage = 3;
1047            goto lbl_rcomm;
1048        lbl_3:
1049            state.optstate.f = state.optstate.f+AP.Math.Sqr(state.w[i]*(state.f-state.tasky[i]));
1050            v = AP.Math.Sqr(state.w[i])*2*(state.f-state.tasky[i]);
1051            for(i_=0; i_<=k-1;i_++)
1052            {
1053                state.optstate.g[i_] = state.optstate.g[i_] + v*state.g[i_];
1054            }
1055            for(j=0; j<=k-1; j++)
1056            {
1057                v = 2*AP.Math.Sqr(state.w[i])*state.g[j];
1058                for(i_=0; i_<=k-1;i_++)
1059                {
1060                    state.optstate.h[j,i_] = state.optstate.h[j,i_] + v*state.g[i_];
1061                }
1062                v = 2*AP.Math.Sqr(state.w[i])*(state.f-state.tasky[i]);
1063                for(i_=0; i_<=k-1;i_++)
1064                {
1065                    state.optstate.h[j,i_] = state.optstate.h[j,i_] + v*state.h[j,i_];
1066                }
1067            }
1068            i = i+1;
1069            goto lbl_24;
1070        lbl_26:
1071            goto lbl_5;
1072        lbl_22:
1073            goto lbl_5;
1074        lbl_6:
1075            minlm.minlmresults(ref state.optstate, ref state.c, ref state.optrep);
1076            state.repterminationtype = state.optrep.terminationtype;
1077           
1078            //
1079            // calculate errors
1080            //
1081            if( state.repterminationtype<=0 )
1082            {
1083                goto lbl_27;
1084            }
1085            state.reprmserror = 0;
1086            state.repavgerror = 0;
1087            state.repavgrelerror = 0;
1088            state.repmaxerror = 0;
1089            relcnt = 0;
1090            i = 0;
1091        lbl_29:
1092            if( i>n-1 )
1093            {
1094                goto lbl_31;
1095            }
1096            for(i_=0; i_<=k-1;i_++)
1097            {
1098                state.c[i_] = state.c[i_];
1099            }
1100            for(i_=0; i_<=m-1;i_++)
1101            {
1102                state.x[i_] = state.taskx[i,i_];
1103            }
1104            state.pointindex = i;
1105            lsfitclearrequestfields(ref state);
1106            state.needf = true;
1107            state.rstate.stage = 4;
1108            goto lbl_rcomm;
1109        lbl_4:
1110            v = state.f;
1111            state.reprmserror = state.reprmserror+AP.Math.Sqr(v-state.tasky[i]);
1112            state.repavgerror = state.repavgerror+Math.Abs(v-state.tasky[i]);
1113            if( (double)(state.tasky[i])!=(double)(0) )
1114            {
1115                state.repavgrelerror = state.repavgrelerror+Math.Abs(v-state.tasky[i])/Math.Abs(state.tasky[i]);
1116                relcnt = relcnt+1;
1117            }
1118            state.repmaxerror = Math.Max(state.repmaxerror, Math.Abs(v-state.tasky[i]));
1119            i = i+1;
1120            goto lbl_29;
1121        lbl_31:
1122            state.reprmserror = Math.Sqrt(state.reprmserror/n);
1123            state.repavgerror = state.repavgerror/n;
1124            if( (double)(relcnt)!=(double)(0) )
1125            {
1126                state.repavgrelerror = state.repavgrelerror/relcnt;
1127            }
1128        lbl_27:
1129            result = false;
1130            return result;
1131           
1132            //
1133            // Saving state
1134            //
1135        lbl_rcomm:
1136            result = true;
1137            state.rstate.ia[0] = n;
1138            state.rstate.ia[1] = m;
1139            state.rstate.ia[2] = k;
1140            state.rstate.ia[3] = i;
1141            state.rstate.ia[4] = j;
1142            state.rstate.ra[0] = v;
1143            state.rstate.ra[1] = relcnt;
1144            return result;
1145        }
1146
1147
1148        /*************************************************************************
1149        Nonlinear least squares fitting results.
1150
1151        Called after LSFitNonlinearIteration() returned False.
1152
1153        INPUT PARAMETERS:
1154            State   -   algorithm state (used by LSFitNonlinearIteration).
1155
1156        OUTPUT PARAMETERS:
1157            Info    -   completetion code:
1158                            * -1    incorrect parameters were specified
1159                            *  1    relative function improvement is no more than
1160                                    EpsF.
1161                            *  2    relative step is no more than EpsX.
1162                            *  4    gradient norm is no more than EpsG
1163                            *  5    MaxIts steps was taken
1164            C       -   array[0..K-1], solution
1165            Rep     -   optimization report. Following fields are set:
1166                        * Rep.TerminationType completetion code:
1167                        * RMSError          rms error on the (X,Y).
1168                        * AvgError          average error on the (X,Y).
1169                        * AvgRelError       average relative error on the non-zero Y
1170                        * MaxError          maximum error
1171                                            NON-WEIGHTED ERRORS ARE CALCULATED
1172
1173
1174          -- ALGLIB --
1175             Copyright 17.08.2009 by Bochkanov Sergey
1176        *************************************************************************/
1177        public static void lsfitnonlinearresults(ref lsfitstate state,
1178            ref int info,
1179            ref double[] c,
1180            ref lsfitreport rep)
1181        {
1182            int i_ = 0;
1183
1184            info = state.repterminationtype;
1185            if( info>0 )
1186            {
1187                c = new double[state.k];
1188                for(i_=0; i_<=state.k-1;i_++)
1189                {
1190                    c[i_] = state.c[i_];
1191                }
1192                rep.rmserror = state.reprmserror;
1193                rep.avgerror = state.repavgerror;
1194                rep.avgrelerror = state.repavgrelerror;
1195                rep.maxerror = state.repmaxerror;
1196            }
1197        }
1198
1199
1200        public static void lsfitscalexy(ref double[] x,
1201            ref double[] y,
1202            int n,
1203            ref double[] xc,
1204            ref double[] yc,
1205            ref int[] dc,
1206            int k,
1207            ref double xa,
1208            ref double xb,
1209            ref double sa,
1210            ref double sb,
1211            ref double[] xoriginal,
1212            ref double[] yoriginal)
1213        {
1214            double xmin = 0;
1215            double xmax = 0;
1216            int i = 0;
1217            int i_ = 0;
1218
1219            System.Diagnostics.Debug.Assert(n>=1, "LSFitScaleXY: incorrect N");
1220            System.Diagnostics.Debug.Assert(k>=0, "LSFitScaleXY: incorrect K");
1221           
1222            //
1223            // Calculate xmin/xmax.
1224            // Force xmin<>xmax.
1225            //
1226            xmin = x[0];
1227            xmax = x[0];
1228            for(i=1; i<=n-1; i++)
1229            {
1230                xmin = Math.Min(xmin, x[i]);
1231                xmax = Math.Max(xmax, x[i]);
1232            }
1233            for(i=0; i<=k-1; i++)
1234            {
1235                xmin = Math.Min(xmin, xc[i]);
1236                xmax = Math.Max(xmax, xc[i]);
1237            }
1238            if( (double)(xmin)==(double)(xmax) )
1239            {
1240                if( (double)(xmin)==(double)(0) )
1241                {
1242                    xmin = -1;
1243                    xmax = +1;
1244                }
1245                else
1246                {
1247                    xmin = 0.5*xmin;
1248                }
1249            }
1250           
1251            //
1252            // Transform abscissas: map [XA,XB] to [0,1]
1253            //
1254            // Store old X[] in XOriginal[] (it will be used
1255            // to calculate relative error).
1256            //
1257            xoriginal = new double[n];
1258            for(i_=0; i_<=n-1;i_++)
1259            {
1260                xoriginal[i_] = x[i_];
1261            }
1262            xa = xmin;
1263            xb = xmax;
1264            for(i=0; i<=n-1; i++)
1265            {
1266                x[i] = 2*(x[i]-0.5*(xa+xb))/(xb-xa);
1267            }
1268            for(i=0; i<=k-1; i++)
1269            {
1270                System.Diagnostics.Debug.Assert(dc[i]>=0, "LSFitScaleXY: internal error!");
1271                xc[i] = 2*(xc[i]-0.5*(xa+xb))/(xb-xa);
1272                yc[i] = yc[i]*Math.Pow(0.5*(xb-xa), dc[i]);
1273            }
1274           
1275            //
1276            // Transform function values: map [SA,SB] to [0,1]
1277            // SA = mean(Y),
1278            // SB = SA+stddev(Y).
1279            //
1280            // Store old Y[] in YOriginal[] (it will be used
1281            // to calculate relative error).
1282            //
1283            yoriginal = new double[n];
1284            for(i_=0; i_<=n-1;i_++)
1285            {
1286                yoriginal[i_] = y[i_];
1287            }
1288            sa = 0;
1289            for(i=0; i<=n-1; i++)
1290            {
1291                sa = sa+y[i];
1292            }
1293            sa = sa/n;
1294            sb = 0;
1295            for(i=0; i<=n-1; i++)
1296            {
1297                sb = sb+AP.Math.Sqr(y[i]-sa);
1298            }
1299            sb = Math.Sqrt(sb/n)+sa;
1300            if( (double)(sb)==(double)(sa) )
1301            {
1302                sb = 2*sa;
1303            }
1304            if( (double)(sb)==(double)(sa) )
1305            {
1306                sb = sa+1;
1307            }
1308            for(i=0; i<=n-1; i++)
1309            {
1310                y[i] = (y[i]-sa)/(sb-sa);
1311            }
1312            for(i=0; i<=k-1; i++)
1313            {
1314                if( dc[i]==0 )
1315                {
1316                    yc[i] = (yc[i]-sa)/(sb-sa);
1317                }
1318                else
1319                {
1320                    yc[i] = yc[i]/(sb-sa);
1321                }
1322            }
1323        }
1324
1325
1326        /*************************************************************************
1327        Internal fitting subroutine
1328        *************************************************************************/
1329        private static void lsfitlinearinternal(ref double[] y,
1330            ref double[] w,
1331            ref double[,] fmatrix,
1332            int n,
1333            int m,
1334            ref int info,
1335            ref double[] c,
1336            ref lsfitreport rep)
1337        {
1338            double threshold = 0;
1339            double[,] ft = new double[0,0];
1340            double[,] q = new double[0,0];
1341            double[,] l = new double[0,0];
1342            double[,] r = new double[0,0];
1343            double[] b = new double[0];
1344            double[] wmod = new double[0];
1345            double[] tau = new double[0];
1346            int i = 0;
1347            int j = 0;
1348            double v = 0;
1349            double[] sv = new double[0];
1350            double[,] u = new double[0,0];
1351            double[,] vt = new double[0,0];
1352            double[] tmp = new double[0];
1353            double[] utb = new double[0];
1354            double[] sutb = new double[0];
1355            int relcnt = 0;
1356            int i_ = 0;
1357
1358            if( n<1 | m<1 )
1359            {
1360                info = -1;
1361                return;
1362            }
1363            info = 1;
1364            threshold = Math.Sqrt(AP.Math.MachineEpsilon);
1365           
1366            //
1367            // Degenerate case, needs special handling
1368            //
1369            if( n<m )
1370            {
1371               
1372                //
1373                // Create design matrix.
1374                //
1375                ft = new double[n, m];
1376                b = new double[n];
1377                wmod = new double[n];
1378                for(j=0; j<=n-1; j++)
1379                {
1380                    v = w[j];
1381                    for(i_=0; i_<=m-1;i_++)
1382                    {
1383                        ft[j,i_] = v*fmatrix[j,i_];
1384                    }
1385                    b[j] = w[j]*y[j];
1386                    wmod[j] = 1;
1387                }
1388               
1389                //
1390                // LQ decomposition and reduction to M=N
1391                //
1392                c = new double[m];
1393                for(i=0; i<=m-1; i++)
1394                {
1395                    c[i] = 0;
1396                }
1397                rep.taskrcond = 0;
1398                ortfac.rmatrixlq(ref ft, n, m, ref tau);
1399                ortfac.rmatrixlqunpackq(ref ft, n, m, ref tau, n, ref q);
1400                ortfac.rmatrixlqunpackl(ref ft, n, m, ref l);
1401                lsfitlinearinternal(ref b, ref wmod, ref l, n, n, ref info, ref tmp, ref rep);
1402                if( info<=0 )
1403                {
1404                    return;
1405                }
1406                for(i=0; i<=n-1; i++)
1407                {
1408                    v = tmp[i];
1409                    for(i_=0; i_<=m-1;i_++)
1410                    {
1411                        c[i_] = c[i_] + v*q[i,i_];
1412                    }
1413                }
1414                return;
1415            }
1416           
1417            //
1418            // N>=M. Generate design matrix and reduce to N=M using
1419            // QR decomposition.
1420            //
1421            ft = new double[n, m];
1422            b = new double[n];
1423            for(j=0; j<=n-1; j++)
1424            {
1425                v = w[j];
1426                for(i_=0; i_<=m-1;i_++)
1427                {
1428                    ft[j,i_] = v*fmatrix[j,i_];
1429                }
1430                b[j] = w[j]*y[j];
1431            }
1432            ortfac.rmatrixqr(ref ft, n, m, ref tau);
1433            ortfac.rmatrixqrunpackq(ref ft, n, m, ref tau, m, ref q);
1434            ortfac.rmatrixqrunpackr(ref ft, n, m, ref r);
1435            tmp = new double[m];
1436            for(i=0; i<=m-1; i++)
1437            {
1438                tmp[i] = 0;
1439            }
1440            for(i=0; i<=n-1; i++)
1441            {
1442                v = b[i];
1443                for(i_=0; i_<=m-1;i_++)
1444                {
1445                    tmp[i_] = tmp[i_] + v*q[i,i_];
1446                }
1447            }
1448            b = new double[m];
1449            for(i_=0; i_<=m-1;i_++)
1450            {
1451                b[i_] = tmp[i_];
1452            }
1453           
1454            //
1455            // R contains reduced MxM design upper triangular matrix,
1456            // B contains reduced Mx1 right part.
1457            //
1458            // Determine system condition number and decide
1459            // should we use triangular solver (faster) or
1460            // SVD-based solver (more stable).
1461            //
1462            // We can use LU-based RCond estimator for this task.
1463            //
1464            rep.taskrcond = rcond.rmatrixlurcondinf(ref r, m);
1465            if( (double)(rep.taskrcond)>(double)(threshold) )
1466            {
1467               
1468                //
1469                // use QR-based solver
1470                //
1471                c = new double[m];
1472                c[m-1] = b[m-1]/r[m-1,m-1];
1473                for(i=m-2; i>=0; i--)
1474                {
1475                    v = 0.0;
1476                    for(i_=i+1; i_<=m-1;i_++)
1477                    {
1478                        v += r[i,i_]*c[i_];
1479                    }
1480                    c[i] = (b[i]-v)/r[i,i];
1481                }
1482            }
1483            else
1484            {
1485               
1486                //
1487                // use SVD-based solver
1488                //
1489                if( !svd.rmatrixsvd(r, m, m, 1, 1, 2, ref sv, ref u, ref vt) )
1490                {
1491                    info = -4;
1492                    return;
1493                }
1494                utb = new double[m];
1495                sutb = new double[m];
1496                for(i=0; i<=m-1; i++)
1497                {
1498                    utb[i] = 0;
1499                }
1500                for(i=0; i<=m-1; i++)
1501                {
1502                    v = b[i];
1503                    for(i_=0; i_<=m-1;i_++)
1504                    {
1505                        utb[i_] = utb[i_] + v*u[i,i_];
1506                    }
1507                }
1508                if( (double)(sv[0])>(double)(0) )
1509                {
1510                    rep.taskrcond = sv[m-1]/sv[0];
1511                    for(i=0; i<=m-1; i++)
1512                    {
1513                        if( (double)(sv[i])>(double)(threshold*sv[0]) )
1514                        {
1515                            sutb[i] = utb[i]/sv[i];
1516                        }
1517                        else
1518                        {
1519                            sutb[i] = 0;
1520                        }
1521                    }
1522                }
1523                else
1524                {
1525                    rep.taskrcond = 0;
1526                    for(i=0; i<=m-1; i++)
1527                    {
1528                        sutb[i] = 0;
1529                    }
1530                }
1531                c = new double[m];
1532                for(i=0; i<=m-1; i++)
1533                {
1534                    c[i] = 0;
1535                }
1536                for(i=0; i<=m-1; i++)
1537                {
1538                    v = sutb[i];
1539                    for(i_=0; i_<=m-1;i_++)
1540                    {
1541                        c[i_] = c[i_] + v*vt[i,i_];
1542                    }
1543                }
1544            }
1545           
1546            //
1547            // calculate errors
1548            //
1549            rep.rmserror = 0;
1550            rep.avgerror = 0;
1551            rep.avgrelerror = 0;
1552            rep.maxerror = 0;
1553            relcnt = 0;
1554            for(i=0; i<=n-1; i++)
1555            {
1556                v = 0.0;
1557                for(i_=0; i_<=m-1;i_++)
1558                {
1559                    v += fmatrix[i,i_]*c[i_];
1560                }
1561                rep.rmserror = rep.rmserror+AP.Math.Sqr(v-y[i]);
1562                rep.avgerror = rep.avgerror+Math.Abs(v-y[i]);
1563                if( (double)(y[i])!=(double)(0) )
1564                {
1565                    rep.avgrelerror = rep.avgrelerror+Math.Abs(v-y[i])/Math.Abs(y[i]);
1566                    relcnt = relcnt+1;
1567                }
1568                rep.maxerror = Math.Max(rep.maxerror, Math.Abs(v-y[i]));
1569            }
1570            rep.rmserror = Math.Sqrt(rep.rmserror/n);
1571            rep.avgerror = rep.avgerror/n;
1572            if( relcnt!=0 )
1573            {
1574                rep.avgrelerror = rep.avgrelerror/relcnt;
1575            }
1576        }
1577
1578
1579        /*************************************************************************
1580        Internal subroutine
1581        *************************************************************************/
1582        private static void lsfitclearrequestfields(ref lsfitstate state)
1583        {
1584            state.needf = false;
1585            state.needfg = false;
1586            state.needfgh = false;
1587        }
1588    }
1589}
Note: See TracBrowser for help on using the repository browser.