1
/* ocamlgsl - OCaml interface to GSL */
2
/* Copyright (�) 2002 - Olivier Andrieu */
3
/* distributed under the terms of the GPL version 2 */
6
#include <gsl/gsl_blas.h>
8
#include "mlgsl_vector_double.h"
9
#include "mlgsl_matrix_double.h"
11
#include "mlgsl_blas.h"
16
value ml_gsl_blas_ddot(value X, value Y)
19
_DECLARE_VECTOR2(X, Y);
20
_CONVERT_VECTOR2(X, Y);
21
gsl_blas_ddot(&v_X, &v_Y, &r);
22
return copy_double(r);
25
value ml_gsl_blas_dnrm2(value X)
29
return copy_double(gsl_blas_dnrm2(&v_X));
32
value ml_gsl_blas_dasum(value X)
36
return copy_double(gsl_blas_dasum(&v_X));
39
value ml_gsl_blas_idamax(value X)
43
return Val_int(gsl_blas_idamax(&v_X));
46
value ml_gsl_blas_dswap(value X, value Y)
48
_DECLARE_VECTOR2(X, Y);
49
_CONVERT_VECTOR2(X, Y);
50
gsl_blas_dswap(&v_X, &v_Y);
54
value ml_gsl_blas_dcopy(value X, value Y)
56
_DECLARE_VECTOR2(X, Y);
57
_CONVERT_VECTOR2(X, Y);
58
gsl_blas_dcopy(&v_X, &v_Y);
62
value ml_gsl_blas_daxpy(value alpha, value X, value Y)
64
_DECLARE_VECTOR2(X, Y);
65
_CONVERT_VECTOR2(X, Y);
66
gsl_blas_daxpy(Double_val(alpha), &v_X, &v_Y);
70
/* FIXME: drotg drotmg drotm */
72
value ml_gsl_blas_drot(value X, value Y, value c, value s)
74
_DECLARE_VECTOR2(X, Y);
75
_CONVERT_VECTOR2(X, Y);
76
gsl_blas_drot(&v_X, &v_Y, Double_val(c), Double_val(s));
80
value ml_gsl_blas_dscal(value alpha, value X)
84
gsl_blas_dscal(Double_val(alpha), &v_X);
91
value ml_gsl_blas_dgemv(value transa, value alpha, value A,
92
value X, value beta, value Y)
95
_DECLARE_VECTOR2(X, Y);
97
_CONVERT_VECTOR2(X, Y);
98
gsl_blas_dgemv(CBLAS_TRANS_val(transa), Double_val(alpha),
99
&m_A, &v_X, Double_val(beta), &v_Y);
103
value ml_gsl_blas_dgemv_bc(value *argv, int argc)
105
return ml_gsl_blas_dgemv(argv[0], argv[1], argv[2],
106
argv[3], argv[4], argv[5]);
109
value ml_gsl_blas_dtrmv(value uplo, value transa, value diag,
116
gsl_blas_dtrmv(CBLAS_UPLO_val(uplo), CBLAS_TRANS_val(transa),
117
CBLAS_DIAG_val(diag), &m_A, &v_X);
121
value ml_gsl_blas_dtrsv(value uplo, value transa, value diag,
128
gsl_blas_dtrsv(CBLAS_UPLO_val(uplo), CBLAS_TRANS_val(transa),
129
CBLAS_DIAG_val(diag), &m_A, &v_X);
133
value ml_gsl_blas_dsymv(value uplo, value alpha, value A,
134
value X, value beta, value Y)
137
_DECLARE_VECTOR2(X, Y);
139
_CONVERT_VECTOR2(X, Y);
140
gsl_blas_dsymv(CBLAS_UPLO_val(uplo), Double_val(alpha),
141
&m_A, &v_X, Double_val(beta), &v_Y);
145
value ml_gsl_blas_dsymv_bc(value *argv, int argc)
147
return ml_gsl_blas_dsymv(argv[0], argv[1], argv[2],
148
argv[3], argv[4], argv[5]);
151
value ml_gsl_blas_dger(value alpha, value X, value Y, value A)
154
_DECLARE_VECTOR2(X, Y);
156
_CONVERT_VECTOR2(X, Y);
157
gsl_blas_dger(Double_val(alpha), &v_X, &v_Y, &m_A);
161
value ml_gsl_blas_dsyr(value uplo ,value alpha, value X, value A)
167
gsl_blas_dsyr(CBLAS_UPLO_val(uplo), Double_val(alpha),
172
value ml_gsl_blas_dsyr2(value uplo ,value alpha, value X, value Y, value A)
175
_DECLARE_VECTOR2(X, Y);
177
_CONVERT_VECTOR2(X, Y);
178
gsl_blas_dsyr2(CBLAS_UPLO_val(uplo), Double_val(alpha),
187
value ml_gsl_blas_dgemm(value transa, value transb,
188
value alpha, value A, value B,
191
_DECLARE_MATRIX3(A, B, C);
192
_CONVERT_MATRIX3(A, B, C);
193
gsl_blas_dgemm(CBLAS_TRANS_val(transa), CBLAS_TRANS_val(transb),
194
Double_val(alpha), &m_A, &m_B, Double_val(beta), &m_C);
198
value ml_gsl_blas_dgemm_bc(value *argv, int argc)
200
return ml_gsl_blas_dgemm(argv[0], argv[1], argv[2],
201
argv[3], argv[4], argv[5], argv[6]);
205
value ml_gsl_blas_dsymm(value side, value uplo,
206
value alpha, value A, value B,
209
_DECLARE_MATRIX3(A, B, C);
210
_CONVERT_MATRIX3(A, B, C);
211
gsl_blas_dsymm(CBLAS_SIDE_val(side), CBLAS_UPLO_val(uplo),
212
Double_val(alpha), &m_A, &m_B, Double_val(beta), &m_C);
216
value ml_gsl_blas_dsymm_bc(value *argv, int argc)
218
return ml_gsl_blas_dsymm(argv[0], argv[1], argv[2],
219
argv[3], argv[4], argv[5], argv[6]);
222
value ml_gsl_blas_dtrmm(value side, value uplo,
223
value transa, value diag,
224
value alpha, value A, value B)
226
_DECLARE_MATRIX2(A, B);
227
_CONVERT_MATRIX2(A, B);
228
gsl_blas_dtrmm(CBLAS_SIDE_val(side), CBLAS_UPLO_val(uplo),
229
CBLAS_TRANS_val(transa), CBLAS_DIAG_val(diag),
230
Double_val(alpha), &m_A, &m_B);
234
value ml_gsl_blas_dtrmm_bc(value *argv, int argc)
236
return ml_gsl_blas_dtrmm(argv[0], argv[1], argv[2],
237
argv[3], argv[4], argv[5], argv[6]);
240
value ml_gsl_blas_dtrsm(value side, value uplo,
241
value transa, value diag,
242
value alpha, value A, value B)
244
_DECLARE_MATRIX2(A, B);
245
_CONVERT_MATRIX2(A, B);
246
gsl_blas_dtrsm(CBLAS_SIDE_val(side), CBLAS_UPLO_val(uplo),
247
CBLAS_TRANS_val(transa), CBLAS_DIAG_val(diag),
248
Double_val(alpha), &m_A, &m_B);
252
value ml_gsl_blas_dtrsm_bc(value *argv, int argc)
254
return ml_gsl_blas_dtrsm(argv[0], argv[1], argv[2],
255
argv[3], argv[4], argv[5], argv[6]);
258
value ml_gsl_blas_dsyrk(value uplo, value trans, value alpha,
259
value A, value beta, value C)
261
_DECLARE_MATRIX2(A, C);
262
_CONVERT_MATRIX2(A, C);
263
gsl_blas_dsyrk(CBLAS_UPLO_val(uplo), CBLAS_TRANS_val(trans),
264
Double_val(alpha), &m_A,
265
Double_val(beta), &m_C);
269
value ml_gsl_blas_dsyrk_bc(value *argv, int argc)
271
return ml_gsl_blas_dsyrk(argv[0], argv[1], argv[2],
272
argv[3], argv[4], argv[5]);
276
value ml_gsl_blas_dsyr2k(value uplo, value trans, value alpha,
277
value A, value B, value beta, value C)
279
_DECLARE_MATRIX3(A, B, C);
280
_CONVERT_MATRIX3(A, B, C);
281
gsl_blas_dsyr2k(CBLAS_UPLO_val(uplo), CBLAS_TRANS_val(trans),
282
Double_val(alpha), &m_A, &m_B,
283
Double_val(beta), &m_C);
287
value ml_gsl_blas_dsyr2k_bc(value *argv, int argc)
289
return ml_gsl_blas_dsyr2k(argv[0], argv[1], argv[2],
290
argv[3], argv[4], argv[5], argv[6]);