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
|
/* ocamlgsl - OCaml interface to GSL */
/* Copyright (©) 2002-2005 - Olivier Andrieu */
/* distributed under the terms of the GPL version 2 */
#include <gsl/gsl_blas.h>
#include "mlgsl_complex.h"
#include "mlgsl_vector_complex.h"
#include "mlgsl_matrix_complex.h"
#include "mlgsl_blas.h"
/* LEVEL1 complex */
CAMLprim value ml_gsl_blas_zdotu(value X, value Y)
{
gsl_complex r;
_DECLARE_VECTOR2(X, Y);
_CONVERT_VECTOR2(X, Y);
gsl_blas_zdotu(&v_X, &v_Y, &r);
return copy_complex(&r);
}
CAMLprim value ml_gsl_blas_zdotc(value X, value Y)
{
gsl_complex r;
_DECLARE_VECTOR2(X, Y);
_CONVERT_VECTOR2(X, Y);
gsl_blas_zdotc(&v_X, &v_Y, &r);
return copy_complex(&r);
}
CAMLprim value ml_gsl_blas_znrm2(value X)
{
_DECLARE_VECTOR(X);
_CONVERT_VECTOR(X);
return copy_double(gsl_blas_dznrm2(&v_X));
}
CAMLprim value ml_gsl_blas_zasum(value X)
{
_DECLARE_VECTOR(X);
_CONVERT_VECTOR(X);
return copy_double(gsl_blas_dzasum(&v_X));
}
CAMLprim value ml_gsl_blas_izamax(value X)
{
_DECLARE_VECTOR(X);
_CONVERT_VECTOR(X);
return Val_int(gsl_blas_izamax(&v_X));
}
CAMLprim value ml_gsl_blas_zswap(value X, value Y)
{
_DECLARE_VECTOR2(X, Y);
_CONVERT_VECTOR2(X, Y);
gsl_blas_zswap(&v_X, &v_Y);
return Val_unit;
}
CAMLprim value ml_gsl_blas_zcopy(value X, value Y)
{
_DECLARE_VECTOR2(X, Y);
_CONVERT_VECTOR2(X, Y);
gsl_blas_zcopy(&v_X, &v_Y);
return Val_unit;
}
CAMLprim value ml_gsl_blas_zaxpy(value alpha, value X, value Y)
{
_DECLARE_VECTOR2(X, Y);
_DECLARE_COMPLEX(alpha);
_CONVERT_VECTOR2(X, Y);
_CONVERT_COMPLEX(alpha);
gsl_blas_zaxpy(z_alpha, &v_X, &v_Y);
return Val_unit;
}
CAMLprim value ml_gsl_blas_zscal(value alpha, value X)
{
_DECLARE_VECTOR(X);
_DECLARE_COMPLEX(alpha);
_CONVERT_VECTOR(X);
_CONVERT_COMPLEX(alpha);
gsl_blas_zscal(z_alpha, &v_X);
return Val_unit;
}
CAMLprim value ml_gsl_blas_zdscal(value alpha, value X)
{
_DECLARE_VECTOR(X);
_CONVERT_VECTOR(X);
gsl_blas_zdscal(Double_val(alpha), &v_X);
return Val_unit;
}
/* LEVEL2 complex */
CAMLprim value ml_gsl_blas_zgemv(value transa, value alpha, value A,
value X, value beta, value Y)
{
_DECLARE_MATRIX(A);
_DECLARE_COMPLEX2(alpha, beta);
_DECLARE_VECTOR2(X, Y);
_CONVERT_MATRIX(A);
_CONVERT_COMPLEX2(alpha, beta);
_CONVERT_VECTOR2(X, Y);
gsl_blas_zgemv(CBLAS_TRANS_val(transa), z_alpha,
&m_A, &v_X, z_beta, &v_Y);
return Val_unit;
}
CAMLprim value ml_gsl_blas_zgemv_bc(value *argv, int argc)
{
return ml_gsl_blas_zgemv(argv[0], argv[1], argv[2],
argv[3], argv[4], argv[5]);
}
CAMLprim value ml_gsl_blas_ztrmv(value uplo, value transa, value diag,
value A, value X)
{
_DECLARE_MATRIX(A);
_DECLARE_VECTOR(X);
_CONVERT_MATRIX(A);
_CONVERT_VECTOR(X);
gsl_blas_ztrmv(CBLAS_UPLO_val(uplo), CBLAS_TRANS_val(transa),
CBLAS_DIAG_val(diag), &m_A, &v_X);
return Val_unit;
}
CAMLprim value ml_gsl_blas_ztrsv(value uplo, value transa, value diag,
value A, value X)
{
_DECLARE_MATRIX(A);
_DECLARE_VECTOR(X);
_CONVERT_MATRIX(A);
_CONVERT_VECTOR(X);
gsl_blas_ztrsv(CBLAS_UPLO_val(uplo), CBLAS_TRANS_val(transa),
CBLAS_DIAG_val(diag), &m_A, &v_X);
return Val_unit;
}
CAMLprim value ml_gsl_blas_zhemv(value uplo, value alpha, value A,
value X, value beta, value Y)
{
_DECLARE_MATRIX(A);
_DECLARE_COMPLEX2(alpha, beta);
_DECLARE_VECTOR2(X,Y);
_CONVERT_MATRIX(A);
_CONVERT_COMPLEX2(alpha, beta);
_CONVERT_VECTOR2(X,Y);
gsl_blas_zhemv(CBLAS_UPLO_val(uplo), z_alpha, &m_A, &v_X,
z_beta, &v_Y);
return Val_unit;
}
CAMLprim value ml_gsl_blas_zhemv_bc(value *argv, int argc)
{
return ml_gsl_blas_zhemv(argv[0], argv[1], argv[2],
argv[3], argv[4], argv[5]);
}
CAMLprim value ml_gsl_blas_zgeru(value alpha, value X, value Y, value A)
{
_DECLARE_MATRIX(A);
_DECLARE_VECTOR2(X,Y);
_DECLARE_COMPLEX(alpha);
_CONVERT_MATRIX(A);
_CONVERT_VECTOR2(X,Y);
_CONVERT_COMPLEX(alpha);
gsl_blas_zgeru(z_alpha, &v_X, &v_Y, &m_A);
return Val_unit;
}
CAMLprim value ml_gsl_blas_zgerc(value alpha, value X, value Y, value A)
{
_DECLARE_MATRIX(A);
_DECLARE_VECTOR2(X,Y);
_DECLARE_COMPLEX(alpha);
_CONVERT_MATRIX(A);
_CONVERT_VECTOR2(X,Y);
_CONVERT_COMPLEX(alpha);
gsl_blas_zgerc(z_alpha, &v_X, &v_Y, &m_A);
return Val_unit;
}
CAMLprim value ml_gsl_blas_zher(value uplo, value alpha, value X, value A)
{
_DECLARE_MATRIX(A);
_DECLARE_VECTOR(X);
_CONVERT_MATRIX(A);
_CONVERT_VECTOR(X);
gsl_blas_zher(CBLAS_UPLO_val(uplo), Double_val(alpha), &v_X, &m_A);
return Val_unit;
}
CAMLprim value ml_gsl_blas_zher2(value uplo, value alpha,
value X, value Y, value A)
{
_DECLARE_MATRIX(A);
_DECLARE_VECTOR2(X,Y);
_DECLARE_COMPLEX(alpha);
_CONVERT_MATRIX(A);
_CONVERT_VECTOR2(X,Y);
_CONVERT_COMPLEX(alpha);
gsl_blas_zher2(CBLAS_UPLO_val(uplo), z_alpha, &v_X, &v_Y, &m_A);
return Val_unit;
}
/* LEVEL3 complex */
CAMLprim value ml_gsl_blas_zgemm(value transa, value transb,
value alpha, value A, value B,
value beta, value C)
{
_DECLARE_MATRIX3(A, B, C);
_DECLARE_COMPLEX2(alpha, beta);
_CONVERT_MATRIX3(A, B, C);
_CONVERT_COMPLEX2(alpha, beta);
gsl_blas_zgemm(CBLAS_TRANS_val(transa), CBLAS_TRANS_val(transb),
z_alpha, &m_A, &m_B, z_beta, &m_C);
return Val_unit;
}
CAMLprim value ml_gsl_blas_zgemm_bc(value *argv, int argc)
{
return ml_gsl_blas_zgemm(argv[0], argv[1], argv[2],
argv[3], argv[4], argv[5], argv[6]);
}
CAMLprim value ml_gsl_blas_zsymm(value side, value uplo,
value alpha, value A, value B,
value beta, value C)
{
_DECLARE_MATRIX3(A, B, C);
_DECLARE_COMPLEX2(alpha, beta);
_CONVERT_MATRIX3(A, B, C);
_CONVERT_COMPLEX2(alpha, beta);
gsl_blas_zsymm(CBLAS_SIDE_val(side), CBLAS_UPLO_val(uplo),
z_alpha, &m_A, &m_B, z_beta, &m_C);
return Val_unit;
}
CAMLprim value ml_gsl_blas_zsymm_bc(value *argv, int argc)
{
return ml_gsl_blas_zsymm(argv[0], argv[1], argv[2],
argv[3], argv[4], argv[5], argv[6]);
}
CAMLprim value ml_gsl_blas_zsyrk(value uplo, value trans, value alpha,
value A, value beta, value C)
{
_DECLARE_MATRIX2(A, C);
_DECLARE_COMPLEX2(alpha, beta);
_CONVERT_MATRIX2(A, C);
_CONVERT_COMPLEX2(alpha, beta);
gsl_blas_zsyrk(CBLAS_UPLO_val(uplo), CBLAS_TRANS_val(trans),
z_alpha, &m_A,
z_beta, &m_C);
return Val_unit;
}
CAMLprim value ml_gsl_blas_zsyrk_bc(value *argv, int argc)
{
return ml_gsl_blas_zsyrk(argv[0], argv[1], argv[2],
argv[3], argv[4], argv[5]);
}
CAMLprim value ml_gsl_blas_zsyr2k(value uplo, value trans, value alpha,
value A, value B, value beta, value C)
{
_DECLARE_MATRIX3(A, B, C);
_DECLARE_COMPLEX2(alpha, beta);
_CONVERT_MATRIX3(A, B, C);
_CONVERT_COMPLEX2(alpha, beta);
gsl_blas_zsyr2k(CBLAS_UPLO_val(uplo), CBLAS_TRANS_val(trans),
z_alpha, &m_A, &m_B,
z_beta, &m_C);
return Val_unit;
}
CAMLprim value ml_gsl_blas_zsyr2k_bc(value *argv, int argc)
{
return ml_gsl_blas_zsyr2k(argv[0], argv[1], argv[2],
argv[3], argv[4], argv[5], argv[6]);
}
CAMLprim value ml_gsl_blas_ztrmm(value side, value uplo,
value transa, value diag,
value alpha, value A, value B)
{
_DECLARE_MATRIX2(A, B);
_DECLARE_COMPLEX(alpha);
_CONVERT_MATRIX2(A, B);
_CONVERT_COMPLEX(alpha);
gsl_blas_ztrmm(CBLAS_SIDE_val(side), CBLAS_UPLO_val(uplo),
CBLAS_TRANS_val(transa), CBLAS_DIAG_val(diag),
z_alpha, &m_A, &m_B);
return Val_unit;
}
CAMLprim value ml_gsl_blas_ztrmm_bc(value *argv, int argc)
{
return ml_gsl_blas_ztrmm(argv[0], argv[1], argv[2],
argv[3], argv[4], argv[5], argv[6]);
}
CAMLprim value ml_gsl_blas_ztrsm(value side, value uplo,
value transa, value diag,
value alpha, value A, value B)
{
_DECLARE_MATRIX2(A, B);
_DECLARE_COMPLEX(alpha);
_CONVERT_MATRIX2(A, B);
_CONVERT_COMPLEX(alpha);
gsl_blas_ztrsm(CBLAS_SIDE_val(side), CBLAS_UPLO_val(uplo),
CBLAS_TRANS_val(transa), CBLAS_DIAG_val(diag),
z_alpha, &m_A, &m_B);
return Val_unit;
}
CAMLprim value ml_gsl_blas_ztrsm_bc(value *argv, int argc)
{
return ml_gsl_blas_ztrsm(argv[0], argv[1], argv[2],
argv[3], argv[4], argv[5], argv[6]);
}
CAMLprim value ml_gsl_blas_zhemm(value side, value uplo,
value alpha, value A, value B,
value beta, value C)
{
_DECLARE_MATRIX3(A, B, C);
_DECLARE_COMPLEX2(alpha, beta);
_CONVERT_MATRIX3(A, B, C);
_CONVERT_COMPLEX2(alpha, beta);
gsl_blas_zhemm(CBLAS_SIDE_val(side), CBLAS_UPLO_val(uplo),
z_alpha, &m_A, &m_B, z_beta, &m_C);
return Val_unit;
}
CAMLprim value ml_gsl_blas_zhemm_bc(value *argv, int argc)
{
return ml_gsl_blas_zhemm(argv[0], argv[1], argv[2],
argv[3], argv[4], argv[5], argv[6]);
}
CAMLprim value ml_gsl_blas_zherk(value uplo, value trans,
value alpha, value A,
value beta, value C)
{
_DECLARE_MATRIX2(A, C);
_CONVERT_MATRIX2(A, C);
gsl_blas_zherk(CBLAS_UPLO_val(uplo), CBLAS_TRANS_val(trans),
Double_val(alpha), &m_A, Double_val(beta), &m_C);
return Val_unit;
}
CAMLprim value ml_gsl_blas_zherk_bc(value *argv, int argc)
{
return ml_gsl_blas_zherk(argv[0], argv[1], argv[2],
argv[3], argv[4], argv[5]);
}
CAMLprim value ml_gsl_blas_zher2k(value uplo, value trans,
value alpha, value A, value B,
value beta, value C)
{
_DECLARE_MATRIX3(A, B, C);
_DECLARE_COMPLEX(alpha);
_CONVERT_MATRIX3(A, B, C);
_CONVERT_COMPLEX(alpha);
gsl_blas_zher2k(CBLAS_UPLO_val(uplo), CBLAS_TRANS_val(trans),
z_alpha, &m_A, &m_B, Double_val(beta), &m_C);
return Val_unit;
}
CAMLprim value ml_gsl_blas_zher2k_bc(value *argv, int argc)
{
return ml_gsl_blas_zher2k(argv[0], argv[1], argv[2],
argv[3], argv[4], argv[5], argv[6]);
}
|