Free cookie consent management tool by TermsFeed Policy Generator

source: branches/2994-AutoDiffForIntervals/HeuristicLab.Problems.DataAnalysis.Regression.Symbolic.Extensions/ConstrainedNLSInternal.cs @ 17200

Last change on this file since 17200 was 17200, checked in by gkronber, 5 years ago

#2994 implemented checking of autodiff with numeric differentials and fixed a bug in IntervalEvaluator.

File size: 22.3 KB
Line 
1using System;
2using System.Collections.Generic;
3using System.Diagnostics;
4using System.Linq;
5using System.Runtime.InteropServices;
6using HeuristicLab.Common;
7using HeuristicLab.Encodings.SymbolicExpressionTreeEncoding;
8
9namespace HeuristicLab.Problems.DataAnalysis.Symbolic.Regression.Extensions {
10  internal class ConstrainedNLSInternal {
11    private readonly int maxIterations;
12    public int MaxIterations => maxIterations;
13
14    private readonly string solver;
15    public string Solver => solver;
16
17    private readonly ISymbolicExpressionTree expr;
18    public ISymbolicExpressionTree Expr => expr;
19
20    private readonly IRegressionProblemData problemData;
21
22    public IRegressionProblemData ProblemData => problemData;
23
24
25    public event Action FunctionEvaluated;
26    public event Action<int, double> ConstraintEvaluated;
27
28    private double bestError = double.MaxValue;
29    public double BestError => bestError;
30
31    private double curError = double.MaxValue;
32    public double CurError => curError;
33
34    private double[] bestSolution;
35    public double[] BestSolution => bestSolution;
36
37    private ISymbolicExpressionTree bestTree;
38    public ISymbolicExpressionTree BestTree => bestTree;
39
40    // begin internal state
41    private IntPtr nlopt;
42    private SymbolicDataAnalysisExpressionTreeLinearInterpreter interpreter;
43    private readonly NLOpt.nlopt_func calculateObjectiveDelegate; // must hold the delegate to prevent GC
44    private readonly IntPtr[] constraintDataPtr; // must hold the objects to prevent GC
45    private readonly NLOpt.nlopt_func[] calculateConstraintDelegates; // must hold the delegates to prevent GC
46    private readonly List<double> thetaValues;
47    private readonly IDictionary<string, Interval> dataIntervals;
48    private readonly int[] trainingRows;
49    private readonly double[] target;
50    private readonly ISymbolicExpressionTree preparedTree;
51    private readonly ISymbolicExpressionTreeNode[] preparedTreeParameterNodes;
52    private readonly List<ConstantTreeNode>[] allThetaNodes;
53    private readonly double[] fi_eval;
54    private readonly double[,] jac_eval;
55    private readonly ISymbolicExpressionTree scaledTree;
56
57    // end internal state
58
59
60    // for data exchange to/from optimizer in native code
61    [StructLayout(LayoutKind.Sequential)]
62    private struct ConstraintData {
63      public int Idx;
64      public ISymbolicExpressionTree Tree;
65      public ISymbolicExpressionTreeNode[] ParameterNodes;
66    }
67
68    internal ConstrainedNLSInternal(string solver, ISymbolicExpressionTree expr, int maxIterations, IRegressionProblemData problemData, double ftol_rel = 0, double ftol_abs = 0, double maxTime = 0) {
69      this.solver = solver;
70      this.expr = expr;
71      this.maxIterations = maxIterations;
72      this.problemData = problemData;
73      this.interpreter = new SymbolicDataAnalysisExpressionTreeLinearInterpreter();
74
75
76      var intervalConstraints = problemData.IntervalConstraints;
77      dataIntervals = problemData.VariableRanges.GetIntervals();
78      trainingRows = problemData.TrainingIndices.ToArray();
79      // buffers
80      target = problemData.TargetVariableTrainingValues.ToArray();
81      var targetStDev = target.StandardDeviationPop();
82      var targetVariance = targetStDev * targetStDev;
83      var targetMean = target.Average();
84      var pred = interpreter.GetSymbolicExpressionTreeValues(expr, problemData.Dataset, trainingRows).ToArray();
85
86      if (pred.Any(pi => double.IsInfinity(pi) || double.IsNaN(pi))) throw new ArgumentException("The expression produces NaN or infinite values.");
87
88      #region linear scaling
89      var predStDev = pred.StandardDeviationPop();
90      if (predStDev == 0) throw new ArgumentException("The expression is constant.");
91      var predMean = pred.Average();
92
93      var scalingFactor = targetStDev / predStDev;
94      var offset = targetMean - predMean * scalingFactor;
95
96      scaledTree = CopyAndScaleTree(expr, scalingFactor, offset);
97      #endregion
98
99      // convert constants to variables named theta...
100      var treeForDerivation = ReplaceConstWithVar(scaledTree, out List<string> thetaNames, out thetaValues); // copies the tree
101
102      // create trees for relevant derivatives
103      Dictionary<string, ISymbolicExpressionTree> derivatives = new Dictionary<string, ISymbolicExpressionTree>();
104      allThetaNodes = thetaNames.Select(_ => new List<ConstantTreeNode>()).ToArray();
105      var constraintTrees = new List<ISymbolicExpressionTree>();
106      foreach (var constraint in intervalConstraints.Constraints) {
107        if (constraint.IsDerivation) {
108          if (!problemData.AllowedInputVariables.Contains(constraint.Variable))
109            throw new ArgumentException($"Invalid constraint: the variable {constraint.Variable} does not exist in the dataset.");
110          var df = DerivativeCalculator.Derive(treeForDerivation, constraint.Variable);
111
112          // NLOpt requires constraint expressions of the form c(x) <= 0
113          // -> we make two expressions, one for the lower bound and one for the upper bound
114
115          if (constraint.Interval.UpperBound < double.PositiveInfinity) {
116            var df_smaller_upper = Subtract((ISymbolicExpressionTree)df.Clone(), CreateConstant(constraint.Interval.UpperBound));
117            // convert variables named theta back to constants
118            var df_prepared = ReplaceVarWithConst(df_smaller_upper, thetaNames, thetaValues, allThetaNodes);
119            constraintTrees.Add(df_prepared);
120          }
121          if (constraint.Interval.LowerBound > double.NegativeInfinity) {
122            var df_larger_lower = Subtract(CreateConstant(constraint.Interval.LowerBound), (ISymbolicExpressionTree)df.Clone());
123            // convert variables named theta back to constants
124            var df_prepared = ReplaceVarWithConst(df_larger_lower, thetaNames, thetaValues, allThetaNodes);
125            constraintTrees.Add(df_prepared);
126          }
127        } else {
128          if (constraint.Interval.UpperBound < double.PositiveInfinity) {
129            var f_smaller_upper = Subtract((ISymbolicExpressionTree)treeForDerivation.Clone(), CreateConstant(constraint.Interval.UpperBound));
130            // convert variables named theta back to constants
131            var df_prepared = ReplaceVarWithConst(f_smaller_upper, thetaNames, thetaValues, allThetaNodes);
132            constraintTrees.Add(df_prepared);
133          }
134          if (constraint.Interval.LowerBound > double.NegativeInfinity) {
135            var f_larger_lower = Subtract(CreateConstant(constraint.Interval.LowerBound), (ISymbolicExpressionTree)treeForDerivation.Clone());
136            // convert variables named theta back to constants
137            var df_prepared = ReplaceVarWithConst(f_larger_lower, thetaNames, thetaValues, allThetaNodes);
138            constraintTrees.Add(df_prepared);
139          }
140        }
141      }
142
143      preparedTree = ReplaceVarWithConst(treeForDerivation, thetaNames, thetaValues, allThetaNodes);
144      preparedTreeParameterNodes = GetParameterNodes(preparedTree, allThetaNodes);
145
146      var dim = thetaValues.Count;
147      fi_eval = new double[target.Length]; // init buffer;
148      jac_eval = new double[target.Length, dim]; // init buffer
149
150
151      var minVal = Math.Min(-1000.0, thetaValues.Min());
152      var maxVal = Math.Max(1000.0, thetaValues.Max());
153      var lb = Enumerable.Repeat(minVal, thetaValues.Count).ToArray();
154      var up = Enumerable.Repeat(maxVal, thetaValues.Count).ToArray();
155      nlopt = NLOpt.nlopt_create(GetSolver(solver), (uint)dim);
156
157      NLOpt.nlopt_set_lower_bounds(nlopt, lb);
158      NLOpt.nlopt_set_upper_bounds(nlopt, up);
159      calculateObjectiveDelegate = new NLOpt.nlopt_func(CalculateObjective); // keep a reference to the delegate (see below)
160      NLOpt.nlopt_set_min_objective(nlopt, calculateObjectiveDelegate, IntPtr.Zero);
161
162
163      constraintDataPtr = new IntPtr[constraintTrees.Count];
164      calculateConstraintDelegates = new NLOpt.nlopt_func[constraintTrees.Count]; // make sure we keep a reference to the delegates (otherwise GC will free delegate objects see https://stackoverflow.com/questions/7302045/callback-delegates-being-collected#7302258)
165      for (int i = 0; i < constraintTrees.Count; i++) {
166        var constraintData = new ConstraintData() { Idx = i, Tree = constraintTrees[i], ParameterNodes = GetParameterNodes(constraintTrees[i], allThetaNodes) };
167        constraintDataPtr[i] = Marshal.AllocHGlobal(Marshal.SizeOf<ConstraintData>());
168        Marshal.StructureToPtr(constraintData, constraintDataPtr[i], fDeleteOld: false);
169        calculateConstraintDelegates[i] = new NLOpt.nlopt_func(CalculateConstraint);
170        NLOpt.nlopt_add_inequality_constraint(nlopt, calculateConstraintDelegates[i], constraintDataPtr[i], 1e-8);
171      }
172
173      NLOpt.nlopt_set_ftol_rel(nlopt, ftol_rel);
174      NLOpt.nlopt_set_ftol_abs(nlopt, ftol_abs);
175      NLOpt.nlopt_set_maxtime(nlopt, maxTime);
176      NLOpt.nlopt_set_maxeval(nlopt, maxIterations);
177    }
178
179    ~ConstrainedNLSInternal() {
180      if (nlopt != IntPtr.Zero)
181        NLOpt.nlopt_destroy(nlopt);
182      if (constraintDataPtr != null) {
183        for (int i = 0; i < constraintDataPtr.Length; i++)
184          if (constraintDataPtr[i] != IntPtr.Zero)
185            Marshal.FreeHGlobal(constraintDataPtr[i]);
186      }
187    }
188
189
190    internal void Optimize() {
191      var x = thetaValues.ToArray();  /* initial guess */
192      double minf = double.MaxValue; /* minimum objective value upon return */
193      var res = NLOpt.nlopt_optimize(nlopt, x, ref minf);
194      bestSolution = x;
195      bestError = minf;
196
197      if (res < 0) {
198        throw new InvalidOperationException($"NLOpt failed {res} {NLOpt.nlopt_get_errmsg(nlopt)}");
199      } else {
200        // update parameters in tree
201        var pIdx = 0;
202        // here we lose the two last parameters (for linear scaling)
203        foreach (var node in scaledTree.IterateNodesPostfix()) {
204          if (node is ConstantTreeNode constTreeNode) {
205            constTreeNode.Value = x[pIdx++];
206          } else if (node is VariableTreeNode varTreeNode) {
207            varTreeNode.Weight = x[pIdx++];
208          }
209        }
210        if (pIdx != x.Length) throw new InvalidProgramException();
211      }
212      bestTree = scaledTree;
213    }
214
215    double CalculateObjective(uint dim, double[] curX, double[] grad, IntPtr data) {
216      UpdateThetaValues(curX);
217      var sse = 0.0;
218
219      if (grad != null) {
220        var autoDiffEval = new VectorAutoDiffEvaluator();
221        autoDiffEval.Evaluate(preparedTree, problemData.Dataset, trainingRows,
222          preparedTreeParameterNodes, fi_eval, jac_eval);
223
224        // calc sum of squared errors and gradient
225        for (int j = 0; j < grad.Length; j++) grad[j] = 0;
226        for (int i = 0; i < target.Length; i++) {
227          var r = target[i] - fi_eval[i];
228          sse += r * r;
229          for (int j = 0; j < grad.Length; j++) {
230            grad[j] -= 2 * r * jac_eval[i, j];
231          }
232        }
233        // average
234        for (int j = 0; j < grad.Length; j++) { grad[j] /= target.Length; }
235
236        #region check gradient
237        var eval = new VectorEvaluator();
238        for (int i = 0; i < dim; i++) {
239          // make two additional evaluations
240          var xForNumDiff = (double[])curX.Clone();
241          double delta = Math.Abs(xForNumDiff[i] * 1e-5);
242          xForNumDiff[i] += delta;
243          UpdateThetaValues(xForNumDiff);
244          var evalHigh = eval.Evaluate(preparedTree, problemData.Dataset, trainingRows);
245          var mseHigh = MSE(target, evalHigh);
246          xForNumDiff[i] = curX[i] - delta;
247          UpdateThetaValues(xForNumDiff);
248          var evalLow = eval.Evaluate(preparedTree, problemData.Dataset, trainingRows);
249          var mseLow = MSE(target, evalLow);
250
251          var numericDiff = (mseHigh - mseLow) / (2 * delta);
252          var autoDiff = grad[i];
253          if ((Math.Abs(autoDiff) < 1e-10 && Math.Abs(numericDiff) > 1e-2)
254            || (Math.Abs(autoDiff) >= 1e-10 && Math.Abs((numericDiff - autoDiff) / numericDiff) > 1e-2))
255            throw new InvalidProgramException();
256        }
257        #endregion
258      } else {
259        var eval = new VectorEvaluator();
260        var prediction = eval.Evaluate(preparedTree, problemData.Dataset, trainingRows);
261
262        // calc sum of squared errors
263        sse = 0.0;
264        for (int i = 0; i < target.Length; i++) {
265          var r = target[i] - prediction[i];
266          sse += r * r;
267        }
268      }
269
270      UpdateBestSolution(sse / target.Length, curX);
271      RaiseFunctionEvaluated();
272
273      if (double.IsNaN(sse)) return double.MaxValue;
274      return sse / target.Length;
275    }
276
277    private double MSE(double[] a, double[] b) {
278      Trace.Assert(a.Length == b.Length);
279      var sse = 0.0;
280      for (int i = 0; i < a.Length; i++) sse += (a[i] - b[i]) * (a[i] - b[i]);
281      return sse / a.Length;
282    }
283
284    private void UpdateBestSolution(double curF, double[] curX) {
285      if (double.IsNaN(curF) || double.IsInfinity(curF)) return;
286      else if (curF < bestError) {
287        bestError = curF;
288        bestSolution = (double[])curX.Clone();
289      }
290    }
291
292    private void UpdateConstraintViolations(int constraintIdx, double value) {
293      if (double.IsNaN(value) || double.IsInfinity(value)) return;
294      RaiseConstraintEvaluated(constraintIdx, value);
295      // else if (curF < bestError) {
296      //   bestError = curF;
297      //   bestSolution = (double[])curX.Clone();
298      // }
299    }
300
301    double CalculateConstraint(uint dim, double[] curX, double[] grad, IntPtr data) {
302      UpdateThetaValues(curX);
303      var intervalEvaluator = new IntervalEvaluator();
304      var constraintData = Marshal.PtrToStructure<ConstraintData>(data);
305
306      if (grad != null) for (int j = 0; j < grad.Length; j++) grad[j] = 0; // clear grad
307
308      var interval = intervalEvaluator.Evaluate(constraintData.Tree, dataIntervals, constraintData.ParameterNodes,
309        out double[] lowerGradient, out double[] upperGradient);
310
311      // we transformed this to a constraint c(x) <= 0, so only the upper bound is relevant for us
312      if (grad != null) for (int j = 0; j < grad.Length; j++) { grad[j] = upperGradient[j]; }
313
314      #region check gradient
315      for (int i = 0; i < dim; i++) {
316        // make two additional evaluations
317        var xForNumDiff = (double[])curX.Clone();
318        double delta = Math.Abs(xForNumDiff[i] * 1e-5);
319        xForNumDiff[i] += delta;
320        UpdateThetaValues(xForNumDiff);
321        var evalHigh = intervalEvaluator.Evaluate(constraintData.Tree, dataIntervals, constraintData.ParameterNodes,
322          out double[] unusedLowerGradientHigh, out double[] unusedUpperGradientHigh);
323
324        xForNumDiff[i] = curX[i] - delta;
325        UpdateThetaValues(xForNumDiff);
326        var evalLow = intervalEvaluator.Evaluate(constraintData.Tree, dataIntervals, constraintData.ParameterNodes,
327          out double[] unusedLowerGradientLow, out double[] unusedUpperGradientLow);
328       
329        var numericDiff = (evalHigh.UpperBound - evalLow.UpperBound) / (2 * delta);
330        var autoDiff = grad[i];
331        // XXX GRADIENT FOR UPPER BOUND FOR FUNCTION IS ZERO FOR THE FIRST SET OF VARIABLES? GRADIENT FOR ADDITION INCORRECT?!
332        if ((Math.Abs(autoDiff) < 1e-10 && Math.Abs(numericDiff) > 1e-2)
333          || (Math.Abs(autoDiff) >= 1e-10 && Math.Abs((numericDiff - autoDiff) / numericDiff) > 1e-2))
334          throw new InvalidProgramException();
335      }
336      #endregion
337
338
339      UpdateConstraintViolations(constraintData.Idx, interval.UpperBound);
340      if (double.IsNaN(interval.UpperBound)) return double.MaxValue;
341      else return interval.UpperBound;
342    }
343
344
345    void UpdateThetaValues(double[] theta) {
346      for (int i = 0; i < theta.Length; ++i) {
347        foreach (var constNode in allThetaNodes[i]) constNode.Value = theta[i];
348      }
349    }
350
351    internal void RequestStop() {
352      NLOpt.nlopt_set_force_stop(nlopt, 1); // hopefully NLOpt is thread safe  , val must be <> 0 otherwise no effect
353    }
354
355    private void RaiseFunctionEvaluated() {
356      FunctionEvaluated?.Invoke();
357    }
358
359    private void RaiseConstraintEvaluated(int idx, double value) {
360      ConstraintEvaluated?.Invoke(idx, value);
361    }
362
363
364    #region helper
365
366    private static ISymbolicExpressionTree CopyAndScaleTree(ISymbolicExpressionTree tree, double scalingFactor, double offset) {
367      var m = (ISymbolicExpressionTree)tree.Clone();
368
369      var add = MakeNode<Addition>(MakeNode<Multiplication>(m.Root.GetSubtree(0).GetSubtree(0), CreateConstant(scalingFactor)), CreateConstant(offset));
370      m.Root.GetSubtree(0).RemoveSubtree(0);
371      m.Root.GetSubtree(0).AddSubtree(add);
372      return m;
373    }
374
375    private static void UpdateConstants(ISymbolicExpressionTreeNode[] nodes, double[] constants) {
376      if (nodes.Length != constants.Length) throw new InvalidOperationException();
377      for (int i = 0; i < nodes.Length; i++) {
378        if (nodes[i] is VariableTreeNode varNode) varNode.Weight = constants[i];
379        else if (nodes[i] is ConstantTreeNode constNode) constNode.Value = constants[i];
380      }
381    }
382
383    private NLOpt.nlopt_algorithm GetSolver(string solver) {
384      if (solver.Contains("MMA")) return NLOpt.nlopt_algorithm.NLOPT_LD_MMA;
385      if (solver.Contains("COBYLA")) return NLOpt.nlopt_algorithm.NLOPT_LN_COBYLA;
386      if (solver.Contains("CCSAQ")) return NLOpt.nlopt_algorithm.NLOPT_LD_CCSAQ;
387      if (solver.Contains("ISRES")) return NLOpt.nlopt_algorithm.NLOPT_GN_ISRES;
388      throw new ArgumentException($"Unknown solver {solver}");
389    }
390
391    private static ISymbolicExpressionTreeNode[] GetParameterNodes(ISymbolicExpressionTree tree, List<ConstantTreeNode>[] allNodes) {
392      // TODO better solution necessary
393      var treeConstNodes = tree.IterateNodesPostfix().OfType<ConstantTreeNode>().ToArray();
394      var paramNodes = new ISymbolicExpressionTreeNode[allNodes.Length];
395      for (int i = 0; i < paramNodes.Length; i++) {
396        paramNodes[i] = allNodes[i].SingleOrDefault(n => treeConstNodes.Contains(n));
397      }
398      return paramNodes;
399    }
400
401    private static ISymbolicExpressionTree ReplaceVarWithConst(ISymbolicExpressionTree tree, List<string> thetaNames, List<double> thetaValues, List<ConstantTreeNode>[] thetaNodes) {
402      var copy = (ISymbolicExpressionTree)tree.Clone();
403      var nodes = copy.IterateNodesPostfix().ToList();
404      for (int i = 0; i < nodes.Count; i++) {
405        var n = nodes[i] as VariableTreeNode;
406        if (n != null) {
407          var thetaIdx = thetaNames.IndexOf(n.VariableName);
408          if (thetaIdx >= 0) {
409            var parent = n.Parent;
410            if (thetaNodes[thetaIdx].Any()) {
411              // HACK: REUSE CONSTANT TREE NODE IN SEVERAL TREES
412              // we use this trick to allow autodiff over thetas when thetas occurr multiple times in the tree (e.g. in derived trees)
413              var constNode = thetaNodes[thetaIdx].First();
414              var childIdx = parent.IndexOfSubtree(n);
415              parent.RemoveSubtree(childIdx);
416              parent.InsertSubtree(childIdx, constNode);
417            } else {
418              var constNode = (ConstantTreeNode)CreateConstant(thetaValues[thetaIdx]);
419              var childIdx = parent.IndexOfSubtree(n);
420              parent.RemoveSubtree(childIdx);
421              parent.InsertSubtree(childIdx, constNode);
422              thetaNodes[thetaIdx].Add(constNode);
423            }
424          }
425        }
426      }
427      return copy;
428    }
429
430    private static ISymbolicExpressionTree ReplaceConstWithVar(ISymbolicExpressionTree tree, out List<string> thetaNames, out List<double> thetaValues) {
431      thetaNames = new List<string>();
432      thetaValues = new List<double>();
433      var copy = (ISymbolicExpressionTree)tree.Clone();
434      var nodes = copy.IterateNodesPostfix().ToList();
435
436      int n = 1;
437      for (int i = 0; i < nodes.Count; ++i) {
438        var node = nodes[i];
439        if (node is ConstantTreeNode constantTreeNode) {
440          var thetaVar = (VariableTreeNode)new Problems.DataAnalysis.Symbolic.Variable().CreateTreeNode();
441          thetaVar.Weight = 1;
442          thetaVar.VariableName = $"θ{n++}";
443
444          thetaNames.Add(thetaVar.VariableName);
445          thetaValues.Add(constantTreeNode.Value);
446
447          var parent = constantTreeNode.Parent;
448          if (parent != null) {
449            var index = constantTreeNode.Parent.IndexOfSubtree(constantTreeNode);
450            parent.RemoveSubtree(index);
451            parent.InsertSubtree(index, thetaVar);
452          }
453        }
454        if (node is VariableTreeNode varTreeNode) {
455          var thetaVar = (VariableTreeNode)new Problems.DataAnalysis.Symbolic.Variable().CreateTreeNode();
456          thetaVar.Weight = 1;
457          thetaVar.VariableName = $"θ{n++}";
458
459          thetaNames.Add(thetaVar.VariableName);
460          thetaValues.Add(varTreeNode.Weight);
461
462          var parent = varTreeNode.Parent;
463          if (parent != null) {
464            var index = varTreeNode.Parent.IndexOfSubtree(varTreeNode);
465            parent.RemoveSubtree(index);
466            var prodNode = MakeNode<Multiplication>();
467            varTreeNode.Weight = 1.0;
468            prodNode.AddSubtree(varTreeNode);
469            prodNode.AddSubtree(thetaVar);
470            parent.InsertSubtree(index, prodNode);
471          }
472        }
473      }
474      return copy;
475    }
476
477    private static ISymbolicExpressionTreeNode CreateConstant(double value) {
478      var constantNode = (ConstantTreeNode)new Constant().CreateTreeNode();
479      constantNode.Value = value;
480      return constantNode;
481    }
482
483    private static ISymbolicExpressionTree Subtract(ISymbolicExpressionTree t, ISymbolicExpressionTreeNode b) {
484      var sub = MakeNode<Subtraction>(t.Root.GetSubtree(0).GetSubtree(0), b);
485      t.Root.GetSubtree(0).RemoveSubtree(0);
486      t.Root.GetSubtree(0).InsertSubtree(0, sub);
487      return t;
488    }
489    private static ISymbolicExpressionTree Subtract(ISymbolicExpressionTreeNode b, ISymbolicExpressionTree t) {
490      var sub = MakeNode<Subtraction>(b, t.Root.GetSubtree(0).GetSubtree(0));
491      t.Root.GetSubtree(0).RemoveSubtree(0);
492      t.Root.GetSubtree(0).InsertSubtree(0, sub);
493      return t;
494    }
495
496    private static ISymbolicExpressionTreeNode MakeNode<T>(params ISymbolicExpressionTreeNode[] fs) where T : ISymbol, new() {
497      var node = new T().CreateTreeNode();
498      foreach (var f in fs) node.AddSubtree(f);
499      return node;
500    }
501    #endregion
502  }
503}
Note: See TracBrowser for help on using the repository browser.