ningshuxia
2022-12-12 e78f5936fee9ab4fff600515bb20a41a28f329c4
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
// <copyright file="FindMinimum.cs" company="Math.NET">
// Math.NET Numerics, part of the Math.NET Project
// http://numerics.mathdotnet.com
// http://github.com/mathnet/mathnet-numerics
//
// Copyright (c) 2009-2017 Math.NET
//
// Permission is hereby granted, free of charge, to any person
// obtaining a copy of this software and associated documentation
// files (the "Software"), to deal in the Software without
// restriction, including without limitation the rights to use,
// copy, modify, merge, publish, distribute, sublicense, and/or sell
// copies of the Software, and to permit persons to whom the
// Software is furnished to do so, subject to the following
// conditions:
//
// The above copyright notice and this permission notice shall be
// included in all copies or substantial portions of the Software.
//
// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND,
// EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES
// OF MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND
// NONINFRINGEMENT. IN NO EVENT SHALL THE AUTHORS OR COPYRIGHT
// HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER LIABILITY,
// WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR
// OTHER DEALINGS IN THE SOFTWARE.
// </copyright>
 
using System;
using IStation.Numerics.LinearAlgebra;
using IStation.Numerics.Optimization;
 
namespace IStation.Numerics
{
    public static class FindMinimum
    {
        /// <summary>
        /// Find value x that minimizes the scalar function f(x), constrained within bounds, using the Golden Section algorithm.
        /// For more options and diagnostics consider to use <see cref="GoldenSectionMinimizer"/> directly.
        /// </summary>
        public static double OfScalarFunctionConstrained(Func<double, double> function, double lowerBound, double upperBound, double tolerance=1e-5, int maxIterations=1000)
        {
            var objective = ObjectiveFunction.ScalarValue(function);
            var result = GoldenSectionMinimizer.Minimum(objective, lowerBound, upperBound, tolerance, maxIterations);
            return result.MinimizingPoint;
        }
 
        /// <summary>
        /// Find vector x that minimizes the function f(x) using the Nelder-Mead Simplex algorithm.
        /// For more options and diagnostics consider to use <see cref="NelderMeadSimplex"/> directly.
        /// </summary>
        public static double OfScalarFunction(Func<double, double> function, double initialGuess, double tolerance = 1e-8, int maxIterations = 1000)
        {
            var objective = ObjectiveFunction.Value(v => function(v[0]));
            var result = NelderMeadSimplex.Minimum(objective, CreateVector.Dense(new[] { initialGuess }), tolerance, maxIterations);
            return result.MinimizingPoint[0];
        }
 
        /// <summary>
        /// Find vector x that minimizes the function f(x) using the Nelder-Mead Simplex algorithm.
        /// For more options and diagnostics consider to use <see cref="NelderMeadSimplex"/> directly.
        /// </summary>
        public static Tuple<double, double> OfFunction(Func<double, double, double> function, double initialGuess0, double initialGuess1, double tolerance = 1e-8, int maxIterations = 1000)
        {
            var objective = ObjectiveFunction.Value(v => function(v[0], v[1]));
            var result = NelderMeadSimplex.Minimum(objective, CreateVector.Dense(new[] { initialGuess0, initialGuess1 }), tolerance, maxIterations);
            return Tuple.Create(result.MinimizingPoint[0], result.MinimizingPoint[1]);
        }
 
        /// <summary>
        /// Find vector x that minimizes the function f(x) using the Nelder-Mead Simplex algorithm.
        /// For more options and diagnostics consider to use <see cref="NelderMeadSimplex"/> directly.
        /// </summary>
        public static Tuple<double, double, double> OfFunction(Func<double, double, double, double> function, double initialGuess0, double initialGuess1, double initialGuess2, double tolerance = 1e-8, int maxIterations = 1000)
        {
            var objective = ObjectiveFunction.Value(v => function(v[0], v[1], v[2]));
            var result = NelderMeadSimplex.Minimum(objective, CreateVector.Dense(new[] { initialGuess0, initialGuess1, initialGuess2 }), tolerance, maxIterations);
            return Tuple.Create(result.MinimizingPoint[0], result.MinimizingPoint[1], result.MinimizingPoint[2]);
        }
 
        /// <summary>
        /// Find vector x that minimizes the function f(x) using the Nelder-Mead Simplex algorithm.
        /// For more options and diagnostics consider to use <see cref="NelderMeadSimplex"/> directly.
        /// </summary>
        public static Vector<double> OfFunction(Func<Vector<double>, double> function, Vector<double> initialGuess, double tolerance=1e-8, int maxIterations=1000)
        {
            var objective = ObjectiveFunction.Value(function);
            var result = NelderMeadSimplex.Minimum(objective, initialGuess, tolerance, maxIterations);
            return result.MinimizingPoint;
        }
 
        /// <summary>
        /// Find vector x that minimizes the function f(x), constrained within bounds, using the Broyden–Fletcher–Goldfarb–Shanno Bounded (BFGS-B) algorithm.
        /// The missing gradient is evaluated numerically (forward difference).
        /// For more options and diagnostics consider to use <see cref="BfgsBMinimizer"/> directly.
        /// </summary>
        public static Vector<double> OfFunctionConstrained(Func<Vector<double>, double> function, Vector<double> lowerBound, Vector<double> upperBound, Vector<double> initialGuess, double gradientTolerance=1e-5, double parameterTolerance=1e-5, double functionProgressTolerance=1e-5, int maxIterations=1000)
        {
            var objective = ObjectiveFunction.Value(function);
            var objectiveWithGradient = new Optimization.ObjectiveFunctions.ForwardDifferenceGradientObjectiveFunction(objective, lowerBound, upperBound);
            var algorithm = new BfgsBMinimizer(gradientTolerance, parameterTolerance, functionProgressTolerance, maxIterations);
            var result = algorithm.FindMinimum(objectiveWithGradient, lowerBound, upperBound, initialGuess);
            return result.MinimizingPoint;
        }
 
        /// <summary>
        /// Find vector x that minimizes the function f(x) using the Broyden–Fletcher–Goldfarb–Shanno (BFGS) algorithm.
        /// For more options and diagnostics consider to use <see cref="BfgsMinimizer"/> directly.
        /// An alternative routine using conjugate gradients (CG) is available in <see cref="ConjugateGradientMinimizer"/>.
        /// </summary>
        public static Vector<double> OfFunctionGradient(Func<Vector<double>, double> function, Func<Vector<double>, Vector<double>> gradient, Vector<double> initialGuess, double gradientTolerance=1e-5, double parameterTolerance=1e-5, double functionProgressTolerance=1e-5, int maxIterations=1000)
        {
            var objective = ObjectiveFunction.Gradient(function, gradient);
            var algorithm = new BfgsMinimizer(gradientTolerance, parameterTolerance, functionProgressTolerance, maxIterations);
            var result = algorithm.FindMinimum(objective, initialGuess);
            return result.MinimizingPoint;
        }
 
        /// <summary>
        /// Find vector x that minimizes the function f(x) using the Broyden–Fletcher–Goldfarb–Shanno (BFGS) algorithm.
        /// For more options and diagnostics consider to use <see cref="BfgsMinimizer"/> directly.
        /// An alternative routine using conjugate gradients (CG) is available in <see cref="ConjugateGradientMinimizer"/>.
        /// </summary>
        public static Vector<double> OfFunctionGradient(Func<Vector<double>, Tuple<double, Vector<double>>> functionGradient, Vector<double> initialGuess, double gradientTolerance=1e-5, double parameterTolerance=1e-5, double functionProgressTolerance=1e-5, int maxIterations=1000)
        {
            var objective = ObjectiveFunction.Gradient(functionGradient);
            var algorithm = new BfgsMinimizer(gradientTolerance, parameterTolerance, functionProgressTolerance, maxIterations);
            var result = algorithm.FindMinimum(objective, initialGuess);
            return result.MinimizingPoint;
        }
 
        /// <summary>
        /// Find vector x that minimizes the function f(x), constrained within bounds, using the Broyden–Fletcher–Goldfarb–Shanno Bounded (BFGS-B) algorithm.
        /// For more options and diagnostics consider to use <see cref="BfgsBMinimizer"/> directly.
        /// </summary>
        public static Vector<double> OfFunctionGradientConstrained(Func<Vector<double>, double> function, Func<Vector<double>, Vector<double>> gradient, Vector<double> lowerBound, Vector<double> upperBound, Vector<double> initialGuess, double gradientTolerance=1e-5, double parameterTolerance=1e-5, double functionProgressTolerance=1e-5, int maxIterations=1000)
        {
            var objective = ObjectiveFunction.Gradient(function, gradient);
            var algorithm = new BfgsBMinimizer(gradientTolerance, parameterTolerance, functionProgressTolerance, maxIterations);
            var result = algorithm.FindMinimum(objective, lowerBound, upperBound, initialGuess);
            return result.MinimizingPoint;
        }
 
        /// <summary>
        /// Find vector x that minimizes the function f(x), constrained within bounds, using the Broyden–Fletcher–Goldfarb–Shanno Bounded (BFGS-B) algorithm.
        /// For more options and diagnostics consider to use <see cref="BfgsBMinimizer"/> directly.
        /// </summary>
        public static Vector<double> OfFunctionGradientConstrained(Func<Vector<double>, Tuple<double, Vector<double>>> functionGradient, Vector<double> lowerBound, Vector<double> upperBound, Vector<double> initialGuess, double gradientTolerance=1e-5, double parameterTolerance=1e-5, double functionProgressTolerance=1e-5, int maxIterations=1000)
        {
            var objective = ObjectiveFunction.Gradient(functionGradient);
            var algorithm = new BfgsBMinimizer(gradientTolerance, parameterTolerance, functionProgressTolerance, maxIterations);
            var result = algorithm.FindMinimum(objective, lowerBound, upperBound, initialGuess);
            return result.MinimizingPoint;
        }
 
        /// <summary>
        /// Find vector x that minimizes the function f(x) using the Newton algorithm.
        /// For more options and diagnostics consider to use <see cref="NewtonMinimizer"/> directly.
        /// </summary>
        public static Vector<double> OfFunctionGradientHessian(Func<Vector<double>, double> function, Func<Vector<double>, Vector<double>> gradient, Func<Vector<double>, Matrix<double>> hessian, Vector<double> initialGuess, double gradientTolerance=1e-8, int maxIterations=1000)
        {
            var objective = ObjectiveFunction.GradientHessian(function, gradient, hessian);
            var result = NewtonMinimizer.Minimum(objective, initialGuess, gradientTolerance, maxIterations);
            return result.MinimizingPoint;
        }
 
        /// <summary>
        /// Find vector x that minimizes the function f(x) using the Newton algorithm.
        /// For more options and diagnostics consider to use <see cref="NewtonMinimizer"/> directly.
        /// </summary>
        public static Vector<double> OfFunctionGradientHessian(Func<Vector<double>, Tuple<double, Vector<double>, Matrix<double>>> functionGradientHessian, Vector<double> initialGuess, double gradientTolerance=1e-8, int maxIterations=1000)
        {
            var objective = ObjectiveFunction.GradientHessian(functionGradientHessian);
            var result = NewtonMinimizer.Minimum(objective, initialGuess, gradientTolerance, maxIterations);
            return result.MinimizingPoint;
        }
    }
}