github.com/gopherd/gonum@v0.0.4/optimize/functions/vlse.go (about)

     1  // Copyright ©2017 The Gonum Authors. All rights reserved.
     2  // Use of this source code is governed by a BSD-style
     3  // license that can be found in the LICENSE file.
     4  
     5  package functions
     6  
     7  import "math"
     8  
     9  // This file implements functions from the Virtual Library of Simulation Experiments.
    10  //  https://www.sfu.ca/~ssurjano/optimization.html
    11  // In many cases gradients and Hessians have been added. In some cases, these
    12  // are not defined at certain points or manifolds. The gradient in these locations
    13  // has been set to 0.
    14  
    15  // Ackley implements the Ackley function, a function of arbitrary dimension that
    16  // has many local minima. It has a single global minimum of 0 at 0. Its typical
    17  // domain is the hypercube of [-32.768, 32.768]^d.
    18  //  f(x) = -20 * exp(-0.2 sqrt(1/d sum_i x_i^2)) - exp(1/d sum_i cos(2π x_i)) + 20 + exp(1)
    19  // where d is the input dimension.
    20  //
    21  // Reference:
    22  //  https://www.sfu.ca/~ssurjano/ackley.html (obtained June 2017)
    23  type Ackley struct{}
    24  
    25  func (Ackley) Func(x []float64) float64 {
    26  	var ss, sc float64
    27  	for _, v := range x {
    28  		ss += v * v
    29  		sc += math.Cos(2 * math.Pi * v)
    30  	}
    31  	id := 1 / float64(len(x))
    32  	return -20*math.Exp(-0.2*math.Sqrt(id*ss)) - math.Exp(id*sc) + 20 + math.E
    33  }
    34  
    35  // Bukin6 implements Bukin's 6th function. The function is two-dimensional, with
    36  // the typical domain as x_0 ∈ [-15, -5], x_1 ∈ [-3, 3]. The function has a unique
    37  // global minimum at [-10, 1], and many local minima.
    38  //  f(x) = 100 * sqrt(|x_1 - 0.01*x_0^2|) + 0.01*|x_0+10|
    39  // Reference:
    40  //  https://www.sfu.ca/~ssurjano/bukin6.html (obtained June 2017)
    41  type Bukin6 struct{}
    42  
    43  func (Bukin6) Func(x []float64) float64 {
    44  	if len(x) != 2 {
    45  		panic(badInputDim)
    46  	}
    47  	return 100*math.Sqrt(math.Abs(x[1]-0.01*x[0]*x[0])) + 0.01*math.Abs(x[0]+10)
    48  }
    49  
    50  // CamelThree implements the three-hump camel function, a two-dimensional function
    51  // with three local minima, one of which is global.
    52  // The function is given by
    53  //  f(x) = 2*x_0^2 - 1.05*x_0^4 + x_0^6/6 + x_0*x_1 + x_1^2
    54  // with the global minimum at
    55  //  x^* = (0, 0)
    56  //  f(x^*) = 0
    57  // The typical domain is x_i ∈ [-5, 5] for all i.
    58  // Reference:
    59  //  https://www.sfu.ca/~ssurjano/camel3.html (obtained December 2017)
    60  type CamelThree struct{}
    61  
    62  func (c CamelThree) Func(x []float64) float64 {
    63  	if len(x) != 2 {
    64  		panic("camelthree: dimension must be 2")
    65  	}
    66  	x0 := x[0]
    67  	x1 := x[1]
    68  	x02 := x0 * x0
    69  	x04 := x02 * x02
    70  	return 2*x02 - 1.05*x04 + x04*x02/6 + x0*x1 + x1*x1
    71  }
    72  
    73  // CamelSix implements the six-hump camel function, a two-dimensional function.
    74  // with six local minima, two of which are global.
    75  // The function is given by
    76  //  f(x) = (4 - 2.1*x_0^2 + x_0^4/3)*x_0^2 + x_0*x_1 + (-4 + 4*x_1^2)*x_1^2
    77  // with the global minima at
    78  //  x^* = (0.0898, -0.7126), (-0.0898, 0.7126)
    79  //  f(x^*) = -1.0316
    80  // The typical domain is x_0 ∈ [-3, 3], x_1 ∈ [-2, 2].
    81  // Reference:
    82  //  https://www.sfu.ca/~ssurjano/camel6.html (obtained December 2017)
    83  type CamelSix struct{}
    84  
    85  func (c CamelSix) Func(x []float64) float64 {
    86  	if len(x) != 2 {
    87  		panic("camelsix: dimension must be 2")
    88  	}
    89  	x0 := x[0]
    90  	x1 := x[1]
    91  	x02 := x0 * x0
    92  	x12 := x1 * x1
    93  	return (4-2.1*x02+x02*x02/3)*x02 + x0*x1 + (-4+4*x12)*x12
    94  }
    95  
    96  // CrossInTray implements the cross-in-tray function. The cross-in-tray function
    97  // is a two-dimensional function with many local minima, and four global minima
    98  // at (±1.3491, ±1.3491). The function is typically evaluated in the square
    99  // [-10,10]^2.
   100  //  f(x) = -0.001(|sin(x_0)sin(x_1)exp(|100-sqrt((x_0^2+x_1^2)/π)|)|+1)^0.1
   101  // Reference:
   102  //  https://www.sfu.ca/~ssurjano/crossit.html (obtained June 2017)
   103  type CrossInTray struct{}
   104  
   105  func (CrossInTray) Func(x []float64) float64 {
   106  	if len(x) != 2 {
   107  		panic(badInputDim)
   108  	}
   109  	x0 := x[0]
   110  	x1 := x[1]
   111  	exp := math.Abs(100 - math.Sqrt((x0*x0+x1*x1)/math.Pi))
   112  	return -0.0001 * math.Pow(math.Abs(math.Sin(x0)*math.Sin(x1)*math.Exp(exp))+1, 0.1)
   113  }
   114  
   115  // DixonPrice implements the DixonPrice function, a function of arbitrary dimension
   116  // Its typical domain is the hypercube of [-10, 10]^d.
   117  // The function is given by
   118  //  f(x) = (x_0-1)^2 + \sum_{i=1}^{d-1} (i+1) * (2*x_i^2-x_{i-1})^2
   119  // where d is the input dimension. There is a single global minimum, which has
   120  // a location and value of
   121  //  x_i^* = 2^{-(2^{i+1}-2)/(2^{i+1})} for i = 0, ..., d-1.
   122  //  f(x^*) = 0
   123  // Reference:
   124  //  https://www.sfu.ca/~ssurjano/dixonpr.html (obtained June 2017)
   125  type DixonPrice struct{}
   126  
   127  func (DixonPrice) Func(x []float64) float64 {
   128  	xp := x[0]
   129  	v := (xp - 1) * (xp - 1)
   130  	for i := 1; i < len(x); i++ {
   131  		xn := x[i]
   132  		tmp := (2*xn*xn - xp)
   133  		v += float64(i+1) * tmp * tmp
   134  		xp = xn
   135  	}
   136  	return v
   137  }
   138  
   139  // DropWave implements the drop-wave function, a two-dimensional function with
   140  // many local minima and one global minimum at 0. The function is typically evaluated
   141  // in the square [-5.12, 5.12]^2.
   142  //  f(x) = - (1+cos(12*sqrt(x0^2+x1^2))) / (0.5*(x0^2+x1^2)+2)
   143  // Reference:
   144  //  https://www.sfu.ca/~ssurjano/drop.html (obtained June 2017)
   145  type DropWave struct{}
   146  
   147  func (DropWave) Func(x []float64) float64 {
   148  	if len(x) != 2 {
   149  		panic(badInputDim)
   150  	}
   151  	x0 := x[0]
   152  	x1 := x[1]
   153  	num := 1 + math.Cos(12*math.Sqrt(x0*x0+x1*x1))
   154  	den := 0.5*(x0*x0+x1*x1) + 2
   155  	return -num / den
   156  }
   157  
   158  // Eggholder implements the Eggholder function, a two-dimensional function with
   159  // many local minima and one global minimum at [512, 404.2319]. The function
   160  // is typically evaluated in the square [-512, 512]^2.
   161  //  f(x) = -(x_1+47)*sin(sqrt(|x_1+x_0/2+47|))-x_1*sin(sqrt(|x_0-(x_1+47)|))
   162  // Reference:
   163  //  https://www.sfu.ca/~ssurjano/egg.html (obtained June 2017)
   164  type Eggholder struct{}
   165  
   166  func (Eggholder) Func(x []float64) float64 {
   167  	if len(x) != 2 {
   168  		panic(badInputDim)
   169  	}
   170  	x0 := x[0]
   171  	x1 := x[1]
   172  	return -(x1+47)*math.Sin(math.Sqrt(math.Abs(x1+x0/2+47))) -
   173  		x0*math.Sin(math.Sqrt(math.Abs(x0-x1-47)))
   174  }
   175  
   176  // GramacyLee implements the Gramacy-Lee function, a one-dimensional function
   177  // with many local minima. The function is typically evaluated on the domain [0.5, 2.5].
   178  //  f(x) = sin(10πx)/(2x) + (x-1)^4
   179  // Reference:
   180  //  https://www.sfu.ca/~ssurjano/grlee12.html (obtained June 2017)
   181  type GramacyLee struct{}
   182  
   183  func (GramacyLee) Func(x []float64) float64 {
   184  	if len(x) != 1 {
   185  		panic(badInputDim)
   186  	}
   187  	x0 := x[0]
   188  	return math.Sin(10*math.Pi*x0)/(2*x0) + math.Pow(x0-1, 4)
   189  }
   190  
   191  // Griewank implements the Griewank function, a function of arbitrary dimension that
   192  // has many local minima. It has a single global minimum of 0 at 0. Its typical
   193  // domain is the hypercube of [-600, 600]^d.
   194  //  f(x) = \sum_i x_i^2/4000 - \prod_i cos(x_i/sqrt(i)) + 1
   195  // where d is the input dimension.
   196  //
   197  // Reference:
   198  //  https://www.sfu.ca/~ssurjano/griewank.html (obtained June 2017)
   199  type Griewank struct{}
   200  
   201  func (Griewank) Func(x []float64) float64 {
   202  	var ss float64
   203  	pc := 1.0
   204  	for i, v := range x {
   205  		ss += v * v
   206  		pc *= math.Cos(v / math.Sqrt(float64(i+1)))
   207  	}
   208  	return ss/4000 - pc + 1
   209  }
   210  
   211  // HolderTable implements the Holder table function. The Holder table function
   212  // is a two-dimensional function with many local minima, and four global minima
   213  // at (±8.05502, ±9.66459). The function is typically evaluated in the square [-10,10]^2.
   214  //  f(x) = -|sin(x_0)cos(x1)exp(|1-sqrt(x_0^2+x1^2)/π|)|
   215  // Reference:
   216  //  https://www.sfu.ca/~ssurjano/holder.html (obtained June 2017)
   217  type HolderTable struct{}
   218  
   219  func (HolderTable) Func(x []float64) float64 {
   220  	if len(x) != 2 {
   221  		panic(badInputDim)
   222  	}
   223  	x0 := x[0]
   224  	x1 := x[1]
   225  	return -math.Abs(math.Sin(x0) * math.Cos(x1) * math.Exp(math.Abs(1-math.Sqrt(x0*x0+x1*x1)/math.Pi)))
   226  }
   227  
   228  // Langermann2 implements the two-dimensional version of the Langermann function.
   229  // The Langermann function has many local minima. The function is typically
   230  // evaluated in the square [0,10]^2.
   231  //  f(x) = \sum_1^5 c_i exp(-(1/π)\sum_{j=1}^2(x_j-A_{ij})^2) * cos(π\sum_{j=1}^2 (x_j - A_{ij})^2)
   232  //  c = [5]float64{1,2,5,2,3}
   233  //  A = [5][2]float64{{3,5},{5,2},{2,1},{1,4},{7,9}}
   234  // Reference:
   235  //  https://www.sfu.ca/~ssurjano/langer.html (obtained June 2017)
   236  type Langermann2 struct{}
   237  
   238  func (Langermann2) Func(x []float64) float64 {
   239  	if len(x) != 2 {
   240  		panic(badInputDim)
   241  	}
   242  	var (
   243  		c = [5]float64{1, 2, 5, 2, 3}
   244  		A = [5][2]float64{{3, 5}, {5, 2}, {2, 1}, {1, 4}, {7, 9}}
   245  	)
   246  	var f float64
   247  	for i, cv := range c {
   248  		var ss float64
   249  		for j, av := range A[i] {
   250  			xja := x[j] - av
   251  			ss += xja * xja
   252  		}
   253  		f += cv * math.Exp(-(1/math.Pi)*ss) * math.Cos(math.Pi*ss)
   254  	}
   255  	return f
   256  }
   257  
   258  // Levy implements the Levy function, a function of arbitrary dimension that
   259  // has many local minima. It has a single global minimum of 0 at 1. Its typical
   260  // domain is the hypercube of [-10, 10]^d.
   261  //  f(x) = sin^2(π*w_0) + \sum_{i=0}^{d-2}(w_i-1)^2*[1+10sin^2(π*w_i+1)] +
   262  //            (w_{d-1}-1)^2*[1+sin^2(2π*w_{d-1})]
   263  //   w_i = 1 + (x_i-1)/4
   264  // where d is the input dimension.
   265  //
   266  // Reference:
   267  //  https://www.sfu.ca/~ssurjano/levy.html (obtained June 2017)
   268  type Levy struct{}
   269  
   270  func (Levy) Func(x []float64) float64 {
   271  	w1 := 1 + (x[0]-1)/4
   272  	s1 := math.Sin(math.Pi * w1)
   273  	sum := s1 * s1
   274  	for i := 0; i < len(x)-1; i++ {
   275  		wi := 1 + (x[i]-1)/4
   276  		s := math.Sin(math.Pi*wi + 1)
   277  		sum += (wi - 1) * (wi - 1) * (1 + 10*s*s)
   278  	}
   279  	wd := 1 + (x[len(x)-1]-1)/4
   280  	sd := math.Sin(2 * math.Pi * wd)
   281  	return sum + (wd-1)*(wd-1)*(1+sd*sd)
   282  }
   283  
   284  // Levy13 implements the Levy-13 function, a two-dimensional function
   285  // with many local minima. It has a single global minimum of 0 at 1. Its typical
   286  // domain is the square [-10, 10]^2.
   287  //  f(x) = sin^2(3π*x_0) + (x_0-1)^2*[1+sin^2(3π*x_1)] + (x_1-1)^2*[1+sin^2(2π*x_1)]
   288  // Reference:
   289  //  https://www.sfu.ca/~ssurjano/levy13.html (obtained June 2017)
   290  type Levy13 struct{}
   291  
   292  func (Levy13) Func(x []float64) float64 {
   293  	if len(x) != 2 {
   294  		panic(badInputDim)
   295  	}
   296  	x0 := x[0]
   297  	x1 := x[1]
   298  	s0 := math.Sin(3 * math.Pi * x0)
   299  	s1 := math.Sin(3 * math.Pi * x1)
   300  	s2 := math.Sin(2 * math.Pi * x1)
   301  	return s0*s0 + (x0-1)*(x0-1)*(1+s1*s1) + (x1-1)*(x1-1)*(1+s2*s2)
   302  }
   303  
   304  // Rastrigin implements the Rastrigen function, a function of arbitrary dimension
   305  // that has many local minima. It has a single global minimum of 0 at 0. Its typical
   306  // domain is the hypercube of [-5.12, 5.12]^d.
   307  //  f(x) = 10d + \sum_i [x_i^2 - 10cos(2π*x_i)]
   308  // where d is the input dimension.
   309  //
   310  // Reference:
   311  //  https://www.sfu.ca/~ssurjano/rastr.html (obtained June 2017)
   312  type Rastrigin struct{}
   313  
   314  func (Rastrigin) Func(x []float64) float64 {
   315  	sum := 10 * float64(len(x))
   316  	for _, v := range x {
   317  		sum += v*v - 10*math.Cos(2*math.Pi*v)
   318  	}
   319  	return sum
   320  }
   321  
   322  // Schaffer2 implements the second Schaffer function, a two-dimensional function
   323  // with many local minima. It has a single global minimum of 0 at 0. Its typical
   324  // domain is the square [-100, 100]^2.
   325  //  f(x) = 0.5 + (sin^2(x_0^2-x_1^2)-0.5) / (1+0.001*(x_0^2+x_1^2))^2
   326  // Reference:
   327  //  https://www.sfu.ca/~ssurjano/schaffer2.html (obtained June 2017)
   328  type Schaffer2 struct{}
   329  
   330  func (Schaffer2) Func(x []float64) float64 {
   331  	if len(x) != 2 {
   332  		panic(badInputDim)
   333  	}
   334  	x0 := x[0]
   335  	x1 := x[1]
   336  	s := math.Sin(x0*x0 - x1*x1)
   337  	den := 1 + 0.001*(x0*x0+x1*x1)
   338  	return 0.5 + (s*s-0.5)/(den*den)
   339  }
   340  
   341  // Schaffer4 implements the fourth Schaffer function, a two-dimensional function
   342  // with many local minima. Its typical domain is the square [-100, 100]^2.
   343  //  f(x) = 0.5 + (cos(sin(|x_0^2-x_1^2|))-0.5) / (1+0.001*(x_0^2+x_1^2))^2
   344  // Reference:
   345  //  https://www.sfu.ca/~ssurjano/schaffer4.html (obtained June 2017)
   346  type Schaffer4 struct{}
   347  
   348  func (Schaffer4) Func(x []float64) float64 {
   349  	if len(x) != 2 {
   350  		panic(badInputDim)
   351  	}
   352  	x0 := x[0]
   353  	x1 := x[1]
   354  	den := 1 + 0.001*(x0*x0+x1*x1)
   355  	return 0.5 + (math.Cos(math.Sin(math.Abs(x0*x0-x1*x1)))-0.5)/(den*den)
   356  }
   357  
   358  // Schwefel implements the Schwefel function, a function of arbitrary dimension
   359  // that has many local minima. Its typical domain is the hypercube of [-500, 500]^d.
   360  //  f(x) = 418.9829*d - \sum_i x_i*sin(sqrt(|x_i|))
   361  // where d is the input dimension.
   362  //
   363  // Reference:
   364  //  https://www.sfu.ca/~ssurjano/schwef.html (obtained June 2017)
   365  type Schwefel struct{}
   366  
   367  func (Schwefel) Func(x []float64) float64 {
   368  	var sum float64
   369  	for _, v := range x {
   370  		sum += v * math.Sin(math.Sqrt(math.Abs(v)))
   371  	}
   372  	return 418.9829*float64(len(x)) - sum
   373  }
   374  
   375  // Shubert implements the Shubert function, a two-dimensional function
   376  // with many local minima and many global minima. Its typical domain is the
   377  // square [-10, 10]^2.
   378  //  f(x) = (sum_{i=1}^5 i cos((i+1)*x_0+i)) * (\sum_{i=1}^5 i cos((i+1)*x_1+i))
   379  // Reference:
   380  //  https://www.sfu.ca/~ssurjano/shubert.html (obtained June 2017)
   381  type Shubert struct{}
   382  
   383  func (Shubert) Func(x []float64) float64 {
   384  	if len(x) != 2 {
   385  		panic(badInputDim)
   386  	}
   387  	x0 := x[0]
   388  	x1 := x[1]
   389  	var s0, s1 float64
   390  	for i := 1.0; i <= 5.0; i++ {
   391  		s0 += i * math.Cos((i+1)*x0+i)
   392  		s1 += i * math.Cos((i+1)*x1+i)
   393  	}
   394  	return s0 * s1
   395  }