double* get_vector_ptr(VALUE ary, size_t *stride, size_t *n) { gsl_vector *v = NULL; gsl_vector_complex *vc = NULL; gsl_matrix *m; if (VECTOR_P(ary)) { Data_Get_Struct(ary, gsl_vector, v); *stride = v->stride; *n = v->size; return v->data; } else if (VECTOR_COMPLEX_P(ary)) { Data_Get_Struct(ary, gsl_vector_complex, vc); *stride = vc->stride; *n = vc->size*2; return vc->data; } else if (MATRIX_P(ary)) { Data_Get_Struct(ary, gsl_matrix, m); *stride = 1; *n = m->size1*m->size2; return m->data; #ifdef HAVE_NARRAY_H } else if (NA_IsNArray(ary)) { VALUE ary2; *n = NA_TOTAL(ary); *stride = 1; ary2 = na_change_type(ary, NA_DFLOAT); return NA_PTR_TYPE(ary2,double*); #endif #ifdef HAVE_NMATRIX_H } else if (NM_IsNMatrix(ary)) {
static VALUE rb_gsl_math_eval(double (*func)(const double), VALUE xx) { VALUE x, ary; size_t i, size; #ifdef HAVE_NARRAY_H struct NARRAY *na; double *ptr1, *ptr2; #endif if (CLASS_OF(xx) == rb_cRange) xx = rb_gsl_range2ary(xx); switch (TYPE(xx)) { case T_FIXNUM: case T_BIGNUM: case T_FLOAT: return rb_float_new((*func)(NUM2DBL(xx))); break; case T_ARRAY: size = RARRAY_LEN(xx); ary = rb_ary_new2(size); for (i = 0; i < size; i++) { x = rb_ary_entry(xx, i); Need_Float(x); rb_ary_store(ary, i, rb_float_new((*func)(RFLOAT_VALUE(x)))); } return ary; break; default: #ifdef HAVE_NARRAY_H if (NA_IsNArray(xx)) { GetNArray(xx, na); ptr1 = (double*) na->ptr; size = na->total; ary = na_make_object(NA_DFLOAT, na->rank, na->shape, CLASS_OF(xx)); ptr2 = NA_PTR_TYPE(ary, double*); for (i = 0; i < size; i++) ptr2[i] = (*func)(ptr1[i]); return ary; } #endif if (VECTOR_P(xx)) { return vector_eval_create(xx, func); } else if (MATRIX_P(xx)) { return matrix_eval_create(xx, func); } else { rb_raise(rb_eTypeError, "wrong argument type %s (Array or Vector or Matrix expected)", rb_class2name(CLASS_OF(xx))); } break; } /* never reach here */ return Qnil; }
static VALUE rb_gsl_deriv_eval(VALUE obj, VALUE xx, VALUE hh, int (*deriv)(const gsl_function *, double, double, double *, double *)) { gsl_function *f = NULL; double result, abserr, h; VALUE x, ary, aerr; gsl_vector *v = NULL, *vnew = NULL, *verr = NULL; gsl_matrix *m = NULL, *mnew = NULL, *merr = NULL; size_t n, i, j; int status; Need_Float(hh); Data_Get_Struct(obj, gsl_function, f); h = NUM2DBL(hh); if (CLASS_OF(xx) == rb_cRange) xx = rb_gsl_range2ary(xx); switch (TYPE(xx)) { case T_FIXNUM: case T_BIGNUM: case T_FLOAT: status = (*deriv)(f, NUM2DBL(xx), h, &result, &abserr); return rb_ary_new3(3, rb_float_new(result), rb_float_new(abserr), INT2FIX(status)); break; case T_ARRAY: // n = RARRAY(xx)->len; n = RARRAY_LEN(xx); ary = rb_ary_new2(n); aerr = rb_ary_new2(n); for (i = 0; i < n; i++) { x = rb_ary_entry(xx, i); Need_Float(x); (*deriv)(f, NUM2DBL(x), h, &result, &abserr); rb_ary_store(ary, i, rb_float_new(result)); rb_ary_store(aerr, i, rb_float_new(abserr)); } return rb_ary_new3(2, ary, aerr); break; default: #ifdef HAVE_NARRAY_H if (NA_IsNArray(xx)) { struct NARRAY *na; double *ptr1, *ptr2, *ptr3; VALUE ary2, ary3; GetNArray(xx, na); n = na->total; ptr1 = (double*) na->ptr; ary2 = na_make_object(NA_DFLOAT, na->rank, na->shape, CLASS_OF(xx)); ary3 = na_make_object(NA_DFLOAT, na->rank, na->shape, CLASS_OF(xx)); ptr2 = NA_PTR_TYPE(ary2, double*); ptr3 = NA_PTR_TYPE(ary3, double*); for (i = 0; i < n; i++) { (*deriv)(f, ptr1[i], h, &result, &abserr); ptr2[i] = result; ptr3[i] = abserr; } return rb_ary_new3(2, ary2, ary3); } #endif if (VECTOR_P(xx)) { Data_Get_Struct(xx, gsl_vector, v); vnew = gsl_vector_alloc(v->size); verr = gsl_vector_alloc(v->size); for (i = 0; i < v->size; i++) { (*deriv)(f, gsl_vector_get(v, i), h, &result, &abserr); gsl_vector_set(vnew, i, result); gsl_vector_set(verr, i, abserr); } return rb_ary_new3(2, Data_Wrap_Struct(cgsl_vector, 0, gsl_vector_free, vnew), Data_Wrap_Struct(cgsl_vector, 0, gsl_vector_free, verr)); } else if (MATRIX_P(xx)) { Data_Get_Struct(xx, gsl_matrix, m); mnew = gsl_matrix_alloc(m->size1, m->size2); merr = gsl_matrix_alloc(m->size1, m->size2); for (i = 0; i < m->size1; i++) { for (j = 0; j < m->size2; j++) { (*deriv)(f, gsl_matrix_get(m, i, j), h, &result, &abserr); gsl_matrix_set(mnew, i, j, result); gsl_matrix_set(merr, i, j, abserr); } } return rb_ary_new3(2, Data_Wrap_Struct(cgsl_matrix, 0, gsl_matrix_free, mnew), Data_Wrap_Struct(cgsl_matrix, 0, gsl_matrix_free, merr)); } else { rb_raise(rb_eTypeError, "wrong argument type"); } break; } return Qnil; /* never reach here */ }
static VALUE rb_gsl_cheb_eval(VALUE obj, VALUE xx) { gsl_cheb_series *p = NULL; VALUE x, ary; size_t i, j, n; gsl_vector *v = NULL, *vnew = NULL; gsl_matrix *m = NULL, *mnew = NULL; Data_Get_Struct(obj, gsl_cheb_series, p); if (CLASS_OF(xx) == rb_cRange) xx = rb_gsl_range2ary(xx); switch (TYPE(xx)) { case T_FIXNUM: case T_BIGNUM: case T_FLOAT: return rb_float_new(gsl_cheb_eval(p, NUM2DBL(xx))); break; case T_ARRAY: // n = RARRAY(xx)->len; n = RARRAY_LEN(xx); ary = rb_ary_new2(n); for (i = 0; i < n; i++) { x = rb_ary_entry(xx, i); Need_Float(xx); rb_ary_store(ary, i, rb_float_new(gsl_cheb_eval(p, NUM2DBL(x)))); } return ary; break; default: #ifdef HAVE_NARRAY_H if (NA_IsNArray(xx)) { struct NARRAY *na; double *ptr1, *ptr2; GetNArray(xx, na); ptr1 = (double*) na->ptr; n = na->total; ary = na_make_object(NA_DFLOAT, na->rank, na->shape, CLASS_OF(xx)); ptr2 = NA_PTR_TYPE(ary,double*); for (i = 0; i < n; i++) ptr2[i] = gsl_cheb_eval(p, ptr1[i]); return ary; } #endif if (VECTOR_P(xx)) { Data_Get_Struct(xx, gsl_vector, v); vnew = gsl_vector_alloc(v->size); for (i = 0; i < v->size; i++) { gsl_vector_set(vnew, i, gsl_cheb_eval(p, gsl_vector_get(v, i))); } return Data_Wrap_Struct(cgsl_vector, 0, gsl_vector_free, vnew); } else if (MATRIX_P(xx)) { Data_Get_Struct(xx, gsl_matrix, m); mnew = gsl_matrix_alloc(m->size1, m->size2); for (i = 0; i < m->size1; i++) { for (j = 0; j < m->size2; j++) { gsl_matrix_set(mnew, i, j, gsl_cheb_eval(p, gsl_matrix_get(m, i, j))); } } return Data_Wrap_Struct(cgsl_matrix, 0, gsl_matrix_free, mnew); } else { rb_raise(rb_eTypeError, "wrong argument type"); } break; } return Qnil; /* never reach here */ }
static VALUE rb_gsl_cheb_eval_n_err(VALUE obj, VALUE nn, VALUE xx) { gsl_cheb_series *p = NULL; double result, err; VALUE x, ary, aerr; size_t n, order, i, j; gsl_vector *v, *vnew, *verr; gsl_matrix *m, *mnew, *merr; CHECK_FIXNUM(nn); order = FIX2INT(nn); Data_Get_Struct(obj, gsl_cheb_series, p); if (CLASS_OF(xx) == rb_cRange) xx = rb_gsl_range2ary(xx); switch (TYPE(xx)) { case T_FIXNUM: case T_BIGNUM: case T_FLOAT: gsl_cheb_eval_n_err(p, order, NUM2DBL(xx), &result, &err); return rb_ary_new3(2, rb_float_new(result), rb_float_new(err)); break; case T_ARRAY: // n = RARRAY(xx)->len; n = RARRAY_LEN(xx); ary = rb_ary_new2(n); aerr = rb_ary_new2(n); for (i = 0; i < n; i++) { x = rb_ary_entry(xx, i); Need_Float(xx); gsl_cheb_eval_n_err(p, order, NUM2DBL(x), &result, &err); rb_ary_store(ary, i, rb_float_new(result)); rb_ary_store(aerr, i, rb_float_new(err)); } return rb_ary_new3(2, ary, aerr); break; default: #ifdef HAVE_NARRAY_H if (NA_IsNArray(xx)) { struct NARRAY *na; double *ptr1, *ptr2, *ptr3; GetNArray(xx, na); ptr1 = (double*) na->ptr; n = na->total; ary = na_make_object(NA_DFLOAT, na->rank, na->shape, CLASS_OF(xx)); aerr = na_make_object(NA_DFLOAT, na->rank, na->shape, CLASS_OF(xx)); ptr2 = NA_PTR_TYPE(ary,double*); ptr3 = NA_PTR_TYPE(aerr,double*); for (i = 0; i < n; i++) { gsl_cheb_eval_n_err(p, order, ptr1[i], &result, &err); ptr2[i] = result; ptr3[i] = err; } return rb_ary_new3(2, ary, aerr); } #endif if (VECTOR_P(xx)) { Data_Get_Struct(xx, gsl_vector, v); vnew = gsl_vector_alloc(v->size); verr = gsl_vector_alloc(v->size); for (i = 0; i < v->size; i++) { gsl_cheb_eval_n_err(p, order, gsl_vector_get(v, i), &result, &err); gsl_vector_set(vnew, i, result); gsl_vector_set(verr, i, err); } return rb_ary_new3(2, Data_Wrap_Struct(cgsl_vector, 0, gsl_vector_free, vnew), Data_Wrap_Struct(cgsl_vector, 0, gsl_vector_free, verr)); } else if (MATRIX_P(xx)) { Data_Get_Struct(xx, gsl_matrix, m); mnew = gsl_matrix_alloc(m->size1, m->size2); merr = gsl_matrix_alloc(m->size1, m->size2); for (i = 0; i < m->size1; i++) { for (j = 0; j < m->size2; j++) { gsl_cheb_eval_n_err(p, order, gsl_matrix_get(m, i, j), &result, &err); gsl_matrix_set(mnew, i, j, result); gsl_matrix_set(merr, i, j, err); } } return rb_ary_new3(2, Data_Wrap_Struct(cgsl_matrix, 0, gsl_matrix_free, mnew), Data_Wrap_Struct(cgsl_matrix, 0, gsl_matrix_free, merr)); } else { rb_raise(rb_eTypeError, "wrong argument type"); } break; } return Qnil; /* never reach here */ }
static VALUE rb_gsl_interp_evaluate(VALUE obj, VALUE xxa, VALUE yya, VALUE xx, double (*eval)(const gsl_interp *, const double [], const double [], double, gsl_interp_accel *)) { rb_gsl_interp *rgi = NULL; double *ptrx = NULL, *ptry = NULL; gsl_vector *v = NULL, *vnew = NULL; gsl_matrix *m = NULL, *mnew = NULL; VALUE ary, x; double val; size_t n, i, j, size, stridex, stridey; #ifdef HAVE_NARRAY_H struct NARRAY *na = NULL; double *ptrz = NULL, *ptr = NULL; #endif Data_Get_Struct(obj, rb_gsl_interp, rgi); ptrx = get_vector_ptr(xxa, &stridex, &size); if (size != rgi->p->size ){ rb_raise(rb_eTypeError, "size mismatch (xa:%d != %d)", (int) size, (int) rgi->p->size); } ptry = get_vector_ptr(yya, &stridey, &size); if (size != rgi->p->size ){ rb_raise(rb_eTypeError, "size mismatch (ya:%d != %d)", (int) size, (int) rgi->p->size); } if (CLASS_OF(xx) == rb_cRange) xx = rb_gsl_range2ary(xx); switch (TYPE(xx)) { case T_FIXNUM: case T_BIGNUM: case T_FLOAT: Need_Float(xx); return rb_float_new((*eval)(rgi->p, ptrx, ptry, NUM2DBL(xx), rgi->a)); break; case T_ARRAY: n = RARRAY(xx)->len; ary = rb_ary_new2(n); for (i = 0; i < n; i++) { x = rb_ary_entry(xx, i); Need_Float(x); val = (*eval)(rgi->p, ptrx, ptry, NUM2DBL(x), rgi->a); rb_ary_store(ary, i, rb_float_new(val)); } return ary; break; default: #ifdef HAVE_NARRAY_H if (NA_IsNArray(xx)) { GetNArray(xx, na); ptrz = (double*) na->ptr; ary = na_make_object(NA_DFLOAT, na->rank, na->shape, CLASS_OF(xx)); ptr = NA_PTR_TYPE(ary, double*); for (i = 0; i < na->total; i++) ptr[i] = (*eval)(rgi->p, ptrx, ptry, ptrz[i], rgi->a); return ary; } #endif if (VECTOR_P(xx)) { Data_Get_Struct(xx, gsl_vector, v); vnew = gsl_vector_alloc(v->size); for (i = 0; i < v->size; i++) { val = (*eval)(rgi->p, ptrx, ptry, gsl_vector_get(v, i), rgi->a); gsl_vector_set(vnew, i, val); } return Data_Wrap_Struct(cgsl_vector, 0, gsl_vector_free, vnew); } else if (MATRIX_P(xx)) { Data_Get_Struct(xx, gsl_matrix, m); mnew = gsl_matrix_alloc(m->size1, m->size2); for (i = 0; i < m->size1; i++) { for (j = 0; j < m->size2; j++) { val = (*eval)(rgi->p, ptrx, ptry, gsl_matrix_get(m, i, j), rgi->a); gsl_matrix_set(mnew, i, j, val); } } return Data_Wrap_Struct(cgsl_matrix, 0, gsl_matrix_free, mnew); } else { rb_raise(rb_eTypeError, "wrong argument type %s", rb_class2name(CLASS_OF(xx))); } break; } /* never reach here */ return Qnil; }
static VALUE rb_gsl_pow_int(VALUE obj, VALUE xx, VALUE nn) { VALUE x, ary, argv[2]; size_t i, j, size; int n; gsl_vector *v = NULL, *vnew = NULL; gsl_matrix *m = NULL, *mnew = NULL; if (CLASS_OF(xx) == rb_cRange) xx = rb_gsl_range2ary(xx); switch (TYPE(xx)) { case T_FIXNUM: case T_BIGNUM: case T_FLOAT: return rb_float_new(gsl_pow_int(NUM2DBL(xx), FIX2INT(nn))); break; case T_ARRAY: CHECK_FIXNUM(nn); n = FIX2INT(nn); size = RARRAY_LEN(xx); ary = rb_ary_new2(size); for (i = 0; i < size; i++) { x = rb_ary_entry(xx, i); Need_Float(x); // rb_ary_store(ary, i, rb_float_new(gsl_pow_int(RFLOAT(x)->value, n))); rb_ary_store(ary, i, rb_float_new(gsl_pow_int(NUM2DBL(x), n))); } return ary; break; default: #ifdef HAVE_NARRAY_H if (NA_IsNArray(xx)) { struct NARRAY *na; double *ptr1, *ptr2; CHECK_FIXNUM(nn); n = FIX2INT(nn); GetNArray(xx, na); ptr1 = (double*) na->ptr; size = na->total; ary = na_make_object(NA_DFLOAT, na->rank, na->shape, CLASS_OF(xx)); ptr2 = NA_PTR_TYPE(ary, double*); for (i = 0; i < size; i++) ptr2[i] = gsl_pow_int(ptr1[i], n); return ary; } #endif if (VECTOR_P(xx)) { CHECK_FIXNUM(nn); n = FIX2INT(nn); Data_Get_Struct(xx, gsl_vector, v); vnew = gsl_vector_alloc(v->size); for (i = 0; i < v->size; i++) { gsl_vector_set(vnew, i, gsl_pow_int(gsl_vector_get(v, i), n)); } return Data_Wrap_Struct(cgsl_vector, 0, gsl_vector_free, vnew); } else if (MATRIX_P(xx)) { CHECK_FIXNUM(nn); n = FIX2INT(nn); Data_Get_Struct(xx, gsl_matrix, m); mnew = gsl_matrix_alloc(m->size1, m->size2); for (i = 0; i < m->size1; i++) { for (j = 0; j < m->size2; j++) { gsl_matrix_set(mnew, i, j, gsl_pow_int(gsl_matrix_get(m, i, j), n)); } } return Data_Wrap_Struct(cgsl_matrix, 0, gsl_matrix_free, mnew); } else if (COMPLEX_P(xx) || VECTOR_COMPLEX_P(xx) || MATRIX_COMPLEX_P(xx)) { argv[0] = xx; argv[1] = nn; return rb_gsl_complex_pow_real(2, argv, obj); } else { rb_raise(rb_eTypeError, "wrong argument type %s (Array or Vector or Matrix expected)", rb_class2name(CLASS_OF(xx))); } break; } /* never reach here */ return Qnil; }
static VALUE rb_gsl_math_eval2(double (*func)(const double, const double), VALUE xx, VALUE yy) { VALUE x, y, ary; size_t i, j, size; gsl_vector *v = NULL, *v2 = NULL, *vnew = NULL; gsl_matrix *m = NULL, *m2 = NULL, *mnew = NULL; if (CLASS_OF(xx) == rb_cRange) xx = rb_gsl_range2ary(xx); switch (TYPE(xx)) { case T_FIXNUM: case T_BIGNUM: case T_FLOAT: Need_Float(yy); return rb_float_new((*func)(NUM2DBL(xx), NUM2DBL(yy))); break; case T_ARRAY: Check_Type(yy, T_ARRAY); size = RARRAY_LEN(xx); // if (size != RARRAY(yy)->len) rb_raise(rb_eRuntimeError, "array sizes are different."); if ((int) size != RARRAY_LEN(yy)) rb_raise(rb_eRuntimeError, "array sizes are different."); ary = rb_ary_new2(size); for (i = 0; i < size; i++) { x = rb_ary_entry(xx, i); y = rb_ary_entry(yy, i); Need_Float(x); Need_Float(y); // rb_ary_store(ary, i, rb_float_new((*func)(RFLOAT(x)->value, RFLOAT(y)->value))); rb_ary_store(ary, i, rb_float_new((*func)(NUM2DBL(x), NUM2DBL(y)))); } return ary; break; default: #ifdef HAVE_NARRAY_H if (NA_IsNArray(xx)) { struct NARRAY *nax, *nay; double *ptr1, *ptr2, *ptr3; GetNArray(xx, nax); GetNArray(yy, nay); ptr1 = (double*) nax->ptr; ptr2 = (double*) nay->ptr; size = nax->total; ary = na_make_object(NA_DFLOAT, nax->rank, nax->shape, CLASS_OF(xx)); ptr3 = NA_PTR_TYPE(ary, double*); for (i = 0; i < size; i++) ptr3[i] = (*func)(ptr1[i], ptr2[i]); return ary; } #endif if (VECTOR_P(xx)) { CHECK_VECTOR(yy); Data_Get_Struct(xx, gsl_vector, v); Data_Get_Struct(yy, gsl_vector, v2); vnew = gsl_vector_alloc(v->size); for (i = 0; i < v->size; i++) { gsl_vector_set(vnew, i, (*func)(gsl_vector_get(v, i), gsl_vector_get(v2, i))); } return Data_Wrap_Struct(cgsl_vector, 0, gsl_vector_free, vnew); } else if (MATRIX_P(xx)) { CHECK_MATRIX(yy); Data_Get_Struct(xx, gsl_matrix, m); Data_Get_Struct(yy, gsl_matrix, m2); mnew = gsl_matrix_alloc(m->size1, m->size2); for (i = 0; i < m->size1; i++) { for (j = 0; j < m->size2; j++) { gsl_matrix_set(mnew, i, j, (*func)(gsl_matrix_get(m, i, j), gsl_matrix_get(m2, i, j))); } } return Data_Wrap_Struct(cgsl_matrix, 0, gsl_matrix_free, mnew); } else { rb_raise(rb_eTypeError, "wrong argument type %s " "(Array or Vector or Matrix expected)", rb_class2name(CLASS_OF(xx))); } break; } /* never reach here */ return Qnil; }
static VALUE rb_gsl_wavelet_transform0(int argc, VALUE *argv, VALUE obj, int sss) { gsl_wavelet *w = NULL; gsl_vector *v = NULL, *vnew; gsl_wavelet_direction dir = forward; gsl_wavelet_workspace *work = NULL; int itmp, flag = 0; // local variable "status" declared and set, but never used //int status; double *ptr1, *ptr2; size_t n, stride; int naflag = 0; VALUE ary, ret; #ifdef HAVE_NARRAY_H struct NARRAY *na1 = NULL; #endif switch (TYPE(obj)) { case T_MODULE: case T_CLASS: case T_OBJECT: if (argc < 2) rb_raise(rb_eArgError, "too few arguments"); CHECK_WAVELET(argv[0]); if (MATRIX_P(argv[1])) { return rb_gsl_wavelet2d(argc, argv, obj, gsl_wavelet2d_transform_matrix, sss); } if (VECTOR_P(argv[1])) { Data_Get_Struct(argv[0], gsl_wavelet, w); Data_Get_Struct(argv[1], gsl_vector, v); ret = argv[1]; ptr1 = v->data; n = v->size; stride = v->stride; #ifdef HAVE_NARRAY_H } else if (NA_IsNArray(argv[1])) { GetNArray(argv[1], na1); ret = argv[1]; ptr1 = (double*) na1->ptr; n = na1->total; naflag = 1; stride = 1; #endif } else { rb_raise(rb_eTypeError, "wrong argument type (Vector expected)"); } itmp = 2; break; default: if (argc < 1) rb_raise(rb_eArgError, "too few arguments"); if (MATRIX_P(argv[0])) { return rb_gsl_wavelet2d(argc, argv, obj, gsl_wavelet2d_transform_matrix, sss); } if (VECTOR_P(obj)) { CHECK_WAVELET(argv[0]); Data_Get_Struct(argv[0], gsl_wavelet, w); Data_Get_Struct(obj, gsl_vector, v); ret = obj; ptr1 = v->data; n = v->size; stride = v->stride; } else if (VECTOR_P(argv[0])) { CHECK_WAVELET(obj); Data_Get_Struct(obj, gsl_wavelet, w); Data_Get_Struct(argv[0], gsl_vector, v); ret = argv[0]; ptr1 = v->data; n = v->size; stride = v->stride; #ifdef HAVE_NARRAY_H } else if (NA_IsNArray(obj)) { CHECK_WAVELET(argv[0]); Data_Get_Struct(argv[0], gsl_wavelet, w); GetNArray(obj, na1); ret = obj; ptr1 = (double*) na1->ptr; n = na1->total; naflag = 1; stride = 1; } else if (NA_IsNArray(argv[0])) { CHECK_WAVELET(obj); Data_Get_Struct(obj, gsl_wavelet, w); GetNArray(argv[0], na1); ret = argv[0]; ptr1 = (double*) na1->ptr; n = na1->total; naflag = 1; stride = 1; #endif } else { rb_raise(rb_eTypeError, "wrong argument type"); } itmp = 1; break; } switch (argc - itmp) { case 2: CHECK_FIXNUM(argv[itmp]); CHECK_WORKSPACE(argv[itmp+1]); dir = FIX2INT(argv[itmp]); Data_Get_Struct(argv[itmp+1], gsl_wavelet_workspace, work); break; case 1: if (TYPE(argv[itmp]) == T_FIXNUM) { dir = FIX2INT(argv[itmp]); work = gsl_wavelet_workspace_alloc(v->size); flag = 1; } else if (rb_obj_is_kind_of(argv[itmp], cgsl_wavelet_workspace)) { Data_Get_Struct(argv[itmp], gsl_wavelet_workspace, work); } else { rb_raise(rb_eTypeError, "wrong argument type"); } break; case 0: work = gsl_wavelet_workspace_alloc(v->size); flag = 1; break; default: rb_raise(rb_eArgError, "too many arguments"); break; } if (naflag == 0) { if (sss == RB_GSL_DWT_COPY) { vnew = gsl_vector_alloc(v->size); gsl_vector_memcpy(vnew, v); ary = Data_Wrap_Struct(cgsl_vector, 0, gsl_vector_free, vnew); ptr2 = vnew->data; } else { ary = ret; ptr2 = ptr1; } } else { #ifdef HAVE_NARRAY_H if (sss == RB_GSL_DWT_COPY) { ary = na_make_object(NA_DFLOAT, na1->rank, na1->shape, cNArray); ptr2 = NA_PTR_TYPE(ary, double*); memcpy(ptr2, ptr1, sizeof(double)*n); } else {
/* * Calculates a function at x, and returns the rusult. */ static VALUE rb_gsl_function_eval(VALUE obj, VALUE x) { gsl_function *F = NULL; VALUE ary, proc, params, result, arynew, x2; gsl_vector *v = NULL, *vnew = NULL; gsl_matrix *m = NULL, *mnew = NULL; size_t i, j, n; Data_Get_Struct(obj, gsl_function, F); ary = (VALUE) F->params; proc = rb_ary_entry(ary, 0); params = rb_ary_entry(ary, 1); if (CLASS_OF(x) == rb_cRange) x = rb_gsl_range2ary(x); switch (TYPE(x)) { case T_FIXNUM: case T_BIGNUM: case T_FLOAT: if (NIL_P(params)) result = rb_funcall(proc, RBGSL_ID_call, 1, x); else result = rb_funcall(proc, RBGSL_ID_call, 2, x, params); return result; break; case T_ARRAY: // n = RARRAY(x)->len; n = RARRAY_LEN(x); arynew = rb_ary_new2(n); for (i = 0; i < n; i++) { x2 = rb_ary_entry(x, i); Need_Float(x2); if (NIL_P(params)) result = rb_funcall(proc, RBGSL_ID_call, 1, x2); else result = rb_funcall(proc, RBGSL_ID_call, 2, x2, params); rb_ary_store(arynew, i, result); } return arynew; break; default: #ifdef HAVE_NARRAY_H if (NA_IsNArray(x)) { double *ptr1, *ptr2; struct NARRAY *na; GetNArray(x, na); ptr1 = (double *) na->ptr; n = na->total; ary = na_make_object(NA_DFLOAT, na->rank, na->shape, CLASS_OF(x)); ptr2 = NA_PTR_TYPE(ary, double*); for (i = 0; i < n; i++) { x2 = rb_float_new(ptr1[i]); if (NIL_P(params)) result = rb_funcall(proc, RBGSL_ID_call, 1, x2); else result = rb_funcall(proc, RBGSL_ID_call, 2, x2, params); ptr2[i] = NUM2DBL(result); } return ary; } #endif #ifdef HAVE_NMATRIX_H if (NM_IsNMatrix(x)) { double *ptr1, *ptr2; NM_DENSE_STORAGE *nm; nm = NM_STORAGE_DENSE(x); ptr1 = (double *) nm->elements; n = NM_DENSE_COUNT(x); ary = rb_nmatrix_dense_create(FLOAT64, nm->shape, nm->dim, nm->elements, n); ptr2 = (double*)NM_DENSE_ELEMENTS(ary); for (i = 0; i < n; i++) { x2 = rb_float_new(ptr1[i]); if (NIL_P(params)) result = rb_funcall(proc, RBGSL_ID_call, 1, x2); else result = rb_funcall(proc, RBGSL_ID_call, 2, x2, params); ptr2[i] = NUM2DBL(result); } return ary; } #endif if (VECTOR_P(x)) { Data_Get_Struct(x, gsl_vector, v); vnew = gsl_vector_alloc(v->size); for (i = 0; i < v->size; i++) { x2 = rb_float_new(gsl_vector_get(v, i)); if (NIL_P(params)) result = rb_funcall(proc, RBGSL_ID_call, 1, x2); else result = rb_funcall(proc, RBGSL_ID_call, 2, x2, params); gsl_vector_set(vnew, i, NUM2DBL(result)); } return Data_Wrap_Struct(cgsl_vector, 0, gsl_vector_free, vnew); } else if (MATRIX_P(x)) { Data_Get_Struct(x, gsl_matrix, m); mnew = gsl_matrix_alloc(m->size1, m->size2); for (i = 0; i < m->size1; i++) { for (j = 0; j < m->size2; j++) { x2 = rb_float_new(gsl_matrix_get(m, i, j)); if (NIL_P(params)) result = rb_funcall(proc, RBGSL_ID_call, 1, x2); else result = rb_funcall(proc, RBGSL_ID_call, 2, x2, params); gsl_matrix_set(mnew, i, j, NUM2DBL(result)); } } return Data_Wrap_Struct(cgsl_matrix, 0, gsl_matrix_free, mnew); } else { rb_raise(rb_eTypeError, "wrong argument type"); } break; } /* never reach here */ return Qnil; }
static VALUE rb_gsl_spline_evaluate(VALUE obj, VALUE xx, double (*eval)(const gsl_spline *, double, gsl_interp_accel *)) { rb_gsl_spline *rgs = NULL; gsl_vector *v = NULL, *vnew = NULL; gsl_matrix *m = NULL, *mnew = NULL; VALUE ary, x; double val; size_t n, i, j; #ifdef HAVE_NARRAY_H double *ptr1 = NULL, *ptr2 = NULL; struct NARRAY *na = NULL; #endif Data_Get_Struct(obj, rb_gsl_spline, rgs); if (CLASS_OF(xx) == rb_cRange) xx = rb_gsl_range2ary(xx); switch (TYPE(xx)) { case T_FIXNUM: case T_BIGNUM: case T_FLOAT: Need_Float(xx); return rb_float_new((*eval)(rgs->s, NUM2DBL(xx), rgs->a)); break; case T_ARRAY: n = RARRAY_LEN(xx); ary = rb_ary_new2(n); for (i = 0; i < n; i++) { x = rb_ary_entry(xx, i); Need_Float(x); val = (*eval)(rgs->s, NUM2DBL(x), rgs->a); rb_ary_store(ary, i, rb_float_new(val)); } return ary; break; default: #ifdef HAVE_NARRAY_H if (NA_IsNArray(xx)) { GetNArray(xx, na); ptr1 = (double *) na->ptr; n = na->total; ary = na_make_object(NA_DFLOAT, na->rank, na->shape, CLASS_OF(xx)); ptr2 = NA_PTR_TYPE(ary, double*); for (i = 0; i < n; i++) ptr2[i] = (*eval)(rgs->s, ptr1[i], rgs->a); return ary; } #endif if (VECTOR_P(xx)) { Data_Get_Struct(xx, gsl_vector, v); vnew = gsl_vector_alloc(v->size); for (i = 0; i < v->size; i++) { val = (*eval)(rgs->s, gsl_vector_get(v, i), rgs->a); gsl_vector_set(vnew, i, val); } return Data_Wrap_Struct(cgsl_vector, 0, gsl_vector_free, vnew); } else if (MATRIX_P(xx)) { Data_Get_Struct(xx, gsl_matrix, m); mnew = gsl_matrix_alloc(m->size1, m->size2); for (i = 0; i < m->size1; i++) { for (j = 0; j < m->size2; j++) { val = (*eval)(rgs->s, gsl_matrix_get(m, i, j), rgs->a); gsl_matrix_set(mnew, i, j, val); } } return Data_Wrap_Struct(cgsl_matrix, 0, gsl_matrix_free, mnew); } else { rb_raise(rb_eTypeError, "wrong argument type %s", rb_class2name(CLASS_OF(xx))); } break; } /* never reach here */ return Qnil; }