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