File: mlgsl_fun.c

package info (click to toggle)
ocamlgsl 1.24.3-1
  • links: PTS, VCS
  • area: main
  • in suites: bullseye
  • size: 1,244 kB
  • sloc: ml: 8,748; ansic: 7,395; makefile: 37
file content (320 lines) | stat: -rw-r--r-- 9,404 bytes parent folder | download | duplicates (3)
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
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
/* gsl-ocaml - OCaml interface to GSL                       */
/* Copyright (©) 2002-2012 - Olivier Andrieu                */
/* Distributed under the terms of the GPL version 3         */

#include <string.h>

#include <gsl/gsl_math.h>

#include <caml/callback.h>
#include <caml/alloc.h>
#include <caml/fail.h>
#include <caml/memory.h>
#include <caml/bigarray.h>

#include "wrappers.h"
#include "mlgsl_fun.h"


/* CALLBACKS */
double gslfun_callback(double x, void *params)
{
  struct callback_params *p=params;
  value res;
  value v_x = copy_double(x);
  res=callback(p->closure, v_x);
  return Double_val(res);
}

/* FDF CALLBACKS */
double gslfun_callback_indir(double x, void *params)
{
  value res;
  value v_x = copy_double(x);
  value *closure = params;
  res=callback(*closure, v_x);
  return Double_val(res);
}
 
double gslfun_callback_f(double x, void *params)
{
  struct callback_params *p=params;
  value res;
  value v_x=copy_double(x);
  res=callback(Field(p->closure, 0), v_x);
  return Double_val(res);
}

double gslfun_callback_df(double x, void *params)
{
  struct callback_params *p=params;
  value res;
  value v_x=copy_double(x);
  res=callback(Field(p->closure, 1), v_x);
  return Double_val(res);
}

void gslfun_callback_fdf(double x, void *params, 
			 double *f, double *df)
{
  struct callback_params *p=params;
  value res;
  value v_x=copy_double(x);
  res=callback(Field(p->closure, 2), v_x);
  *f =Double_val(Field(res, 0));
  *df=Double_val(Field(res, 1));
}


/* MONTE CALLBACKS */
double gsl_monte_callback(double *x_arr, size_t dim, void *params)
{
  struct callback_params *p=params;
  value res;

  memcpy(Double_array_val(p->dbl), x_arr, dim*sizeof(double));
  res=callback(p->closure, p->dbl);
  return Double_val(res);
}

double gsl_monte_callback_fast(double *x_arr, size_t dim, void *params)
{
  struct callback_params *p=params;
  value res;

  res=callback(p->closure, (value)x_arr);
  return Double_val(res);
}



/* MULTIROOT CALLBACKS */
int gsl_multiroot_callback(const gsl_vector *x, void *params, gsl_vector *F)
{
  int barr_flags = BIGARRAY_FLOAT64 | BIGARRAY_C_LAYOUT;
  struct callback_params *p=params;
  value x_barr, f_barr;
  int len = x->size;
  gsl_vector_view x_v, f_v;

  x_barr = alloc_bigarray_dims(barr_flags, 1, NULL, len);
  f_barr = alloc_bigarray_dims(barr_flags, 1, NULL, len);
  x_v = gsl_vector_view_array(Data_bigarray_val(x_barr), len);
  f_v = gsl_vector_view_array(Data_bigarray_val(f_barr), len);

  gsl_vector_memcpy(&x_v.vector, x);
  callback2(p->closure, x_barr, f_barr);
  gsl_vector_memcpy(F, &f_v.vector);
  return GSL_SUCCESS;
}

int gsl_multiroot_callback_f(const gsl_vector *x, void *params, gsl_vector *F)
{
  int barr_flags = BIGARRAY_FLOAT64 | BIGARRAY_C_LAYOUT;
  struct callback_params *p=params;
  value x_barr, f_barr;
  int len = x->size;
  gsl_vector_view x_v, f_v;

  x_barr = alloc_bigarray_dims(barr_flags, 1, NULL, len);
  f_barr = alloc_bigarray_dims(barr_flags, 1, NULL, len);
  x_v = gsl_vector_view_array(Data_bigarray_val(x_barr), len);
  f_v = gsl_vector_view_array(Data_bigarray_val(f_barr), len);

  gsl_vector_memcpy(&x_v.vector, x);
  callback2(Field(p->closure, 0), x_barr, f_barr);
  gsl_vector_memcpy(F, &f_v.vector);
  return GSL_SUCCESS;
}

int gsl_multiroot_callback_df(const gsl_vector *x, void *params, gsl_matrix *J)
{
  int barr_flags = BIGARRAY_FLOAT64 | BIGARRAY_C_LAYOUT;
  struct callback_params *p=params;
  value x_barr, j_barr;
  int len = x->size;
  gsl_vector_view x_v;
  gsl_matrix_view j_v;

  x_barr = alloc_bigarray_dims(barr_flags, 1, NULL, len);
  j_barr = alloc_bigarray_dims(barr_flags, 2, NULL, len, len);
  x_v = gsl_vector_view_array(Data_bigarray_val(x_barr), len);
  j_v = gsl_matrix_view_array(Data_bigarray_val(j_barr), len, len);

  gsl_vector_memcpy(&x_v.vector, x);
  callback2(Field(p->closure, 1), x_barr, j_barr);
  gsl_matrix_memcpy(J, &j_v.matrix);
  return GSL_SUCCESS;
}

int gsl_multiroot_callback_fdf(const gsl_vector *x, void *params, 
			   gsl_vector *F, gsl_matrix *J)
{
  int barr_flags = BIGARRAY_FLOAT64 | BIGARRAY_C_LAYOUT;
  struct callback_params *p=params;
  value x_barr, f_barr, j_barr;
  int len = x->size;
  gsl_vector_view x_v, f_v;
  gsl_matrix_view j_v;
  
  x_barr = alloc_bigarray_dims(barr_flags, 1, NULL, len);
  f_barr = alloc_bigarray_dims(barr_flags, 1, NULL, len);
  j_barr = alloc_bigarray_dims(barr_flags, 2, NULL, len, len);
  x_v = gsl_vector_view_array(Data_bigarray_val(x_barr), len);
  f_v = gsl_vector_view_array(Data_bigarray_val(f_barr), len);
  j_v = gsl_matrix_view_array(Data_bigarray_val(j_barr), len, len);

  gsl_vector_memcpy(&x_v.vector, x);
  callback3(Field(p->closure, 2), x_barr, f_barr, j_barr);
  gsl_vector_memcpy(F, &f_v.vector);
  gsl_matrix_memcpy(J, &j_v.matrix);
  return GSL_SUCCESS;
}



/* MULTIMIN CALLBACKS */
double gsl_multimin_callback(const gsl_vector *x, void *params)
{
  int barr_flags = BIGARRAY_FLOAT64 | BIGARRAY_C_LAYOUT;
  struct callback_params *p=params;
  value x_barr;
  int len = x->size;
  gsl_vector_view x_v;
  value res;

  x_barr = alloc_bigarray_dims(barr_flags, 1, NULL, len);
  x_v = gsl_vector_view_array(Data_bigarray_val(x_barr), len);

  gsl_vector_memcpy(&x_v.vector, x);
  res=callback(p->closure, x_barr);
  return Double_val(res);
}

double gsl_multimin_callback_f(const gsl_vector *x, void *params)
{
  int barr_flags = BIGARRAY_FLOAT64 | BIGARRAY_C_LAYOUT;
  struct callback_params *p=params;
  value x_barr;
  int len = x->size;
  gsl_vector_view x_v;
  value res;

  x_barr = alloc_bigarray_dims(barr_flags, 1, NULL, len);
  x_v = gsl_vector_view_array(Data_bigarray_val(x_barr), len);

  gsl_vector_memcpy(&x_v.vector, x);
  res=callback(Field(p->closure, 0), x_barr);
  return Double_val(res);
}

void gsl_multimin_callback_df(const gsl_vector *x, void *params, gsl_vector *G)
{
  int barr_flags = BIGARRAY_FLOAT64 | BIGARRAY_C_LAYOUT;
  struct callback_params *p=params;
  value x_barr, g_barr;
  int len = x->size;
  gsl_vector_view x_v, g_v;

  x_barr = alloc_bigarray_dims(barr_flags, 1, NULL, len);
  g_barr = alloc_bigarray_dims(barr_flags, 1, NULL, len);
  x_v = gsl_vector_view_array(Data_bigarray_val(x_barr), len);
  g_v = gsl_vector_view_array(Data_bigarray_val(g_barr), len);

  gsl_vector_memcpy(&x_v.vector, x);
  callback2(Field(p->closure, 1), x_barr, g_barr);
  gsl_vector_memcpy(G, &g_v.vector);
}

void gsl_multimin_callback_fdf(const gsl_vector *x, void *params, 
			       double *f, gsl_vector *G)
{
  int barr_flags = BIGARRAY_FLOAT64 | BIGARRAY_C_LAYOUT;
  struct callback_params *p=params;
  value x_barr, g_barr;
  int len = x->size;
  gsl_vector_view x_v, g_v;
  value res;
  
  x_barr = alloc_bigarray_dims(barr_flags, 1, NULL, len);
  g_barr = alloc_bigarray_dims(barr_flags, 1, NULL, len);
  x_v = gsl_vector_view_array(Data_bigarray_val(x_barr), len);
  g_v = gsl_vector_view_array(Data_bigarray_val(g_barr), len);

  gsl_vector_memcpy(&x_v.vector, x);
  res=callback2(Field(p->closure, 2), x_barr, g_barr);
  gsl_vector_memcpy(G, &g_v.vector);
  *f=Double_val(res);
}



/* MULTIFIT CALLBACKS */
int gsl_multifit_callback_f(const gsl_vector *X, void *params, gsl_vector *F)
{
  int barr_flags = BIGARRAY_FLOAT64 | BIGARRAY_C_LAYOUT;
  struct callback_params *parms=params;
  value x_barr, f_barr;
  size_t p = X->size;
  size_t n = F->size;
  gsl_vector_view x_v, f_v;

  x_barr = alloc_bigarray_dims(barr_flags, 1, NULL, p);
  f_barr = alloc_bigarray_dims(barr_flags, 1, NULL, n);
  x_v = gsl_vector_view_array(Data_bigarray_val(x_barr), p);
  f_v = gsl_vector_view_array(Data_bigarray_val(f_barr), n);

  gsl_vector_memcpy(&x_v.vector, X);
  callback2(Field(parms->closure, 0), x_barr, f_barr);
  gsl_vector_memcpy(F, &f_v.vector);
  return GSL_SUCCESS;
}

int gsl_multifit_callback_df(const gsl_vector *X, void *params, gsl_matrix *J)
{
  int barr_flags = BIGARRAY_FLOAT64 | BIGARRAY_C_LAYOUT;
  struct callback_params *parms=params;
  value x_barr, j_barr;
  size_t p = X->size;
  size_t n = J->size1;
  gsl_vector_view x_v;
  gsl_matrix_view j_v;
  value res;

  x_barr = alloc_bigarray_dims(barr_flags, 1, NULL, p);
  j_barr = alloc_bigarray_dims(barr_flags, 2, NULL, n, p);
  x_v = gsl_vector_view_array(Data_bigarray_val(x_barr), p);
  j_v = gsl_matrix_view_array(Data_bigarray_val(j_barr), n, p);

  gsl_vector_memcpy(&x_v.vector, X);
  res=callback2(Field(parms->closure, 1), x_barr, j_barr);
  if(Is_exception_result(res))
    return GSL_FAILURE;
  gsl_matrix_memcpy(J, &j_v.matrix);
  return GSL_SUCCESS;
}

int gsl_multifit_callback_fdf(const gsl_vector *X, void *params, 
			      gsl_vector *F, gsl_matrix *J)
{
  int barr_flags = BIGARRAY_FLOAT64 | BIGARRAY_C_LAYOUT;
  struct callback_params *parms=params;
  value x_barr, f_barr, j_barr;
  size_t p = X->size;
  size_t n = F->size;
  gsl_vector_view x_v, f_v;
  gsl_matrix_view j_v;
  
  x_barr = alloc_bigarray_dims(barr_flags, 1, NULL, p);
  f_barr = alloc_bigarray_dims(barr_flags, 1, NULL, n);
  j_barr = alloc_bigarray_dims(barr_flags, 2, NULL, n, p);
  x_v = gsl_vector_view_array(Data_bigarray_val(x_barr), p);
  f_v = gsl_vector_view_array(Data_bigarray_val(f_barr), n);
  j_v = gsl_matrix_view_array(Data_bigarray_val(j_barr), n, p);

  gsl_vector_memcpy(&x_v.vector, X);
  callback3(Field(parms->closure, 2), x_barr, f_barr, j_barr);
  gsl_vector_memcpy(F, &f_v.vector);
  gsl_matrix_memcpy(J, &j_v.matrix);
  return GSL_SUCCESS;
}