1 /** D bindings for GSL.
2     Authors:    Chibisi Chima-Okereke
3     Copyright:  Copyright (c) 2016, Chibisi Chima-Okereke. All rights reserved.
4     License:    Boost License 1.0
5 */
6 
7 module gsl.multimin;
8 import gsl.math;
9 import gsl.vector;
10 import gsl.matrix;
11 import gsl.min;
12 
13 extern (C):
14 
15 alias gsl_multimin_function_struct gsl_multimin_function;
16 alias gsl_multimin_function_fdf_struct gsl_multimin_function_fdf;
17 
18 extern __gshared const(gsl_multimin_fdfminimizer_type)* gsl_multimin_fdfminimizer_steepest_descent;
19 extern __gshared const(gsl_multimin_fdfminimizer_type)* gsl_multimin_fdfminimizer_conjugate_pr;
20 extern __gshared const(gsl_multimin_fdfminimizer_type)* gsl_multimin_fdfminimizer_conjugate_fr;
21 extern __gshared const(gsl_multimin_fdfminimizer_type)* gsl_multimin_fdfminimizer_vector_bfgs;
22 extern __gshared const(gsl_multimin_fdfminimizer_type)* gsl_multimin_fdfminimizer_vector_bfgs2;
23 extern __gshared const(gsl_multimin_fminimizer_type)* gsl_multimin_fminimizer_nmsimplex;
24 extern __gshared const(gsl_multimin_fminimizer_type)* gsl_multimin_fminimizer_nmsimplex2;
25 extern __gshared const(gsl_multimin_fminimizer_type)* gsl_multimin_fminimizer_nmsimplex2rand;
26 
27 struct gsl_multimin_function_struct
28 {
29     double function (const(gsl_vector)*, void*) f;
30     size_t n;
31     void* params;
32 }
33 
34 struct gsl_multimin_function_fdf_struct
35 {
36     double function (const(gsl_vector)*, void*) f;
37     void function (const(gsl_vector)*, void*, gsl_vector*) df;
38     void function (const(gsl_vector)*, void*, double*, gsl_vector*) fdf;
39     size_t n;
40     void* params;
41 }
42 
43 struct gsl_multimin_fminimizer_type
44 {
45     const(char)* name;
46     size_t size;
47     int function (void*, size_t) alloc;
48     int function (void*, gsl_multimin_function*, const(gsl_vector)*, double*, const(gsl_vector)*) set;
49     int function (void*, gsl_multimin_function*, gsl_vector*, double*, double*) iterate;
50     void function (void*) free;
51 }
52 
53 struct gsl_multimin_fminimizer
54 {
55     const(gsl_multimin_fminimizer_type)* type;
56     gsl_multimin_function* f;
57     double fval;
58     gsl_vector* x;
59     double size;
60     void* state;
61 }
62 
63 struct gsl_multimin_fdfminimizer_type
64 {
65     const(char)* name;
66     size_t size;
67     int function (void*, size_t) alloc;
68     int function (void*, gsl_multimin_function_fdf*, const(gsl_vector)*, double*, gsl_vector*, double, double) set;
69     int function (void*, gsl_multimin_function_fdf*, gsl_vector*, double*, gsl_vector*, gsl_vector*) iterate;
70     int function (void*) restart;
71     void function (void*) free;
72 }
73 
74 struct gsl_multimin_fdfminimizer
75 {
76     const(gsl_multimin_fdfminimizer_type)* type;
77     gsl_multimin_function_fdf* fdf;
78     double f;
79     gsl_vector* x;
80     gsl_vector* gradient;
81     gsl_vector* dx;
82     void* state;
83 }
84 
85 int gsl_multimin_diff (const(gsl_multimin_function)* f, const(gsl_vector)* x, gsl_vector* g);
86 gsl_multimin_fminimizer* gsl_multimin_fminimizer_alloc (const(gsl_multimin_fminimizer_type)* T, size_t n);
87 int gsl_multimin_fminimizer_set (gsl_multimin_fminimizer* s, gsl_multimin_function* f, const(gsl_vector)* x, const(gsl_vector)* step_size);
88 void gsl_multimin_fminimizer_free (gsl_multimin_fminimizer* s);
89 const(char)* gsl_multimin_fminimizer_name (const(gsl_multimin_fminimizer)* s);
90 int gsl_multimin_fminimizer_iterate (gsl_multimin_fminimizer* s);
91 gsl_vector* gsl_multimin_fminimizer_x (const(gsl_multimin_fminimizer)* s);
92 double gsl_multimin_fminimizer_minimum (const(gsl_multimin_fminimizer)* s);
93 double gsl_multimin_fminimizer_size (const(gsl_multimin_fminimizer)* s);
94 int gsl_multimin_test_gradient (const(gsl_vector)* g, double epsabs);
95 int gsl_multimin_test_size (const double size, double epsabs);
96 gsl_multimin_fdfminimizer* gsl_multimin_fdfminimizer_alloc (const(gsl_multimin_fdfminimizer_type)* T, size_t n);
97 int gsl_multimin_fdfminimizer_set (gsl_multimin_fdfminimizer* s, gsl_multimin_function_fdf* fdf, const(gsl_vector)* x, double step_size, double tol);
98 void gsl_multimin_fdfminimizer_free (gsl_multimin_fdfminimizer* s);
99 const(char)* gsl_multimin_fdfminimizer_name (const(gsl_multimin_fdfminimizer)* s);
100 int gsl_multimin_fdfminimizer_iterate (gsl_multimin_fdfminimizer* s);
101 int gsl_multimin_fdfminimizer_restart (gsl_multimin_fdfminimizer* s);
102 gsl_vector* gsl_multimin_fdfminimizer_x (const(gsl_multimin_fdfminimizer)* s);
103 gsl_vector* gsl_multimin_fdfminimizer_dx (const(gsl_multimin_fdfminimizer)* s);
104 gsl_vector* gsl_multimin_fdfminimizer_gradient (const(gsl_multimin_fdfminimizer)* s);
105 double gsl_multimin_fdfminimizer_minimum (const(gsl_multimin_fdfminimizer)* s);