Free cookie consent management tool by TermsFeed Policy Generator

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

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

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

File size: 35.2 KB
Line 
1/*************************************************************************
2Copyright (c) 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 minlm
26    {
27        public struct lmstate
28        {
29            public bool wrongparams;
30            public int n;
31            public int m;
32            public double epsf;
33            public double epsx;
34            public int maxits;
35            public int flags;
36            public int usermode;
37            public double[] x;
38            public double f;
39            public double[] fi;
40            public double[,] j;
41            public double[,] h;
42            public double[] g;
43            public bool needf;
44            public bool needfg;
45            public bool needfgh;
46            public bool needfij;
47            public bool xupdated;
48            public lbfgs.lbfgsstate internalstate;
49            public lbfgs.lbfgsreport internalrep;
50            public double[] xprec;
51            public double[] xbase;
52            public double[] xdir;
53            public double[] gbase;
54            public double[,] rawmodel;
55            public double[,] model;
56            public double[] work;
57            public AP.rcommstate rstate;
58            public int repiterationscount;
59            public int repterminationtype;
60            public int repnfunc;
61            public int repnjac;
62            public int repngrad;
63            public int repnhess;
64            public int repncholesky;
65        };
66
67
68        public struct lmreport
69        {
70            public int iterationscount;
71            public int terminationtype;
72            public int nfunc;
73            public int njac;
74            public int ngrad;
75            public int nhess;
76            public int ncholesky;
77        };
78
79
80
81
82        public const int lmmodefj = 0;
83        public const int lmmodefgj = 1;
84        public const int lmmodefgh = 2;
85        public const int lmflagnoprelbfgs = 1;
86        public const int lmflagnointlbfgs = 2;
87        public const int lmprelbfgsm = 5;
88        public const int lmintlbfgsits = 5;
89        public const int lbfgsnorealloc = 1;
90
91
92        /*************************************************************************
93            LEVENBERG-MARQUARDT-LIKE METHOD FOR NON-LINEAR OPTIMIZATION
94
95        Optimization using function gradient and Hessian.  Algorithm -  Levenberg-
96        Marquardt   modification   with   L-BFGS   pre-optimization  and  internal
97        pre-conditioned L-BFGS optimization after each Levenberg-Marquardt step.
98
99        Function F has general form (not "sum-of-squares"):
100
101            F = F(x[0], ..., x[n-1])
102
103        EXAMPLE
104
105        See HTML-documentation.
106
107        INPUT PARAMETERS:
108            N       -   dimension, N>1
109            X       -   initial solution, array[0..N-1]
110            EpsF    -   stopping criterion. Algorithm stops if
111                        |F(k+1)-F(k)| <= EpsF*max{|F(k)|, |F(k+1)|, 1}
112            EpsX    -   stopping criterion. Algorithm stops if
113                        |X(k+1)-X(k)| <= EpsX*(1+|X(k)|)
114            MaxIts  -   stopping criterion. Algorithm stops after MaxIts iterations.
115                        MaxIts=0 means no stopping criterion.
116
117        OUTPUT PARAMETERS:
118            State   -   structure which stores algorithm state between subsequent
119                        calls of MinLMIteration. Used for reverse communication.
120                        This structure should be passed to MinLMIteration subroutine.
121
122        See also MinLMIteration, MinLMResults.
123
124        NOTE
125
126        Passing EpsF=0, EpsX=0 and MaxIts=0 (simultaneously) will lead to automatic
127        stopping criterion selection (small EpsX).
128
129          -- ALGLIB --
130             Copyright 30.03.2009 by Bochkanov Sergey
131        *************************************************************************/
132        public static void minlmfgh(int n,
133            ref double[] x,
134            double epsf,
135            double epsx,
136            int maxits,
137            ref lmstate state)
138        {
139            int i_ = 0;
140
141           
142            //
143            // Prepare RComm
144            //
145            state.rstate.ia = new int[3+1];
146            state.rstate.ba = new bool[0+1];
147            state.rstate.ra = new double[8+1];
148            state.rstate.stage = -1;
149           
150            //
151            // prepare internal structures
152            //
153            lmprepare(n, 0, true, ref state);
154           
155            //
156            // initialize, check parameters
157            //
158            state.xupdated = false;
159            state.n = n;
160            state.m = 0;
161            state.epsf = epsf;
162            state.epsx = epsx;
163            state.maxits = maxits;
164            state.flags = 0;
165            if( (double)(state.epsf)==(double)(0) & (double)(state.epsx)==(double)(0) & state.maxits==0 )
166            {
167                state.epsx = 1.0E-6;
168            }
169            state.usermode = lmmodefgh;
170            state.wrongparams = false;
171            if( n<1 | (double)(epsf)<(double)(0) | (double)(epsx)<(double)(0) | maxits<0 )
172            {
173                state.wrongparams = true;
174                return;
175            }
176            for(i_=0; i_<=n-1;i_++)
177            {
178                state.x[i_] = x[i_];
179            }
180        }
181
182
183        /*************************************************************************
184            LEVENBERG-MARQUARDT-LIKE METHOD FOR NON-LINEAR OPTIMIZATION
185
186        Optimization using function gradient and Jacobian.  Algorithm -  Levenberg-
187        Marquardt   modification   with   L-BFGS   pre-optimization  and  internal
188        pre-conditioned L-BFGS optimization after each Levenberg-Marquardt step.
189
190        Function F is represented as sum of squares:
191
192            F = f[0]^2(x[0],...,x[n-1]) + ... + f[m-1]^2(x[0],...,x[n-1])
193
194        EXAMPLE
195
196        See HTML-documentation.
197
198        INPUT PARAMETERS:
199            N       -   dimension, N>1
200            M       -   number of functions f[i]
201            X       -   initial solution, array[0..N-1]
202            EpsF    -   stopping criterion. Algorithm stops if
203                        |F(k+1)-F(k)| <= EpsF*max{|F(k)|, |F(k+1)|, 1}
204            EpsX    -   stopping criterion. Algorithm stops if
205                        |X(k+1)-X(k)| <= EpsX*(1+|X(k)|)
206            MaxIts  -   stopping criterion. Algorithm stops after MaxIts iterations.
207                        MaxIts=0 means no stopping criterion.
208
209        OUTPUT PARAMETERS:
210            State   -   structure which stores algorithm state between subsequent
211                        calls of MinLMIteration. Used for reverse communication.
212                        This structure should be passed to MinLMIteration subroutine.
213
214        See also MinLMIteration, MinLMResults.
215
216        NOTE
217
218        Passing EpsF=0, EpsX=0 and MaxIts=0 (simultaneously) will lead to automatic
219        stopping criterion selection (small EpsX).
220
221          -- ALGLIB --
222             Copyright 30.03.2009 by Bochkanov Sergey
223        *************************************************************************/
224        public static void minlmfgj(int n,
225            int m,
226            ref double[] x,
227            double epsf,
228            double epsx,
229            int maxits,
230            ref lmstate state)
231        {
232            int i_ = 0;
233
234           
235            //
236            // Prepare RComm
237            //
238            state.rstate.ia = new int[3+1];
239            state.rstate.ba = new bool[0+1];
240            state.rstate.ra = new double[8+1];
241            state.rstate.stage = -1;
242           
243            //
244            // prepare internal structures
245            //
246            lmprepare(n, m, true, ref state);
247           
248            //
249            // initialize, check parameters
250            //
251            state.xupdated = false;
252            state.n = n;
253            state.m = m;
254            state.epsf = epsf;
255            state.epsx = epsx;
256            state.maxits = maxits;
257            state.flags = 0;
258            if( (double)(state.epsf)==(double)(0) & (double)(state.epsx)==(double)(0) & state.maxits==0 )
259            {
260                state.epsx = 1.0E-6;
261            }
262            state.usermode = lmmodefgj;
263            state.wrongparams = false;
264            if( n<1 | m<1 | (double)(epsf)<(double)(0) | (double)(epsx)<(double)(0) | maxits<0 )
265            {
266                state.wrongparams = true;
267                return;
268            }
269            for(i_=0; i_<=n-1;i_++)
270            {
271                state.x[i_] = x[i_];
272            }
273        }
274
275
276        /*************************************************************************
277            CLASSIC LEVENBERG-MARQUARDT METHOD FOR NON-LINEAR OPTIMIZATION
278
279        Optimization using Jacobi matrix. Algorithm  -  classic Levenberg-Marquardt
280        method.
281
282        Function F is represented as sum of squares:
283
284            F = f[0]^2(x[0],...,x[n-1]) + ... + f[m-1]^2(x[0],...,x[n-1])
285
286        EXAMPLE
287
288        See HTML-documentation.
289
290        INPUT PARAMETERS:
291            N       -   dimension, N>1
292            M       -   number of functions f[i]
293            X       -   initial solution, array[0..N-1]
294            EpsF    -   stopping criterion. Algorithm stops if
295                        |F(k+1)-F(k)| <= EpsF*max{|F(k)|, |F(k+1)|, 1}
296            EpsX    -   stopping criterion. Algorithm stops if
297                        |X(k+1)-X(k)| <= EpsX*(1+|X(k)|)
298            MaxIts  -   stopping criterion. Algorithm stops after MaxIts iterations.
299                        MaxIts=0 means no stopping criterion.
300
301        OUTPUT PARAMETERS:
302            State   -   structure which stores algorithm state between subsequent
303                        calls of MinLMIteration. Used for reverse communication.
304                        This structure should be passed to MinLMIteration subroutine.
305
306        See also MinLMIteration, MinLMResults.
307
308        NOTE
309
310        Passing EpsF=0, EpsX=0 and MaxIts=0 (simultaneously) will lead to automatic
311        stopping criterion selection (small EpsX).
312
313          -- ALGLIB --
314             Copyright 30.03.2009 by Bochkanov Sergey
315        *************************************************************************/
316        public static void minlmfj(int n,
317            int m,
318            ref double[] x,
319            double epsf,
320            double epsx,
321            int maxits,
322            ref lmstate state)
323        {
324            int i_ = 0;
325
326           
327            //
328            // Prepare RComm
329            //
330            state.rstate.ia = new int[3+1];
331            state.rstate.ba = new bool[0+1];
332            state.rstate.ra = new double[8+1];
333            state.rstate.stage = -1;
334           
335            //
336            // prepare internal structures
337            //
338            lmprepare(n, m, true, ref state);
339           
340            //
341            // initialize, check parameters
342            //
343            state.xupdated = false;
344            state.n = n;
345            state.m = m;
346            state.epsf = epsf;
347            state.epsx = epsx;
348            state.maxits = maxits;
349            state.flags = 0;
350            if( (double)(state.epsf)==(double)(0) & (double)(state.epsx)==(double)(0) & state.maxits==0 )
351            {
352                state.epsx = 1.0E-6;
353            }
354            state.usermode = lmmodefj;
355            state.wrongparams = false;
356            if( n<1 | m<1 | (double)(epsf)<(double)(0) | (double)(epsx)<(double)(0) | maxits<0 )
357            {
358                state.wrongparams = true;
359                return;
360            }
361            for(i_=0; i_<=n-1;i_++)
362            {
363                state.x[i_] = x[i_];
364            }
365        }
366
367
368        /*************************************************************************
369        One Levenberg-Marquardt iteration.
370
371        Called after inialization of State structure with MinLMXXX subroutine.
372        See HTML docs for examples.
373
374        Input parameters:
375            State   -   structure which stores algorithm state between subsequent
376                        calls and which is used for reverse communication. Must be
377                        initialized with MinLMXXX call first.
378
379        If subroutine returned False, iterative algorithm has converged.
380
381        If subroutine returned True, then:
382        * if State.NeedF=True,      -   function value F at State.X[0..N-1]
383                                        is required
384        * if State.NeedFG=True      -   function value F and gradient G
385                                        are required
386        * if State.NeedFiJ=True     -   function vector f[i] and Jacobi matrix J
387                                        are required
388        * if State.NeedFGH=True     -   function value F, gradient G and Hesian H
389                                        are required
390
391        One and only one of this fields can be set at time.
392
393        Results are stored:
394        * function value            -   in LMState.F
395        * gradient                  -   in LMState.G[0..N-1]
396        * Jacobi matrix             -   in LMState.J[0..M-1,0..N-1]
397        * Hessian                   -   in LMState.H[0..N-1,0..N-1]
398
399          -- ALGLIB --
400             Copyright 10.03.2009 by Bochkanov Sergey
401        *************************************************************************/
402        public static bool minlmiteration(ref lmstate state)
403        {
404            bool result = new bool();
405            int n = 0;
406            int m = 0;
407            int i = 0;
408            double xnorm = 0;
409            double stepnorm = 0;
410            bool spd = new bool();
411            double fbase = 0;
412            double fnew = 0;
413            double lambda = 0;
414            double nu = 0;
415            double lambdaup = 0;
416            double lambdadown = 0;
417            int lbfgsflags = 0;
418            double v = 0;
419            int i_ = 0;
420
421           
422            //
423            // Reverse communication preparations
424            // I know it looks ugly, but it works the same way
425            // anywhere from C++ to Python.
426            //
427            // This code initializes locals by:
428            // * random values determined during code
429            //   generation - on first subroutine call
430            // * values from previous call - on subsequent calls
431            //
432            if( state.rstate.stage>=0 )
433            {
434                n = state.rstate.ia[0];
435                m = state.rstate.ia[1];
436                i = state.rstate.ia[2];
437                lbfgsflags = state.rstate.ia[3];
438                spd = state.rstate.ba[0];
439                xnorm = state.rstate.ra[0];
440                stepnorm = state.rstate.ra[1];
441                fbase = state.rstate.ra[2];
442                fnew = state.rstate.ra[3];
443                lambda = state.rstate.ra[4];
444                nu = state.rstate.ra[5];
445                lambdaup = state.rstate.ra[6];
446                lambdadown = state.rstate.ra[7];
447                v = state.rstate.ra[8];
448            }
449            else
450            {
451                n = -983;
452                m = -989;
453                i = -834;
454                lbfgsflags = 900;
455                spd = true;
456                xnorm = 364;
457                stepnorm = 214;
458                fbase = -338;
459                fnew = -686;
460                lambda = 912;
461                nu = 585;
462                lambdaup = 497;
463                lambdadown = -271;
464                v = -581;
465            }
466            if( state.rstate.stage==0 )
467            {
468                goto lbl_0;
469            }
470            if( state.rstate.stage==1 )
471            {
472                goto lbl_1;
473            }
474            if( state.rstate.stage==2 )
475            {
476                goto lbl_2;
477            }
478            if( state.rstate.stage==3 )
479            {
480                goto lbl_3;
481            }
482            if( state.rstate.stage==4 )
483            {
484                goto lbl_4;
485            }
486            if( state.rstate.stage==5 )
487            {
488                goto lbl_5;
489            }
490            if( state.rstate.stage==6 )
491            {
492                goto lbl_6;
493            }
494           
495            //
496            // Routine body
497            //
498            System.Diagnostics.Debug.Assert(state.usermode==lmmodefj | state.usermode==lmmodefgj | state.usermode==lmmodefgh, "LM: internal error");
499            if( state.wrongparams )
500            {
501                state.repterminationtype = -1;
502                result = false;
503                return result;
504            }
505           
506            //
507            // prepare params
508            //
509            n = state.n;
510            m = state.m;
511            lambdaup = 10;
512            lambdadown = 0.3;
513            nu = 2;
514            lbfgsflags = 0;
515           
516            //
517            // if we have F and G
518            //
519            if( ! ((state.usermode==lmmodefgj | state.usermode==lmmodefgh) & state.flags/lmflagnoprelbfgs%2==0) )
520            {
521                goto lbl_7;
522            }
523           
524            //
525            // First stage of the hybrid algorithm: LBFGS
526            //
527            lbfgs.minlbfgs(n, Math.Min(n, lmprelbfgsm), ref state.x, 0.0, 0.0, 0.0, Math.Max(25, n), 0, ref state.internalstate);
528        lbl_9:
529            if( ! lbfgs.minlbfgsiteration(ref state.internalstate) )
530            {
531                goto lbl_10;
532            }
533           
534            //
535            // RComm
536            //
537            for(i_=0; i_<=n-1;i_++)
538            {
539                state.x[i_] = state.internalstate.x[i_];
540            }
541            lmclearrequestfields(ref state);
542            state.needfg = true;
543            state.rstate.stage = 0;
544            goto lbl_rcomm;
545        lbl_0:
546            state.repnfunc = state.repnfunc+1;
547            state.repngrad = state.repngrad+1;
548           
549            //
550            // Call LBFGS
551            //
552            state.internalstate.f = state.f;
553            for(i_=0; i_<=n-1;i_++)
554            {
555                state.internalstate.g[i_] = state.g[i_];
556            }
557            goto lbl_9;
558        lbl_10:
559            lbfgs.minlbfgsresults(ref state.internalstate, ref state.x, ref state.internalrep);
560        lbl_7:
561           
562            //
563            // Second stage of the hybrid algorithm: LM
564            // Initialize quadratic model.
565            //
566            if( state.usermode!=lmmodefgh )
567            {
568                goto lbl_11;
569            }
570           
571            //
572            // RComm
573            //
574            lmclearrequestfields(ref state);
575            state.needfgh = true;
576            state.rstate.stage = 1;
577            goto lbl_rcomm;
578        lbl_1:
579            state.repnfunc = state.repnfunc+1;
580            state.repngrad = state.repngrad+1;
581            state.repnhess = state.repnhess+1;
582           
583            //
584            // generate raw quadratic model
585            //
586            for(i=0; i<=n-1; i++)
587            {
588                for(i_=0; i_<=n-1;i_++)
589                {
590                    state.rawmodel[i,i_] = state.h[i,i_];
591                }
592            }
593            for(i_=0; i_<=n-1;i_++)
594            {
595                state.gbase[i_] = state.g[i_];
596            }
597            fbase = state.f;
598        lbl_11:
599            if( ! (state.usermode==lmmodefgj | state.usermode==lmmodefj) )
600            {
601                goto lbl_13;
602            }
603           
604            //
605            // RComm
606            //
607            lmclearrequestfields(ref state);
608            state.needfij = true;
609            state.rstate.stage = 2;
610            goto lbl_rcomm;
611        lbl_2:
612            state.repnfunc = state.repnfunc+1;
613            state.repnjac = state.repnjac+1;
614           
615            //
616            // generate raw quadratic model
617            //
618            blas.matrixmatrixmultiply(ref state.j, 0, m-1, 0, n-1, true, ref state.j, 0, m-1, 0, n-1, false, 1.0, ref state.rawmodel, 0, n-1, 0, n-1, 0.0, ref state.work);
619            blas.matrixvectormultiply(ref state.j, 0, m-1, 0, n-1, true, ref state.fi, 0, m-1, 1.0, ref state.gbase, 0, n-1, 0.0);
620            fbase = 0.0;
621            for(i_=0; i_<=m-1;i_++)
622            {
623                fbase += state.fi[i_]*state.fi[i_];
624            }
625        lbl_13:
626            lambda = 0.001;
627        lbl_15:
628            if( false )
629            {
630                goto lbl_16;
631            }
632           
633            //
634            // 1. Model = RawModel+lambda*I
635            // 2. Try to solve (RawModel+Lambda*I)*dx = -g.
636            //    Increase lambda if left part is not positive definite.
637            //
638            for(i=0; i<=n-1; i++)
639            {
640                for(i_=0; i_<=n-1;i_++)
641                {
642                    state.model[i,i_] = state.rawmodel[i,i_];
643                }
644                state.model[i,i] = state.model[i,i]+lambda;
645            }
646            spd = cholesky.spdmatrixcholesky(ref state.model, n, true);
647            state.repncholesky = state.repncholesky+1;
648            if( !spd )
649            {
650                lambda = lambda*lambdaup*nu;
651                nu = nu*2;
652                goto lbl_15;
653            }
654            if( !spdsolve.spdmatrixcholeskysolve(ref state.model, state.gbase, n, true, ref state.xdir) )
655            {
656                lambda = lambda*lambdaup*nu;
657                nu = nu*2;
658                goto lbl_15;
659            }
660            for(i_=0; i_<=n-1;i_++)
661            {
662                state.xdir[i_] = -1*state.xdir[i_];
663            }
664           
665            //
666            // Candidate lambda found.
667            // 1. Save old w in WBase
668            // 1. Test some stopping criterions
669            // 2. If error(w+wdir)>error(w), increase lambda
670            //
671            for(i_=0; i_<=n-1;i_++)
672            {
673                state.xbase[i_] = state.x[i_];
674            }
675            for(i_=0; i_<=n-1;i_++)
676            {
677                state.x[i_] = state.x[i_] + state.xdir[i_];
678            }
679            xnorm = 0.0;
680            for(i_=0; i_<=n-1;i_++)
681            {
682                xnorm += state.xbase[i_]*state.xbase[i_];
683            }
684            stepnorm = 0.0;
685            for(i_=0; i_<=n-1;i_++)
686            {
687                stepnorm += state.xdir[i_]*state.xdir[i_];
688            }
689            xnorm = Math.Sqrt(xnorm);
690            stepnorm = Math.Sqrt(stepnorm);
691            if( (double)(stepnorm)<=(double)(state.epsx*(1+xnorm)) )
692            {
693               
694                //
695                // step size if small enough
696                //
697                state.repterminationtype = 2;
698                goto lbl_16;
699            }
700            lmclearrequestfields(ref state);
701            state.needf = true;
702            state.rstate.stage = 3;
703            goto lbl_rcomm;
704        lbl_3:
705            state.repnfunc = state.repnfunc+1;
706            fnew = state.f;
707            if( (double)(Math.Abs(fnew-fbase))<=(double)(state.epsf*Math.Max(1, Math.Max(Math.Abs(fbase), Math.Abs(fnew)))) )
708            {
709               
710                //
711                // function change is small enough
712                //
713                state.repterminationtype = 1;
714                goto lbl_16;
715            }
716            if( (double)(fnew)>(double)(fbase) )
717            {
718               
719                //
720                // restore state and continue out search for lambda
721                //
722                for(i_=0; i_<=n-1;i_++)
723                {
724                    state.x[i_] = state.xbase[i_];
725                }
726                lambda = lambda*lambdaup*nu;
727                nu = nu*2;
728                goto lbl_15;
729            }
730            if( ! ((state.usermode==lmmodefgj | state.usermode==lmmodefgh) & state.flags/lmflagnointlbfgs%2==0) )
731            {
732                goto lbl_17;
733            }
734            System.Diagnostics.Debug.Assert(state.usermode==lmmodefgh);
735           
736            //
737            // Optimize using inv(cholesky(H)) as preconditioner
738            //
739            if( ! trinverse.rmatrixtrinverse(ref state.model, n, true, false) )
740            {
741                goto lbl_19;
742            }
743           
744            //
745            // if matrix can be inverted use it.
746            // just silently move to next iteration otherwise.
747            // (will be very, very rare, mostly for specially
748            // designed near-degenerate tasks)
749            //
750            for(i_=0; i_<=n-1;i_++)
751            {
752                state.xbase[i_] = state.x[i_];
753            }
754            for(i=0; i<=n-1; i++)
755            {
756                state.xprec[i] = 0;
757            }
758            lbfgs.minlbfgs(n, Math.Min(n, lmintlbfgsits), ref state.xprec, 0.0, 0.0, 0.0, lmintlbfgsits, lbfgsflags, ref state.internalstate);
759        lbl_21:
760            if( ! lbfgs.minlbfgsiteration(ref state.internalstate) )
761            {
762                goto lbl_22;
763            }
764           
765            //
766            // convert XPrec to unpreconditioned form, then call RComm.
767            //
768            for(i=0; i<=n-1; i++)
769            {
770                v = 0.0;
771                for(i_=i; i_<=n-1;i_++)
772                {
773                    v += state.internalstate.x[i_]*state.model[i,i_];
774                }
775                state.x[i] = state.xbase[i]+v;
776            }
777            lmclearrequestfields(ref state);
778            state.needfg = true;
779            state.rstate.stage = 4;
780            goto lbl_rcomm;
781        lbl_4:
782            state.repnfunc = state.repnfunc+1;
783            state.repngrad = state.repngrad+1;
784           
785            //
786            // 1. pass State.F to State.InternalState.F
787            // 2. convert gradient back to preconditioned form
788            //
789            state.internalstate.f = state.f;
790            for(i=0; i<=n-1; i++)
791            {
792                state.internalstate.g[i] = 0;
793            }
794            for(i=0; i<=n-1; i++)
795            {
796                v = state.g[i];
797                for(i_=i; i_<=n-1;i_++)
798                {
799                    state.internalstate.g[i_] = state.internalstate.g[i_] + v*state.model[i,i_];
800                }
801            }
802           
803            //
804            // next iteration
805            //
806            goto lbl_21;
807        lbl_22:
808           
809            //
810            // change LBFGS flags to NoRealloc.
811            // L-BFGS subroutine will use memory allocated from previous run.
812            // it is possible since all subsequent calls will be with same N/M.
813            //
814            lbfgsflags = lbfgsnorealloc;
815           
816            //
817            // back to unpreconditioned X
818            //
819            lbfgs.minlbfgsresults(ref state.internalstate, ref state.xprec, ref state.internalrep);
820            for(i=0; i<=n-1; i++)
821            {
822                v = 0.0;
823                for(i_=i; i_<=n-1;i_++)
824                {
825                    v += state.xprec[i_]*state.model[i,i_];
826                }
827                state.x[i] = state.xbase[i]+v;
828            }
829        lbl_19:
830        lbl_17:
831           
832            //
833            // Accept new position.
834            // Calculate Hessian
835            //
836            if( state.usermode!=lmmodefgh )
837            {
838                goto lbl_23;
839            }
840           
841            //
842            // RComm
843            //
844            lmclearrequestfields(ref state);
845            state.needfgh = true;
846            state.rstate.stage = 5;
847            goto lbl_rcomm;
848        lbl_5:
849            state.repnfunc = state.repnfunc+1;
850            state.repngrad = state.repngrad+1;
851            state.repnhess = state.repnhess+1;
852           
853            //
854            // Update raw quadratic model
855            //
856            for(i=0; i<=n-1; i++)
857            {
858                for(i_=0; i_<=n-1;i_++)
859                {
860                    state.rawmodel[i,i_] = state.h[i,i_];
861                }
862            }
863            for(i_=0; i_<=n-1;i_++)
864            {
865                state.gbase[i_] = state.g[i_];
866            }
867            fbase = state.f;
868        lbl_23:
869            if( ! (state.usermode==lmmodefgj | state.usermode==lmmodefj) )
870            {
871                goto lbl_25;
872            }
873           
874            //
875            // RComm
876            //
877            lmclearrequestfields(ref state);
878            state.needfij = true;
879            state.rstate.stage = 6;
880            goto lbl_rcomm;
881        lbl_6:
882            state.repnfunc = state.repnfunc+1;
883            state.repnjac = state.repnjac+1;
884           
885            //
886            // generate raw quadratic model
887            //
888            blas.matrixmatrixmultiply(ref state.j, 0, m-1, 0, n-1, true, ref state.j, 0, m-1, 0, n-1, false, 1.0, ref state.rawmodel, 0, n-1, 0, n-1, 0.0, ref state.work);
889            blas.matrixvectormultiply(ref state.j, 0, m-1, 0, n-1, true, ref state.fi, 0, m-1, 1.0, ref state.gbase, 0, n-1, 0.0);
890            fbase = 0.0;
891            for(i_=0; i_<=m-1;i_++)
892            {
893                fbase += state.fi[i_]*state.fi[i_];
894            }
895        lbl_25:
896            state.repiterationscount = state.repiterationscount+1;
897            if( state.repiterationscount>=state.maxits & state.maxits>0 )
898            {
899                state.repterminationtype = 5;
900                goto lbl_16;
901            }
902           
903            //
904            // Update lambda
905            //
906            lambda = lambda*lambdadown;
907            nu = 2;
908            goto lbl_15;
909        lbl_16:
910            result = false;
911            return result;
912           
913            //
914            // Saving state
915            //
916        lbl_rcomm:
917            result = true;
918            state.rstate.ia[0] = n;
919            state.rstate.ia[1] = m;
920            state.rstate.ia[2] = i;
921            state.rstate.ia[3] = lbfgsflags;
922            state.rstate.ba[0] = spd;
923            state.rstate.ra[0] = xnorm;
924            state.rstate.ra[1] = stepnorm;
925            state.rstate.ra[2] = fbase;
926            state.rstate.ra[3] = fnew;
927            state.rstate.ra[4] = lambda;
928            state.rstate.ra[5] = nu;
929            state.rstate.ra[6] = lambdaup;
930            state.rstate.ra[7] = lambdadown;
931            state.rstate.ra[8] = v;
932            return result;
933        }
934
935
936        /*************************************************************************
937        Levenberg-Marquardt algorithm results
938
939        Called after MinLMIteration returned False.
940
941        Input parameters:
942            State   -   algorithm state (used by MinLMIteration).
943
944        Output parameters:
945            X       -   array[0..N-1], solution
946            Rep     -   optimization report:
947                        * Rep.TerminationType completetion code:
948                            * -1    incorrect parameters were specified
949                            *  1    relative function improvement is no more than
950                                    EpsF.
951                            *  2    relative step is no more than EpsX.
952                            *  4    gradient norm is no more than EpsG
953                            *  5    MaxIts steps was taken
954                        * Rep.IterationsCount contains iterations count
955                        * Rep.NFunc     - number of function calculations
956                        * Rep.NJac      - number of Jacobi matrix calculations
957                        * Rep.NGrad     - number of gradient calculations
958                        * Rep.NHess     - number of Hessian calculations
959                        * Rep.NCholesky - number of Cholesky decomposition calculations
960
961          -- ALGLIB --
962             Copyright 10.03.2009 by Bochkanov Sergey
963        *************************************************************************/
964        public static void minlmresults(ref lmstate state,
965            ref double[] x,
966            ref lmreport rep)
967        {
968            int i_ = 0;
969
970            x = new double[state.n-1+1];
971            for(i_=0; i_<=state.n-1;i_++)
972            {
973                x[i_] = state.x[i_];
974            }
975            rep.iterationscount = state.repiterationscount;
976            rep.terminationtype = state.repterminationtype;
977            rep.nfunc = state.repnfunc;
978            rep.njac = state.repnjac;
979            rep.ngrad = state.repngrad;
980            rep.nhess = state.repnhess;
981            rep.ncholesky = state.repncholesky;
982        }
983
984
985        /*************************************************************************
986        Prepare internal structures (except for RComm).
987
988        Note: M must be zero for FGH mode, non-zero for FJ/FGJ mode.
989        *************************************************************************/
990        private static void lmprepare(int n,
991            int m,
992            bool havegrad,
993            ref lmstate state)
994        {
995            state.repiterationscount = 0;
996            state.repterminationtype = 0;
997            state.repnfunc = 0;
998            state.repnjac = 0;
999            state.repngrad = 0;
1000            state.repnhess = 0;
1001            state.repncholesky = 0;
1002            if( n<0 | m<0 )
1003            {
1004                return;
1005            }
1006            if( havegrad )
1007            {
1008                state.g = new double[n-1+1];
1009            }
1010            if( m!=0 )
1011            {
1012                state.j = new double[m-1+1, n-1+1];
1013                state.fi = new double[m-1+1];
1014                state.h = new double[0+1, 0+1];
1015            }
1016            else
1017            {
1018                state.j = new double[0+1, 0+1];
1019                state.fi = new double[0+1];
1020                state.h = new double[n-1+1, n-1+1];
1021            }
1022            state.x = new double[n-1+1];
1023            state.rawmodel = new double[n-1+1, n-1+1];
1024            state.model = new double[n-1+1, n-1+1];
1025            state.xbase = new double[n-1+1];
1026            state.xprec = new double[n-1+1];
1027            state.gbase = new double[n-1+1];
1028            state.xdir = new double[n-1+1];
1029            state.work = new double[Math.Max(n, m)+1];
1030        }
1031
1032
1033        /*************************************************************************
1034        Clears request fileds (to be sure that we don't forgot to clear something)
1035        *************************************************************************/
1036        private static void lmclearrequestfields(ref lmstate state)
1037        {
1038            state.needf = false;
1039            state.needfg = false;
1040            state.needfgh = false;
1041            state.needfij = false;
1042        }
1043    }
1044}
Note: See TracBrowser for help on using the repository browser.