gsl-ocaml / mlgsl_fun.c

  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
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
/* ocamlgsl - OCaml interface to GSL                        */
/* Copyright (©) 2002-2005 - Olivier Andrieu                */
/* distributed under the terms of the GPL version 2         */

#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_exn(p->closure, v_x);
  if(Is_exception_result(res))
    return GSL_NAN;
  return Double_val(res);
}

double gslfun_callback_indir(double x, void *params)
{
  value res;
  value v_x = copy_double(x);
  value *closure = params;
  res=callback_exn(*closure, v_x);
  if(Is_exception_result(res))
    return GSL_NAN;
  return Double_val(res);
}
 
/* FDF CALLBACKS */
double gslfun_callback_f(double x, void *params)
{
  struct callback_params *p=params;
  value res;
  value v_x=copy_double(x);
  res=callback_exn(Field(p->closure, 0), v_x);
  if(Is_exception_result(res))
    return GSL_NAN;
  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_exn(Field(p->closure, 1), v_x);
  if(Is_exception_result(res))
    return GSL_NAN;
  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_exn(Field(p->closure, 2), v_x);
  if(Is_exception_result(res)){
    *f=GSL_NAN; *df=GSL_NAN;
    return;
  }
  *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_exn(p->closure, p->dbl);
  if(Is_exception_result(res))
    return GSL_NAN;
  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_exn(p->closure, (value)x_arr);
  if(Is_exception_result(res))
    return GSL_NAN;
  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 | BIGARRAY_EXTERNAL;
  struct callback_params *p=params;
  value x_barr, f_barr;
  int len = x->size;
  LOCALARRAY(double, x_arr, len); 
  LOCALARRAY(double, f_arr, len); 
  gsl_vector_view x_v, f_v;
  value res;

  x_v = gsl_vector_view_array(x_arr, len);
  f_v = gsl_vector_view_array(f_arr, len);
  gsl_vector_memcpy(&x_v.vector, x);
  x_barr = alloc_bigarray_dims(barr_flags, 1, x_arr, len);
  f_barr = alloc_bigarray_dims(barr_flags, 1, f_arr, len);
  res=callback2_exn(p->closure, x_barr, f_barr);
  if(Is_exception_result(res))
    return GSL_FAILURE;
  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 | BIGARRAY_EXTERNAL;
  struct callback_params *p=params;
  value x_barr, f_barr;
  int len = x->size;
  LOCALARRAY(double, x_arr, len); 
  LOCALARRAY(double, f_arr, len); 
  gsl_vector_view x_v, f_v;
  value res;

  x_v = gsl_vector_view_array(x_arr, len);
  f_v = gsl_vector_view_array(f_arr, len);
  gsl_vector_memcpy(&x_v.vector, x);
  x_barr = alloc_bigarray_dims(barr_flags, 1, x_arr, len);
  f_barr = alloc_bigarray_dims(barr_flags, 1, f_arr, len);
  res=callback2_exn(Field(p->closure, 0), x_barr, f_barr);
  if(Is_exception_result(res))
    return GSL_FAILURE;
  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 | BIGARRAY_EXTERNAL;
  struct callback_params *p=params;
  value x_barr, j_barr;
  int len = x->size;
  LOCALARRAY(double, x_arr, len); 
  LOCALARRAY(double, j_arr, len*len); 
  gsl_vector_view x_v;
  gsl_matrix_view j_v;
  value res;

  x_v = gsl_vector_view_array(x_arr, len);
  j_v = gsl_matrix_view_array(j_arr, len, len);
  gsl_vector_memcpy(&x_v.vector, x);
  x_barr = alloc_bigarray_dims(barr_flags, 1, x_arr, len);
  j_barr = alloc_bigarray_dims(barr_flags, 2, j_arr, len, len);
  res=callback2_exn(Field(p->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_multiroot_callback_fdf(const gsl_vector *x, void *params, 
			   gsl_vector *F, gsl_matrix *J)
{
  int barr_flags = BIGARRAY_FLOAT64 | BIGARRAY_C_LAYOUT | BIGARRAY_EXTERNAL;
  struct callback_params *p=params;
  value x_barr, f_barr, j_barr;
  int len = x->size;
  LOCALARRAY(double, x_arr, len); 
  LOCALARRAY(double, f_arr, len); 
  LOCALARRAY(double, j_arr, len*len); 
  gsl_vector_view x_v, f_v;
  gsl_matrix_view j_v;
  value res;
  
  x_v = gsl_vector_view_array(x_arr, len);
  f_v = gsl_vector_view_array(f_arr, len);
  j_v = gsl_matrix_view_array(j_arr, len, len);
  gsl_vector_memcpy(&x_v.vector, x);
  x_barr = alloc_bigarray_dims(barr_flags, 1, x_arr, len);
  f_barr = alloc_bigarray_dims(barr_flags, 1, f_arr, len);
  j_barr = alloc_bigarray_dims(barr_flags, 2, j_arr, len, len);
  res=callback3_exn(Field(p->closure, 2), x_barr, f_barr, j_barr);
  if(Is_exception_result(res))
    return GSL_FAILURE;
  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 | BIGARRAY_EXTERNAL;
  struct callback_params *p=params;
  value x_barr;
  int len = x->size;

/* CR mmottl: Stack allocation may segfault with large lengths, hence
   malloc.  Note that the OCaml callbacks may put the bigarrays into
   references.  This is evil, and these bindings should really get fixed
   to avoid problems of that sort. */
  double *x_arr = malloc(sizeof(double) * len);
#if 0
  LOCALARRAY(double, x_arr, len);
#endif

  gsl_vector_view x_v;
  value res;

  x_v = gsl_vector_view_array(x_arr, len);
  gsl_vector_memcpy(&x_v.vector, x);
  x_barr = alloc_bigarray_dims(barr_flags, 1, x_arr, len);
  res=callback_exn(p->closure, x_barr);

  /* CR mmottl: need to free malloced memory now */
  free(x_arr);

  if(Is_exception_result(res)) {
    return GSL_NAN;
  }

  return Double_val(res);
}

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

/* CR mmottl: Stack allocation may segfault with large lengths, hence
   malloc.  Note that the OCaml callbacks may put the bigarrays into
   references.  This is evil, and these bindings should really get fixed
   to avoid problems of that sort. */
  double *x_arr = malloc(sizeof(double) * len);
#if 0
  LOCALARRAY(double, x_arr, len);
#endif

  gsl_vector_view x_v;
  value res;

  x_v = gsl_vector_view_array(x_arr, len);
  gsl_vector_memcpy(&x_v.vector, x);
  x_barr = alloc_bigarray_dims(barr_flags, 1, x_arr, len);
  res=callback_exn(Field(p->closure, 0), x_barr);

  /* CR mmottl: need to free malloced memory now */
  free(x_arr);

  if(Is_exception_result(res))
    return GSL_NAN;

  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 | BIGARRAY_EXTERNAL;
  struct callback_params *p=params;
  value x_barr, g_barr;
  int len = x->size;

/* CR mmottl: Stack allocation may segfault with large lengths, hence
   malloc.  Note that the OCaml callbacks may put the bigarrays into
   references.  This is evil, and these bindings should really get fixed
   to avoid problems of that sort. */
  double *x_arr = malloc(sizeof(double) * len);
  double *g_arr = malloc(sizeof(double) * len);
#if 0
  LOCALARRAY(double, x_arr, len);
  LOCALARRAY(double, g_arr, len); 
#endif

  gsl_vector_view x_v, g_v;
  value res;

  x_v = gsl_vector_view_array(x_arr, len);
  g_v = gsl_vector_view_array(g_arr, len);
  gsl_vector_memcpy(&x_v.vector, x);
  x_barr = alloc_bigarray_dims(barr_flags, 1, x_arr, len);
  g_barr = alloc_bigarray_dims(barr_flags, 1, g_arr, len);
  res=callback2_exn(Field(p->closure, 1), x_barr, g_barr);

  if(Is_exception_result(res)){
    /* CR mmottl: need to free malloced memory now */
    free(x_arr);
    free(g_arr);

    /* the caml functions raised an exception but there's no way we can
       indicate this to GSL since the return type is void.
       So we set the out param G to NaN. */
    gsl_vector_set_all(G, GSL_NAN);
    return;
  }

  gsl_vector_memcpy(G, &g_v.vector);

  /* CR mmottl: need to free malloced memory now */
  free(x_arr);
  free(g_arr);

}

void gsl_multimin_callback_fdf(const gsl_vector *x, void *params, 
			       double *f, gsl_vector *G)
{
  int barr_flags = BIGARRAY_FLOAT64 | BIGARRAY_C_LAYOUT | BIGARRAY_EXTERNAL;
  struct callback_params *p=params;
  value x_barr, g_barr;
  int len = x->size;

/* CR mmottl: Stack allocation may segfault with large lengths, hence
   malloc.  Note that the OCaml callbacks may put the bigarrays into
   references.  This is evil, and these bindings should really get fixed
   to avoid problems of that sort. */
  double *x_arr = malloc(sizeof(double) * len);
  double *g_arr = malloc(sizeof(double) * len);
#if 0
  LOCALARRAY(double, x_arr, len);
  LOCALARRAY(double, g_arr, len); 
#endif

  gsl_vector_view x_v, g_v;
  value res;
  
  x_v = gsl_vector_view_array(x_arr, len);
  g_v = gsl_vector_view_array(g_arr, len);
  gsl_vector_memcpy(&x_v.vector, x);
  x_barr = alloc_bigarray_dims(barr_flags, 1, x_arr, len);
  g_barr = alloc_bigarray_dims(barr_flags, 1, g_arr, len);
  res=callback2_exn(Field(p->closure, 2), x_barr, g_barr);

  if(Is_exception_result(res)){
    /* CR mmottl: need to free malloced memory now */
    free(x_arr);
    free(g_arr);

    *f=GSL_NAN;
    return;
  }
  gsl_vector_memcpy(G, &g_v.vector);

  /* CR mmottl: need to free malloced memory now */
  free(x_arr);
  free(g_arr);

  *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 | BIGARRAY_EXTERNAL;
  struct callback_params *parms=params;
  value x_barr, f_barr;
  size_t p = X->size;
  size_t n = F->size;
  LOCALARRAY(double, x_arr, p); 
  LOCALARRAY(double, f_arr, n); 
  gsl_vector_view x_v, f_v;
  value res;

  x_v = gsl_vector_view_array(x_arr, p);
  f_v = gsl_vector_view_array(f_arr, n);
  gsl_vector_memcpy(&x_v.vector, X);
  x_barr = alloc_bigarray_dims(barr_flags, 1, x_arr, p);
  f_barr = alloc_bigarray_dims(barr_flags, 1, f_arr, n);
  res=callback2_exn(Field(parms->closure, 0), x_barr, f_barr);
  if(Is_exception_result(res))
    return GSL_FAILURE;
  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 | BIGARRAY_EXTERNAL;
  struct callback_params *parms=params;
  value x_barr, j_barr;
  size_t p = X->size;
  size_t n = J->size1;
  LOCALARRAY(double, x_arr, p); 
  LOCALARRAY(double, j_arr, n*p); 
  gsl_vector_view x_v;
  gsl_matrix_view j_v;
  value res;

  x_v = gsl_vector_view_array(x_arr, p);
  j_v = gsl_matrix_view_array(j_arr, n, p);
  gsl_vector_memcpy(&x_v.vector, X);
  x_barr = alloc_bigarray_dims(barr_flags, 1, x_arr, p);
  j_barr = alloc_bigarray_dims(barr_flags, 2, j_arr, n, p);
  res=callback2_exn(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 | BIGARRAY_EXTERNAL;
  struct callback_params *parms=params;
  value x_barr, f_barr, j_barr;
  size_t p = X->size;
  size_t n = F->size;
  LOCALARRAY(double, x_arr, p); 
  LOCALARRAY(double, f_arr, n); 
  LOCALARRAY(double, j_arr, n*p); 
  gsl_vector_view x_v, f_v;
  gsl_matrix_view j_v;
  value res;
  
  x_v = gsl_vector_view_array(x_arr, p);
  f_v = gsl_vector_view_array(f_arr, n);
  j_v = gsl_matrix_view_array(j_arr, n, p);
  gsl_vector_memcpy(&x_v.vector, X);
  x_barr = alloc_bigarray_dims(barr_flags, 1, x_arr, p);
  f_barr = alloc_bigarray_dims(barr_flags, 1, f_arr, n);
  j_barr = alloc_bigarray_dims(barr_flags, 2, j_arr, n, p);
  res=callback3_exn(Field(parms->closure, 2), x_barr, f_barr, j_barr);
  if(Is_exception_result(res))
    return GSL_FAILURE;
  gsl_vector_memcpy(F, &f_v.vector);
  gsl_matrix_memcpy(J, &j_v.matrix);
  return GSL_SUCCESS;
}
Tip: Filter by directory path e.g. /media app.js to search for public/media/app.js.
Tip: Use camelCasing e.g. ProjME to search for ProjectModifiedEvent.java.
Tip: Filter by extension type e.g. /repo .js to search for all .js files in the /repo directory.
Tip: Separate your search with spaces e.g. /ssh pom.xml to search for src/ssh/pom.xml.
Tip: Use ↑ and ↓ arrow keys to navigate and return to view the file.
Tip: You can also navigate files with Ctrl+j (next) and Ctrl+k (previous) and view the file with Ctrl+o.
Tip: You can also navigate files with Alt+j (next) and Alt+k (previous) and view the file with Alt+o.