Free cookie consent management tool by TermsFeed Policy Generator

source: branches/MathNetNumerics-Exploration-2789/HeuristicLab.Algorithms.DataAnalysis.Experimental/Splines.cs @ 15469

Last change on this file since 15469 was 15469, checked in by gkronber, 6 years ago

#2789 trying to get SBART to work correctly.

File size: 30.3 KB
Line 
1#region License Information
2/* HeuristicLab
3 * Copyright (C) 2002-2016 Heuristic and Evolutionary Algorithms Laboratory (HEAL)
4 *
5 * This file is part of HeuristicLab.
6 *
7 * HeuristicLab is free software: you can redistribute it and/or modify
8 * it under the terms of the GNU General Public License as published by
9 * the Free Software Foundation, either version 3 of the License, or
10 * (at your option) any later version.
11 *
12 * HeuristicLab is distributed in the hope that it will be useful,
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
15 * GNU General Public License for more details.
16 *
17 * You should have received a copy of the GNU General Public License
18 * along with HeuristicLab. If not, see <http://www.gnu.org/licenses/>.
19 */
20#endregion
21
22using System;
23using System.Collections.Generic;
24using System.Diagnostics;
25using System.Linq;
26using System.Threading;
27using HeuristicLab.Common;
28using HeuristicLab.Core;
29using HeuristicLab.Data;
30using HeuristicLab.Optimization;
31using HeuristicLab.Parameters;
32using HeuristicLab.Persistence.Default.CompositeSerializers.Storable;
33using HeuristicLab.Problems.DataAnalysis;
34using HeuristicLab.Random;
35using MathNet.Numerics.Interpolation;
36using MathNet.Numerics.LinearAlgebra.Double;
37using MathNet.Numerics.LinearAlgebra.Double.Solvers;
38
39
40/*
41 * Experimenting with various spline implementations.
42 *
43 * There are many different variants of using splines.
44 *  - Spline Interpolation
45 *  - Natural Spline Interpolation
46 *  - Smoothing Splines (Reinsch)
47 *  - P-Splines
48 *  - B-Splines
49 *  - Regression Splines
50 *  - Penalized Regression Splines
51 *  -
52 *
53 *
54 *  Generally, a spline is a piecewise function whereby each piece connects two knots.
55 *  In almost all cases splines are univariate. There are extensions to more dimensions
56 *  but this is much less convenient.
57 *  A linear spline is a piecewise linear function. Usually the following constaints are enforced
58 *  for each knot point k and the two piecewise functions fl and fr
59 *  fl(k) = fr(k), fl(k)' = fr(k)', fl(k)''=0, fr(k)''=0, additionally constraints for
60 *  the boundary knots are often enforced.
61
62 *  Spline Interpolation solves the problem of finding a smooth interpolating function
63 *  which goes through all knots exactly.
64 *
65 *  In spline smoothing one tries to find a function minimizing the squared residuals
66 *  and maximizing smoothness. For the smoothing spline g the optimization objective is
67 *  || y - g(x) ||² + \lambda integrate (xmin : xmax) g(x)'' dx.
68 *  The roughness penality is the integral over the second derivative of the smoothing spline
69 *  over the whole domain of x. Lambda is used to balance the amount of smoothing.
70 *  If lambda is set to zero, the solution is an interpolating spline. If lambda is set to
71 *  infinity the solution is a maximally smooth function (= line).
72 *  Interestingly, the estimation of a smoothing spline is a linear operation. As a consequence
73 *  a lot of theory from linear models can be reused (e.g. cross-validation and generalized
74 *  cross-validation) with minimal additional effort. GCV can also be used to tune the
75 *  smoothing parameter automatically. Tuned algorithms are available to solve smoothing
76 *  spline estimation efficiently O(n) or O(m²n) where m is the number of knots and n the
77 *  number of data points.
78 *  It is possible to calculate confidence intervals for the smoothing spline and given
79 *  an error model subsequently also confidence intervals for predictions.
80 *
81 *  We are primarily interested in spline smoothing / regression.
82 *  Alglib and Math.NET provide several different implementations for interpolation
83 *  including polynomial (cubic), Hermite, Akima, Linear, Catmull-Rom, Rational, ...
84 *
85 *  For spline smoothing or regression there are also many differen variants.
86 *
87 *  Smoothing Splines by Reinsch (Reinsch, Smoothing by Spline Functions, 1967).
88 *  Pseudo-code for the algorithm is given in the paper. A small change in the algorithm
89 *  to improve convergence is discussed in (Reinsch, Smoothing by Spline Functions II, 1970)
90 *  The algorithm uses the properties of the smoothing matrix (bandedness) to calculate the
91 *  smoothing spline in O(n). Two parameters have to be set (S and rho). If the noise in the
92 *  data is known then both parameter can be set directly, otherwise S should be set to n-1
93 *  and rho should be found via cross-validation.
94 *  It is not easy to calculate a the CV (PRESS) or GCV in the approach by Reinsch and
95 *  therefore the algorithm is not ideal.
96 *  Wahba and colleagues (1990) introduced CV and GCV estimation.
97 *  Hutchinson and colleages improved the algorithm (Fortran implementation available as
98 *  algorithm 642.f from ACM).
99 *
100 *  CUBGCV (Algorithm 642 Hutchingson)
101 *  Several improvments over the Reinsch algorithm. In particular efficient calculation
102 *  of the
103 
104 *  Finbarr O'Sullivan BART
105 *  ...
106 *
107 *  B-Splines
108 *  ...
109 *
110 *  P-Splines (Eilers and Marx), use the B-Spline basis in combination with a penality
111 *  for large variations of subsequent B-Spline coefficients. Conceptually, if the
112 *  coefficients of neighbouring B-Splines are almost the same then the resulting smoothing
113 *  spline is less wiggly. It has been critizized that this is rather subjective.
114 *  P-Splines also allow calculation of CV and GCV for a given smoothing parameter which
115 *  can be used for optimization of lambda.
116 *  It is easy to use P-Splines within GAMs (no backfitting necessary).
117 *  The paper contains MATLAB code also for the calculation of the B-Spline basis.
118 *  The computational effort is O(?). The authors argue that it is efficient.
119 *  The number of knots is not critical it is possible to use every data point as knot.
120 *
121 *  Alglib: penalized regression spline
122 *  ...
123 *
124 *  Alglib: penalized regression spline with support for weights
125 *  ...
126 */
127
128namespace HeuristicLab.Algorithms.DataAnalysis.Experimental {
129  [Item("Splines", "")]
130  [Creatable(CreatableAttribute.Categories.DataAnalysisRegression, Priority = 102)]
131  [StorableClass]
132  public sealed class Splines : FixedDataAnalysisAlgorithm<IRegressionProblem> {
133    [StorableConstructor]
134    private Splines(bool deserializing) : base(deserializing) { }
135    [StorableHook(HookType.AfterDeserialization)]
136    private void AfterDeserialization() {
137    }
138
139    private Splines(Splines original, Cloner cloner)
140      : base(original, cloner) {
141    }
142    public override IDeepCloneable Clone(Cloner cloner) {
143      return new Splines(this, cloner);
144    }
145
146    public Splines()
147      : base() {
148      var validTypes = new ItemSet<StringValue>(
149        new[] {
150          "Smoothing Spline (Reinsch)",
151          "Smoothing Spline (Reinsch with automatic tolerance determination using LOOCV)",
152          "B-Spline Smoothing",
153          "Penalized Regression Spline (alglib)",
154          "CUBGCV",
155          "Finnbar O'Sullivan (sbart)"
156      }.Select(s => new StringValue(s)));
157
158      Parameters.Add(new ConstrainedValueParameter<StringValue>("Type", "The type of spline (as supported by alglib)", validTypes, validTypes.Last()));
159      Parameters.Add(new ValueParameter<DoubleValue>("Lambda", "Regularization parameter for smoothing splines (0..+inf)", new DoubleValue(100)));
160      Parameters.Add(new ValueParameter<DoubleValue>("StdDev (noise)", "Known error in y values. Used only be Reinsch Smoothing Splines", new DoubleValue(0.1)));
161      Problem = new RegressionProblem();
162    }
163
164
165    protected override void Run(CancellationToken cancellationToken) {
166      double[,] inputMatrix = Problem.ProblemData.Dataset.ToArray(Problem.ProblemData.AllowedInputVariables.Concat(new string[] { Problem.ProblemData.TargetVariable }), Problem.ProblemData.TrainingIndices);
167      if (inputMatrix.Cast<double>().Any(x => double.IsNaN(x) || double.IsInfinity(x)))
168        throw new NotSupportedException("Splines does not support NaN or infinity values in the input dataset.");
169
170      var inputVars = Problem.ProblemData.AllowedInputVariables.ToArray();
171      if (inputVars.Length > 3) throw new ArgumentException();
172
173      var y = Problem.ProblemData.TargetVariableTrainingValues.ToArray();
174      if (inputVars.Length == 1) {
175
176        var x = Problem.ProblemData.Dataset.GetDoubleValues(inputVars.First(), Problem.ProblemData.TrainingIndices).ToArray();
177        alglib.spline1dinterpolant c;
178        var type = ((StringValue)Parameters["Type"].ActualValue).Value;
179        switch (type) {
180          case "Smoothing Spline (Reinsch)": {
181              CalculateSmoothingSplineReinsch(x, y, inputVars);
182              break;
183            }
184          case "Smoothing Spline (Reinsch with automatic tolerance determination using LOOCV)": {
185              double optTol, looRMSE;
186              var model = CalculateSmoothingSplineReinsch(x, y, inputVars, Problem.ProblemData.TargetVariable, out optTol, out looRMSE);
187              Results.Add(new Result("Solution", model.CreateRegressionSolution((IRegressionProblemData)Problem.ProblemData.Clone())));
188              Results.Add(new Result("Optimal tolerance", new DoubleValue(optTol)));
189              Results.Add(new Result("RMSE (LOO-CV)", new DoubleValue(looRMSE)));
190              break;
191            }
192          case "CUBGCV": {
193              CubicSplineGCV.CubGcvReport report;
194              var model =
195                CubicSplineGCV.CalculateCubicSpline(x, y,
196                Problem.ProblemData.TargetVariable, inputVars, out report);
197              var targetVar = Problem.ProblemData.TargetVariable;
198              var problemData = (IRegressionProblemData)Problem.ProblemData.Clone();
199              Results.Add(new Result("Solution", model.CreateRegressionSolution(problemData)));
200              Results.Add(new Result("GCV", new DoubleValue(report.generalizedCrossValidation)));
201              Results.Add(new Result("Estimated error variance", new DoubleValue(report.estimatedErrorVariance)));
202              Results.Add(new Result("Estimated RSS degrees of freedom", new DoubleValue(report.estimatedRSSDegreesOfFreedom)));
203              Results.Add(new Result("Estimated true mean squared error at data points", new DoubleValue(report.estimatedTrueMeanSquaredErrorAtDataPoints)));
204              Results.Add(new Result("Mean square of DF", new DoubleValue(report.meanSquareOfDf)));
205              Results.Add(new Result("Mean square residual", new DoubleValue(report.meanSquareResiudal)));
206              Results.Add(new Result("Smoothing parameter (optimized for GCV)", new DoubleValue(report.smoothingParameter)));
207              break;
208            }
209          case "Finnbar O'Sullivan (sbart)": {
210              SBART.SBART_Report report;
211           
212              var model =
213                SBART.CalculateSBART(x, y,
214                Problem.ProblemData.TargetVariable, inputVars, out report);
215              var targetVar = Problem.ProblemData.TargetVariable;
216              var problemData = (IRegressionProblemData)Problem.ProblemData.Clone();
217              Results.Add(new Result("Solution", model.CreateRegressionSolution(problemData)));
218              Results.Add(new Result("GCV", new DoubleValue(report.gcv)));
219              Results.Add(new Result("Smoothing parameter (optimized for GCV)", new DoubleValue(report.smoothingParameter)));
220              break;
221
222            }
223          case "Penalized Regression Spline (alglib)": {
224              var lambda = ((IValueParameter<DoubleValue>)Parameters["Lambda"]).Value.Value; // controls extent of smoothing
225              var model = CalculatePenalizedRegressionSpline(x, y, lambda, x.Length, Problem.ProblemData.TargetVariable, inputVars);
226              var targetVar = Problem.ProblemData.TargetVariable;
227              var problemData = (IRegressionProblemData)Problem.ProblemData.Clone();
228              Results.Add(new Result("Solution", model.CreateRegressionSolution(problemData)));
229
230              break;
231            }
232          default: throw new NotSupportedException();
233        }
234
235      }
236    }
237
238
239    public static IRegressionModel CalculatePenalizedRegressionSpline(double[] x, double[] y, double lambda, int numKnots, string targetVar, string[] inputVars) {
240      int info;
241      alglib.spline1dinterpolant interpolant;
242      alglib.spline1dfitreport rep;
243
244      alglib.spline1dfitpenalized(x, y, numKnots, lambda, out info, out interpolant, out rep);
245      return new Spline1dModel(interpolant, targetVar, inputVars);
246    }
247
248    public static IRegressionEnsembleModel CalculatePenalizedRegressionSpline(double[] x, double[] y, double lambda, string targetVar, string[] inputVars,
249      out double avgTrainRMSE,
250      out double cvRMSE, out double[] predictions) {
251      int info;
252      alglib.spline1dinterpolant interpolant;
253      alglib.spline1dfitreport rep;
254      int n = x.Length;
255      n = (int)Math.Max(50, 3 * Math.Sqrt(n));
256
257      int folds = 10;
258      int foldSize = x.Length <= 30 ? 1 : x.Length / folds;
259      if (foldSize == 1) folds = x.Length;
260
261      double[] x_loo = new double[(folds - 1) * foldSize];
262      double[] y_loo = new double[(folds - 1) * foldSize];
263      double[] w_loo = Enumerable.Repeat(1.0, x_loo.Length).ToArray();
264
265      var models = new List<IRegressionModel>();
266      var sumTrainSE = 0.0;
267      int nTrain = 0;
268      var sumTestSE = 0.0;
269      int nTest = 0;
270
271      predictions = new double[y.Length];
272
273      for (int i = 0; i < folds; i++) {
274        if (i < folds - 1) {
275          //Console.WriteLine("Copy from {0} to {1} n = {2}", (i + 1) * foldSize, i * foldSize, (folds - i - 1) * foldSize);
276          Array.Copy(x, (i + 1) * foldSize, x_loo, i * foldSize, (folds - i - 1) * foldSize);
277          Array.Copy(y, (i + 1) * foldSize, y_loo, i * foldSize, (folds - i - 1) * foldSize);
278        }
279        alglib.spline1dfitpenalizedw(x_loo, y_loo, w_loo, n, lambda, out info, out interpolant, out rep);
280        //Debug.Assert(y_loo.All(yi => yi != 0.0));
281        //Debug.Assert(x_loo.All(xi => xi != 0.0));
282
283        // training error
284        for (int j = 0; j < x_loo.Length; j++) {
285          var y_pred = alglib.spline1dcalc(interpolant, x[j]);
286          double res = y[j] - y_pred;
287          sumTrainSE += res * res;
288          nTrain++;
289        }
290        // test
291        for (int j = i * foldSize; j < (i + 1) * foldSize; j++) {
292          predictions[j] = alglib.spline1dcalc(interpolant, x[j]);
293          double res = y[j] - predictions[j];
294          sumTestSE += res * res;
295          nTest++;
296        }
297
298        models.Add(new Spline1dModel(interpolant, targetVar, inputVars));
299
300        if (i < folds - 1) {
301          //Console.WriteLine("Copy from {0} to {1} n = {2}", i * foldSize, i * foldSize, foldSize);
302          // add current fold for next iteration
303          Array.Copy(x, i * foldSize, x_loo, i * foldSize, foldSize);
304          Array.Copy(y, i * foldSize, y_loo, i * foldSize, foldSize);
305        }
306      }
307
308      cvRMSE = Math.Sqrt(sumTestSE / nTest);
309      avgTrainRMSE = Math.Sqrt(sumTrainSE / nTrain);
310
311      return new RegressionEnsembleModel(models);
312    }
313
314    private struct SplineData {
315      public double a, b, c, d, x, y;
316    };
317
318
319    /*
320     * The procedure Quincunx, takes as arguments
321       the vectors u, v and w which are respectively the diagonal, the first supradiagonal
322       and the second supradiagonal of the banded matrix. The vector on
323       the LHS of the equation (80) is placed in q which contains the solution on the
324       completion of the procedure.
325    */
326    private void Quincunx(int n, MyArray<double> u, MyArray<double> v, MyArray<double> w, MyArray<double> q) {
327      // { factorisation}
328      u[-1] = 0;
329      u[0] = 0;
330      for (int j = 1; j <= n - 1; j++) {
331        u[j] = u[j] - u[j - 2] * Math.Pow(w[j - 2], 2) - u[j - 1] * Math.Pow(v[j - 1], 2);
332        v[j] = (v[j] - u[j - 1] * v[j - 1] * w[j - 1]) / u[j];
333        w[j] = w[j] / u[j];
334      }
335      // { forward substitution}
336      for (int j = 1; j <= n - 1; j++) {
337        q[j] = q[j] - v[j - 1] * q[j - 1] - w[j - 2] * q[j - 2];
338      }
339      for (int j = 1; j < n - 1; j++) {
340        q[j] = q[j] / u[j];
341      }
342      // { back substitution}
343      q[n + 1] = 0;
344      q[n] = 0;
345      for (int j = n - 1; j >= 1; j--) {
346        q[j] = q[j] - v[j] * q[j + 1] - w[j] * q[j + 2];
347      }
348    }
349
350    public void CalculateSmoothingSplineReinsch(double[] x, double[] y, string[] inputVars) {
351      double s = x.Length + 1;
352      double stdTol = ((IValueParameter<DoubleValue>)Parameters["StdDev (noise)"]).Value.Value; // controls extent of smoothing
353
354      var model = CalculateSmoothingSplineReinsch(x, y, inputVars, stdTol, s, Problem.ProblemData.TargetVariable);
355
356      Results.Add(new Result("Solution", model.CreateRegressionSolution((IRegressionProblemData)Problem.ProblemData.Clone())));
357    }
358
359
360    public static IRegressionEnsembleModel CalculateSmoothingSplineReinsch(double[] x, double[] y, double tol, string targetVar, string[] inputVars,
361      out double avgTrainRMSE,
362      out double cvRMSE) {
363      int info;
364
365      int folds = 10;
366      int foldSize = x.Length <= 30 ? 1 : x.Length / folds;
367      if (foldSize == 1) folds = x.Length;
368
369      double[] x_loo = new double[(folds - 1) * foldSize];
370      double[] y_loo = new double[(folds - 1) * foldSize];
371
372      var models = new List<IRegressionModel>();
373      var sumTrainSE = 0.0;
374      int nTrain = 0;
375      var sumTestSE = 0.0;
376      int nTest = 0;
377
378
379      for (int i = 0; i < folds; i++) {
380        if (i < folds - 1) {
381          Array.Copy(x, (i + 1) * foldSize, x_loo, i * foldSize, x.Length - ((i + 1) * foldSize) - 1);
382          Array.Copy(y, (i + 1) * foldSize, y_loo, i * foldSize, y.Length - ((i + 1) * foldSize) - 1);
383        }
384        var model = CalculateSmoothingSplineReinsch(x_loo, y_loo, inputVars, tol, ((folds - 1) * foldSize) - 1, targetVar);
385
386        // training error
387        for(int j = 0;j<x_loo.Length;j++) {
388          var y_pred = model.GetEstimatedValue(x[j]);
389          double res = y[j] - y_pred;
390          sumTrainSE += res * res;
391          nTrain++;
392        }
393        // test
394        for (int j = i * foldSize; j < (i + 1) * foldSize; j++) {
395          var y_pred = model.GetEstimatedValue(x[j]);
396          double res = y[j] - y_pred;
397          sumTestSE += res * res;
398          nTest++;
399        }
400
401        models.Add(model);
402
403        if (i < folds - 1) {
404          // add current fold for next iteration
405          Array.Copy(x, i * foldSize, x_loo, i * foldSize, foldSize);
406          Array.Copy(y, i * foldSize, y_loo, i * foldSize, foldSize);
407        }
408      }
409
410      cvRMSE = Math.Sqrt(sumTestSE / nTest);
411      avgTrainRMSE = Math.Sqrt(sumTrainSE / nTrain);
412
413      return new RegressionEnsembleModel(models);
414    }
415
416
417    // automatically determines tolerance using cv, SLOW!
418    public static IRegressionModel CalculateSmoothingSplineReinsch(double[] x, double[] y, string[] inputVars, string targetVar,
419      out double optTolerance, out double cvRMSE) {
420      var maxTolerance = y.StandardDeviation();
421      var minTolerance = maxTolerance * 1E-6;
422      optTolerance = maxTolerance;
423      var tol = maxTolerance;
424      cvRMSE = double.PositiveInfinity;
425      IRegressionModel bestModel = new ConstantModel(y.Average(), targetVar);
426      while (tol > minTolerance) {
427        double curCVRMSE;
428        double avgTrainRmse;
429        var model = Splines.CalculateSmoothingSplineReinsch(
430          x, y, tol, targetVar, inputVars, out avgTrainRmse, out curCVRMSE);
431
432        if (curCVRMSE < cvRMSE) {
433          cvRMSE = curCVRMSE;
434          optTolerance = tol;
435          bestModel = model;
436        }
437        tol *= 0.85;
438      }
439      return bestModel;
440    }
441
442    public static ReinschSmoothingSplineModel CalculateSmoothingSplineReinsch(double[] xOrig, double[] yOrig, string[] inputVars, double stdTol, double s, string targetVar) {
443      var minX = xOrig.Min();
444      var maxX = xOrig.Max();
445      var range = maxX - minX;
446
447      double[] w = Enumerable.Repeat(stdTol, xOrig.Length).ToArray();
448      SortAndBin((double[])xOrig.Clone(), (double[])yOrig.Clone(), w, out xOrig, out yOrig, out w, scaling: false);
449
450      // See Smoothing by Spline Functions, Reinsch, 1967
451      // move x and y into an array indexed with 1..n to match indexing in Reinsch paper
452      int n = xOrig.Length;
453      var x = new MyArray<double>(1, xOrig);
454      var y = new MyArray<double>(1, yOrig);
455      var inv_dy = new MyArray<double>(1, w);
456
457      int n1 = 1;
458      int n2 = n;
459
460      // results
461      var a = new MyArray<double>(n1, n);
462      var b = new MyArray<double>(n1, n);
463      var c = new MyArray<double>(n1, n);
464      var d = new MyArray<double>(n1, n);
465
466      // smooth
467      {
468        int i, m1, m2; double e, f, f2, g = 0.0, h, p;
469        MyArray<double>
470          r = new MyArray<double>(n1 - 1, n + 2),
471          r1 = new MyArray<double>(n1 - 1, n + 2),
472          r2 = new MyArray<double>(n1 - 1, n + 2),
473          t = new MyArray<double>(n1 - 1, n + 2),
474          t1 = new MyArray<double>(n1 - 1, n + 2),
475          u = new MyArray<double>(n1 - 1, n + 2),
476          v = new MyArray<double>(n1 - 1, n + 2);
477        m1 = n1 - 1;
478        m2 = n2 + 1;
479        r[m1] = r[n1] = r1[n2] = r2[n2] = r2[m2] =
480          u[m1] = u[n1] = u[n2] = u[m2] = p = 0.0;
481        m1 = n1 + 1;
482        m2 = n2 - 1;
483        h = x[m1] - x[n1]; f = (y[m1] - y[n1]) / h;
484        for (i = m1; i <= m2; i++) {
485          g = h; h = x[i + 1] - x[i];
486          e = f; f = (y[i + 1] - y[i]) / h;
487          a[i] = f - e; t[i] = 2 * (g + h) / 3; t1[i] = h / 3;
488          r2[i] = inv_dy[i - 1] / g; r[i] = inv_dy[i + 1] / h;
489          r1[i] = -inv_dy[i] / g - inv_dy[i] / h;
490        }
491        for (i = m1; i <= m2; i++) {
492          b[i] = r[i] * r[i] + r1[i] * r1[i] + r2[i] * r2[i];
493          c[i] = r[i] * r1[i + 1] + r1[i] * r2[i + 1];
494          d[i] = r[i] * r2[i + 2];
495        }
496        f2 = -s;
497        next:
498        for (i = m1; i <= m2; i++) {
499          r1[i - 1] = f * r[i - 1]; r2[i - 2] = g * r[i - 2];
500          r[i] = 1 / (p * b[i] + t[i] - f * r1[i - 1] - g * r2[i - 2]);
501          u[i] = a[i] - r1[i - 1] * u[i - 1] - r2[i - 2] * u[i - 2];
502          f = p * c[i] + t1[i] - h * r1[i - 1]; g = h; h = d[i] * p;
503        }
504        for (i = m2; i >= m1; i--) {
505          u[i] = r[i] * u[i] - r1[i] * u[i + 1] - r2[i] * u[i + 2];
506        }
507        e = h = 0;
508        for (i = n1; i <= m2; i++) {
509          g = h; h = (u[i + 1] - u[i]) / (x[i + 1] - x[i]);
510          v[i] = (h - g) * inv_dy[i] * inv_dy[i]; e = e + v[i] * (h - g);
511        }
512        g = v[n2] = -h * inv_dy[n2] * inv_dy[n2]; e = e - g * h;
513        g = f2; f2 = e * p * p;
514        if (f2 >= s || f2 <= g) goto fin;
515        f = 0; h = (v[m1] - v[n1]) / (x[m1] - x[n1]);
516        for (i = m1; i <= m2; i++) {
517          g = h; h = (v[i + 1] - v[i]) / (x[i + 1] - x[i]);
518          g = h - g - r1[i - 1] * r[i - 1] - r2[i - 2] * r[i - 2];
519          f = f + g * r[i] * g; r[i] = g;
520        }
521        h = e - p * f; if (h <= 0) goto fin;
522        p = p + (s - f2) / ((Math.Sqrt(s / e) + p) * h); goto next;
523
524        fin:
525        for (i = n1; i <= n2; i++) {
526          a[i] = y[i] - p * v[i];
527          c[i] = u[i];
528        }
529        for (i = n1; i <= m2; i++) {
530          h = x[i + 1] - x[i];
531          d[i] = (c[i + 1] - c[i]) / (3 * h);
532          b[i] = (a[i + 1] - a[i]) / h - (h * d[i] + c[i]) * h;
533        }
534      }
535
536      // extrapolate for xx > x[n2]
537      b[b.Length] = b[b.Length - 1];
538      d[1] = 0;
539
540      return new ReinschSmoothingSplineModel(a, b, c, d, x, targetVar, inputVars);
541    }
542
543    public static void SortAndBin(double[] x, double[] y, double[] w, out double[] x2, out double[] y2, out double[] w2, bool scaling = false) {
544      var sortedIdx = Enumerable.Range(0, x.Length).ToArray();
545      // sort by x
546      Array.Sort(x, sortedIdx);
547
548      var xl = new List<double>();
549      var yl = new List<double>();
550      var wl = new List<double>();
551
552      int n = x.Length;
553      var range = x[n - 1] - x[0];
554      var binSize = range / n / 10;
555      {
556        // binning
557        int i = 0;
558        while (i < n) {
559          int k = 0;
560          int j = i;
561          double sumX = 0.0;
562          double sumY = 0.0;
563          double sumW = 0.0;
564          while (j < n && x[j] - x[i] <= binSize) {
565            k++;
566            sumX += x[j];
567            sumY += y[sortedIdx[j]];
568            sumW += w[sortedIdx[j]];
569            j++;
570          }
571          var avgX = sumX / k;
572          if (scaling) avgX = (avgX - x[0]) / range;
573          xl.Add(avgX);
574          yl.Add(sumY / k);
575          wl.Add(sumW);
576          i += k;
577        }
578      }
579
580      x2 = xl.ToArray();
581      y2 = yl.ToArray();
582      w2 = wl.ToArray();
583    }
584  }
585
586
587  // array with non-zero lower index
588  internal class MyArray<T> {
589    private T[] arr;
590    private int lowerBound;
591
592    public int Length { get { return arr.Length; } }
593
594    public T this[int key] {
595      get {
596        return arr[key - lowerBound];
597      }
598      set {
599        arr[key - lowerBound] = value;
600      }
601    }
602
603    public MyArray(int lowerBound, int numElements) {
604      this.lowerBound = lowerBound;
605      arr = new T[numElements];
606    }
607    public MyArray(int lowerBound, T[] source) : this(lowerBound, source.Length) {
608      Array.Copy(source, arr, source.Length);
609    }
610
611    public T[] ToArray() {
612      var res = new T[arr.Length];
613      Array.Copy(arr, res, res.Length);
614      return res;
615    }
616  }
617
618
619  // UNFINISHED
620  public class ReinschSmoothingSplineModel : NamedItem, IRegressionModel {
621    private MyArray<double> a;
622    private MyArray<double> b;
623    private MyArray<double> c;
624    private MyArray<double> d;
625    private MyArray<double> x;
626    private double offset;
627    private double scale;
628
629    public string TargetVariable { get; set; }
630
631    public IEnumerable<string> VariablesUsedForPrediction { get; private set; }
632
633    public event EventHandler TargetVariableChanged;
634
635    public ReinschSmoothingSplineModel(ReinschSmoothingSplineModel orig, Cloner cloner) : base(orig, cloner) {
636      this.TargetVariable = orig.TargetVariable;
637      this.VariablesUsedForPrediction = orig.VariablesUsedForPrediction.ToArray();
638      this.a = orig.a;
639      this.b = orig.b;
640      this.c = orig.c;
641      this.d = orig.d;
642      this.x = orig.x;
643      this.scale = orig.scale;
644      this.offset = orig.offset;
645    }
646    internal ReinschSmoothingSplineModel(
647      MyArray<double> a,
648      MyArray<double> b,
649      MyArray<double> c,
650      MyArray<double> d,
651      MyArray<double> x,
652      string targetVar, string[] inputs, double offset = 0, double scale = 1) : base("SplineModel", "SplineModel") {
653      this.a = a;
654      this.b = b;
655      this.c = c;
656      this.d = d;
657      this.x = x;
658      this.TargetVariable = targetVar;
659      this.VariablesUsedForPrediction = inputs;
660      this.scale = scale;
661      this.offset = offset;
662    }
663
664    public override IDeepCloneable Clone(Cloner cloner) {
665      return new ReinschSmoothingSplineModel(this, cloner);
666    }
667
668    public IRegressionSolution CreateRegressionSolution(IRegressionProblemData problemData) {
669      return new RegressionSolution(this, (IRegressionProblemData)problemData.Clone());
670    }
671
672    public double GetEstimatedValue(double xx) {
673      xx = (xx - offset) * scale;
674      int n = a.Length;
675      if (xx <= x[1]) {
676        double h = xx - x[1];
677        return a[1] + h * (b[1] + h * (c[1] + h * d[1]));
678      } else if (xx >= x[n]) {
679        double h = xx - x[n];
680        return a[n] + h * (b[n] + h * (c[n] + h * d[n]));
681      } else {
682        // binary search
683        int lower = 1;
684        int upper = n;
685        while (true) {
686          if (upper < lower) throw new InvalidProgramException();
687          int i = lower + (upper - lower) / 2;
688          if (x[i] <= xx && xx < x[i + 1]) {
689            double h = xx - x[i];
690            return a[i] + h * (b[i] + h * (c[i] + h * d[i]));
691          } else if (xx < x[i]) {
692            upper = i - 1;
693          } else {
694            lower = i + 1;
695          }
696        }
697      }
698    }
699
700    public IEnumerable<double> GetEstimatedValues(IDataset dataset, IEnumerable<int> rows) {
701      int n = x.Length;
702      var product = dataset.GetDoubleValues(VariablesUsedForPrediction.First(), rows).ToArray();
703      foreach (var factor in VariablesUsedForPrediction.Skip(1)) {
704        product = product.Zip(dataset.GetDoubleValues(factor, rows), (pi, fi) => pi * fi).ToArray();
705      }
706      foreach (var xx in product) {
707        yield return GetEstimatedValue(xx);
708      }
709    }
710  }
711
712  // UNFINISHED
713  public class Spline1dModel : NamedItem, IRegressionModel {
714    alglib.spline1dinterpolant interpolant;
715
716    public string TargetVariable { get; set; }
717
718    public IEnumerable<string> VariablesUsedForPrediction { get; private set; }
719
720    public event EventHandler TargetVariableChanged;
721
722    public Spline1dModel(Spline1dModel orig, Cloner cloner) : base(orig, cloner) {
723      this.TargetVariable = orig.TargetVariable;
724      this.VariablesUsedForPrediction = orig.VariablesUsedForPrediction.ToArray();
725
726      this.interpolant = (alglib.spline1dinterpolant)orig.interpolant.make_copy();
727    }
728    public Spline1dModel(
729      alglib.spline1dinterpolant interpolant,
730      string targetVar, string[] inputs, double offset = 0, double scale = 1) : base("SplineModel", "SplineModel") {
731      this.interpolant = (alglib.spline1dinterpolant)interpolant.make_copy();
732      this.TargetVariable = targetVar;
733      this.VariablesUsedForPrediction = inputs;
734    }
735
736    public override IDeepCloneable Clone(Cloner cloner) {
737      return new Spline1dModel(this, cloner);
738    }
739
740    public IRegressionSolution CreateRegressionSolution(IRegressionProblemData problemData) {
741      return new RegressionSolution(this, (IRegressionProblemData)problemData.Clone());
742    }
743
744    public double GetEstimatedValue(double xx) {
745      return alglib.spline1dcalc(interpolant, xx);
746    }
747
748    public IEnumerable<double> GetEstimatedValues(IDataset dataset, IEnumerable<int> rows) {
749      var product = dataset.GetDoubleValues(VariablesUsedForPrediction.First(), rows).ToArray();
750      foreach (var factor in VariablesUsedForPrediction.Skip(1)) {
751        product = product.Zip(dataset.GetDoubleValues(factor, rows), (pi, fi) => pi * fi).ToArray();
752      }
753      foreach (var xx in product) {
754        yield return GetEstimatedValue(xx);
755      }
756    }
757  }
758}
Note: See TracBrowser for help on using the repository browser.