/* This returns 7 elements array */ static VALUE rb_gsl_fit_linear(int argc, VALUE *argv, VALUE obj) { double *ptrx, *ptry; double c0, c1, cov00, cov01, cov11, sumsq; int status; size_t n, stridex, stridey; switch (argc) { case 2: ptrx = get_vector_ptr(argv[0], &stridex, &n); ptry = get_vector_ptr(argv[1], &stridey, &n); break; case 3: CHECK_FIXNUM(argv[2]); ptrx = get_vector_ptr(argv[0], &stridex, &n); ptry = get_vector_ptr(argv[1], &stridey, &n); n = FIX2INT(argv[2]); break; default: rb_raise(rb_eArgError, "wrong number of arguments (%d for 2 or 3)", argc); break; } status = gsl_fit_linear(ptrx, stridex, ptry, stridey, n, &c0, &c1, &cov00, &cov01, &cov11, &sumsq); return rb_ary_new3(7, rb_float_new(c0), rb_float_new(c1), rb_float_new(cov00), rb_float_new(cov01), rb_float_new(cov11), rb_float_new(sumsq), INT2FIX(status)); }
static VALUE rb_gsl_interp_eval_integ_e(VALUE obj, VALUE xxa, VALUE yya, VALUE aa, VALUE bb) { rb_gsl_interp *rgi = NULL; double *ptr1 = NULL, *ptr2 = NULL; size_t size, stridex, stridey; double y, a, b; int status; Need_Float(aa); Need_Float(bb); Data_Get_Struct(obj, rb_gsl_interp, rgi); ptr1 = get_vector_ptr(xxa, &stridex, &size); ptr2 = get_vector_ptr(yya, &stridey, &size); a = NUM2DBL(aa); b = NUM2DBL(bb); status = gsl_interp_eval_integ_e(rgi->p, ptr1, ptr2, a, b, rgi->a, &y); switch (status) { case GSL_EDOM: rb_gsl_error_handler("gsl_interp_eval_integ_e error", __FILE__, __LINE__, status); break; default: return rb_float_new(y); break; } return Qnil; }
static VALUE rb_gsl_interp_new(int argc, VALUE *argv, VALUE klass) { rb_gsl_interp *sp = NULL; const gsl_interp_type *T = NULL; double *ptrx = NULL, *ptry = NULL; size_t sizex = 0, sizey = 0, size = 0, stride = 1; int i; for (i = 0; i < argc; i++) { switch (TYPE(argv[i])) { case T_STRING: T = get_interp_type(argv[i]); break; case T_FIXNUM: if (T) size = FIX2INT(argv[i]); else T = get_interp_type(argv[i]); break; default: if (ptrx == NULL) { ptrx = get_vector_ptr(argv[i], &stride, &sizex); } else { ptry = get_vector_ptr(argv[i], &stride, &sizey); size = GSL_MIN_INT(sizex, sizey); } break; } } if (size == 0) rb_raise(rb_eRuntimeError, "interp size is not given."); sp = ALLOC(rb_gsl_interp); if (T == NULL) T = gsl_interp_cspline; sp->p = gsl_interp_alloc(T, size); sp->a = gsl_interp_accel_alloc(); if (ptrx && ptry) gsl_interp_init(sp->p, ptrx, ptry, size); return Data_Wrap_Struct(klass, 0, rb_gsl_interp_free, sp); }
static VALUE rb_gsl_interp_init(VALUE obj, VALUE xxa, VALUE yya) { rb_gsl_interp *rgi = NULL; double *ptrx = NULL, *ptry = NULL; size_t size, stride; ptrx = get_vector_ptr(xxa, &stride, &size); ptry = get_vector_ptr(yya, &stride, &size); Data_Get_Struct(obj, rb_gsl_interp, rgi); gsl_interp_init(rgi->p, ptrx, ptry, size); return obj; }
static VALUE rb_gsl_interp_eval_integ(VALUE obj, VALUE xxa, VALUE yya, VALUE aa, VALUE bb) { rb_gsl_interp *rgi = NULL; double *ptr1 = NULL, *ptr2 = NULL; size_t size, stridex, stridey; double a, b; Need_Float(aa); Need_Float(bb); Data_Get_Struct(obj, rb_gsl_interp, rgi); ptr1 = get_vector_ptr(xxa, &stridex, &size); ptr2 = get_vector_ptr(yya, &stridey, &size); a = NUM2DBL(aa); b = NUM2DBL(bb); return rb_float_new(gsl_interp_eval_integ(rgi->p, ptr1, ptr2, a, b, rgi->a)); }
static VALUE rb_gsl_spline_find(VALUE obj, VALUE vv, VALUE xx) { rb_gsl_spline *sp = NULL; double *ptr = NULL, x; size_t size, stride; Data_Get_Struct(obj, rb_gsl_spline, sp); ptr = get_vector_ptr(vv, &stride, &size); x = RFLOAT_VALUE(xx); return INT2FIX(gsl_interp_accel_find(sp->a, ptr, size, x)); }
static VALUE rb_gsl_interp_find(VALUE obj, VALUE vv, VALUE xx) { rb_gsl_interp *rgi = NULL; double *ptr = NULL, x; size_t size, stride; Need_Float(xx); Data_Get_Struct(obj, rb_gsl_interp, rgi); ptr = get_vector_ptr(vv, &stride, &size); x = NUM2DBL(xx); return INT2FIX(gsl_interp_accel_find(rgi->a, ptr, size, x)); }
static VALUE rb_gsl_interp_accel_find(VALUE obj, VALUE vv, VALUE xx) { gsl_interp_accel *a = NULL; double x, *ptr = NULL; size_t size, stride; Need_Float(xx); Data_Get_Struct(obj, gsl_interp_accel, a); ptr = get_vector_ptr(vv, &stride, &size); Need_Float(xx); x = NUM2DBL(xx); return INT2FIX(gsl_interp_accel_find(a, ptr, size, x)); }
VALUE vector_eval_create(VALUE obj, double (*func)(double)) { gsl_vector *vnew; size_t i, size, stride; double *ptr; ptr = get_vector_ptr(obj, &stride, &size); vnew = gsl_vector_alloc(size); for (i = 0; i < size; i++) { gsl_vector_set(vnew, i, (*func)(ptr[i*stride])); } return Data_Wrap_Struct(cgsl_vector, 0, gsl_vector_free, vnew); }
static VALUE rb_gsl_utrunc_accel(VALUE obj) { gsl_sum_levin_utrunc_workspace *w = NULL; double sum, err, sum_plain, *ptr; size_t terms_used, n, stride; ptr = get_vector_ptr(obj, &stride, &n); w = gsl_sum_levin_utrunc_alloc(n); gsl_sum_levin_utrunc_accel(ptr, n, w, &sum, &err); sum_plain = w->sum_plain; terms_used = w->terms_used; gsl_sum_levin_utrunc_free(w); return rb_ary_new3(4, rb_float_new(sum), rb_float_new(err), rb_float_new(sum_plain), INT2FIX(terms_used)); }
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; }