Annotation of OpenXM/src/ox_gsl/call_gsl.c, Revision 1.3
1.3 ! takayama 1: /* $OpenXM: OpenXM/src/ox_gsl/call_gsl.c,v 1.2 2018/03/30 04:43:16 takayama Exp $
1.2 takayama 2: */
1.1 takayama 3: //#include <gsl/gsl_types.h>
4: //#include <gsl/gsl_sys.h>
1.2 takayama 5: #include <unistd.h>
1.1 takayama 6: #include <gsl/gsl_sf_result.h>
1.2 takayama 7: #include <gsl/gsl_errno.h>
1.1 takayama 8: #include <gsl/gsl_sf_gamma.h>
1.3 ! takayama 9: #include <gsl/gsl_integration.h>
1.1 takayama 10: #include "ox_gsl.h"
11: extern int Debug;
1.2 takayama 12: // local prototype declarations
13:
1.1 takayama 14: void call_gsl_sf_lngamma_complex_e() {
1.3 ! takayama 15: int argc;
1.1 takayama 16: double zr;
17: double zi;
18: gsl_sf_result lnr;
19: gsl_sf_result arg;
20: int status;
21: cmo *r[3];
22: cmo *ans;
1.2 takayama 23: // gsl_set_error_handler_off();
24: gsl_set_error_handler((gsl_error_handler_t *)myhandler);
1.3 ! takayama 25: argc = get_i(); // number of args
! 26: if (argc != 2) {
! 27: pops(argc);
! 28: push(make_error2("The argc must be 2 for gsl_sf_lngamma_complex_e.",NULL,0,-1));
! 29: return;
! 30: }
1.1 takayama 31: zr=get_double();
32: zi=get_double();
33: if (Debug) printf("gsl_sf_lngamma_complex_e(zr=%lg,zi=%lg)\n",zr,zi);
34: status = gsl_sf_lngamma_complex_e(zr,zi,&lnr,&arg);
35: r[0] = (cmo *)new_cmo_double(lnr.val);
36: r[1] = (cmo *)new_cmo_double(arg.val);
37: r[2] = (cmo *)new_cmo_int32(status);
1.2 takayama 38: ans = (cmo *)new_cmo_list_array((void *)r,3);
1.1 takayama 39: push(ans);
40: }
1.3 ! takayama 41:
! 42: cmo *Func_x=NULL;
! 43: double f_x(double x,void *params) {
! 44: double d;
! 45: replace(1,"x",x);
! 46: if (eval_cmo(Func_x,&d)==0) GSL_ERROR("eval_cmo fails in f_x",GSL_ETOL);
! 47: return(d);
! 48: }
! 49: void call_gsl_integration_qags() {
! 50: int argc;
! 51: double a;
! 52: double b;
! 53: int status;
! 54: cmo *r[3];
! 55: cmo *ans;
! 56: gsl_integration_workspace * w
! 57: = gsl_integration_workspace_alloc (1000);
! 58: double result, error;
! 59: gsl_function F;
! 60: double epsabs=0;
! 61: double epsrel=1e-7;
! 62: int limit=1000;
! 63:
! 64: // gsl_set_error_handler_off();
! 65: gsl_set_error_handler((gsl_error_handler_t *)myhandler);
! 66: argc = get_i(); // number of args
! 67: if (argc != 3) {
! 68: pops(argc);
! 69: push(make_error2("The argc must be 3 for gsl_integration_qags.",NULL,0,-1));
! 70: return;
! 71: }
! 72: Func_x = pop();
! 73: a = get_double();
! 74: b = get_double();
! 75:
! 76: F.function = &f_x;
! 77:
! 78: status=gsl_integration_qags (&F, a, b, epsabs, epsrel, limit,
! 79: w, &result, &error);
! 80:
! 81: // printf ("result = % .18f\n", result);
! 82: // printf ("estimated error = % .18f\n", error);
! 83: // printf ("intervals = %zu\n", w->size);
! 84:
! 85: gsl_integration_workspace_free(w);
! 86:
! 87: r[0] = (cmo *)new_cmo_double(result);
! 88: r[1] = (cmo *)new_cmo_double(error);
! 89: r[2] = (cmo *)new_cmo_int32(status);
! 90: ans = (cmo *)new_cmo_list_array((void *)r,3);
! 91: push(ans);
! 92: }
FreeBSD-CVSweb <freebsd-cvsweb@FreeBSD.org>