summaryrefslogtreecommitdiff
path: root/lib/GSL
diff options
context:
space:
mode:
authorAlberto Ruiz <aruiz@um.es>2007-09-11 17:34:24 +0000
committerAlberto Ruiz <aruiz@um.es>2007-09-11 17:34:24 +0000
commitcd937c2be2900b8f13506d9ae7c731ad43d74e05 (patch)
tree8a4161c6510aac58b3ab85041145fb19ea6c5615 /lib/GSL
parent834b4837799611fd7fbaa9609ea587e041cb0ca1 (diff)
allow setting off GSL default error handler
Diffstat (limited to 'lib/GSL')
-rw-r--r--lib/GSL/Differentiation.hs2
-rw-r--r--lib/GSL/Integration.hs2
-rw-r--r--lib/GSL/Special.hs52
-rw-r--r--lib/GSL/Special/Internal.hs10
-rw-r--r--lib/GSL/gsl-aux.c5
-rw-r--r--lib/GSL/gsl-aux.h9
6 files changed, 18 insertions, 62 deletions
diff --git a/lib/GSL/Differentiation.hs b/lib/GSL/Differentiation.hs
index 6fd7565..e8e22d2 100644
--- a/lib/GSL/Differentiation.hs
+++ b/lib/GSL/Differentiation.hs
@@ -23,7 +23,7 @@ module GSL.Differentiation (
23) where 23) where
24 24
25import Foreign 25import Foreign
26import Data.Packed.Internal.Common(mkfun,check,(//)) 26import Data.Packed.Internal(mkfun,check,(//))
27 27
28derivGen :: 28derivGen ::
29 Int -- ^ type: 0 central, 1 forward, 2 backward 29 Int -- ^ type: 0 central, 1 forward, 2 backward
diff --git a/lib/GSL/Integration.hs b/lib/GSL/Integration.hs
index dd76be8..4152c9d 100644
--- a/lib/GSL/Integration.hs
+++ b/lib/GSL/Integration.hs
@@ -21,7 +21,7 @@ module GSL.Integration (
21) where 21) where
22 22
23import Foreign 23import Foreign
24import Data.Packed.Internal.Common(mkfun,check,(//)) 24import Data.Packed.Internal(mkfun,check,(//))
25 25
26 26
27-------------------------------------------------------------------- 27--------------------------------------------------------------------
diff --git a/lib/GSL/Special.hs b/lib/GSL/Special.hs
index 9eaf733..bdee1b2 100644
--- a/lib/GSL/Special.hs
+++ b/lib/GSL/Special.hs
@@ -41,6 +41,7 @@ module GSL.Special (
41, module GSL.Special.Synchrotron 41, module GSL.Special.Synchrotron
42, module GSL.Special.Trig 42, module GSL.Special.Trig
43, module GSL.Special.Zeta 43, module GSL.Special.Zeta
44, setErrorHandlerOff
44) 45)
45where 46where
46 47
@@ -74,51 +75,6 @@ import GSL.Special.Trig
74import GSL.Special.Zeta 75import GSL.Special.Zeta
75 76
76 77
77-------------------- simple functions -------------------------- 78-- | This action removes the GSL default error handler which aborts the program, so that
78 79-- GSL errors can be handled by Haskell (using Control.Exception) and ghci doesn't abort.
79{- | The error function (/gsl_sf_erf/), defined as 2\/ \\sqrt \\pi * \int\_0\^t \\exp -t\^2 dt. 80foreign import ccall "gsl-aux.h no_abort_on_error" setErrorHandlerOff :: IO ()
80
81@> erf 1.5
820.9661051464753108@
83
84-}
85foreign import ccall "gsl-aux.h gsl_sf_erf" erf :: Double -> Double
86
87{- | The Gaussian probability density function (/gsl_sf_erf_Z/), defined as (1\/\\sqrt\{2\\pi\}) \\exp(-x\^2\/2).
88
89>> erf_Z 1.5
90>0.12951759566589172
91
92-}
93foreign import ccall "gsl-aux.h gsl_sf_erf_Z" erf_Z :: Double -> Double
94
95-------------------- functions returning sf_result -------------
96
97{- | The regular cylindrical Bessel function of zeroth order, J_0(x). This is
98 the example in the GSL manual, returning the value of the function and
99 the error term:
100
101@\> bessel_J0_e 5.0
102(-0.1775967713143383,1.9302109579684196e-16)@
103
104-}
105bessel_J0_e :: Double -> (Double,Double)
106bessel_J0_e x = createSFR "bessel_J0_e" (gsl_sf_bessel_J0_e x)
107foreign import ccall "gsl-aux.h gsl_sf_bessel_J0_e" gsl_sf_bessel_J0_e :: Double -> Ptr Double -> IO Int
108
109-------------------- functions returning sf_result_e10 -------------
110
111{- | (From the GSL manual) \"This function computes the exponential \exp(x) using the @gsl_sf_result_e10@ type to return a result with extended range. This function may be useful if the value of \exp(x) would overflow the numeric range of double\".
112
113For example:
114
115@\> exp_e10_e 30.0
116(1.0686474581524432,1.4711818964275088e-14,13)@
117
118@\> exp 30.0
1191.0686474581524463e13@
120
121-}
122exp_e10_e :: Double -> (Double,Int,Double)
123exp_e10_e x = createSFR_E10 "exp_e10_e" (gsl_sf_exp_e10_e x)
124foreign import ccall "gsl-aux.h gsl_sf_exp_e10_e" gsl_sf_exp_e10_e :: Double -> Ptr Double -> IO Int
diff --git a/lib/GSL/Special/Internal.hs b/lib/GSL/Special/Internal.hs
index c7455d9..ad1aaa2 100644
--- a/lib/GSL/Special/Internal.hs
+++ b/lib/GSL/Special/Internal.hs
@@ -25,7 +25,7 @@ module GSL.Special.Internal (
25where 25where
26 26
27import Foreign 27import Foreign
28import Data.Packed.Internal.Common(check,(//)) 28import Data.Packed.Internal(check,(//))
29 29
30 30
31data Precision = PrecDouble | PrecSingle | PrecApprox 31data Precision = PrecDouble | PrecSingle | PrecApprox
@@ -39,10 +39,10 @@ type Gsl_mode_t = Int
39 39
40---------------------------------------------------------------- 40----------------------------------------------------------------
41-- | access to a sf_result 41-- | access to a sf_result
42createSFR :: Storable a => t -> (Ptr a -> IO Int) -> (a, a) 42createSFR :: Storable a => String -> (Ptr a -> IO Int) -> (a, a)
43createSFR s f = unsafePerformIO $ do 43createSFR s f = unsafePerformIO $ do
44 p <- mallocArray 2 44 p <- mallocArray 2
45 f p // check "createSFR" [] 45 f p // check s []
46 [val,err] <- peekArray 2 p 46 [val,err] <- peekArray 2 p
47 free p 47 free p
48 return (val,err) 48 return (val,err)
@@ -52,12 +52,12 @@ createSFR s f = unsafePerformIO $ do
52-- the sf_result_e10 contains two doubles and the exponent 52-- the sf_result_e10 contains two doubles and the exponent
53 53
54-- | acces to sf_result_e10 54-- | acces to sf_result_e10
55createSFR_E10 :: (Storable t2, Storable t3, Storable t1) => t -> (Ptr a -> IO Int) -> (t1, t2, t3) 55createSFR_E10 :: (Storable t2, Storable t3, Storable t1) => String -> (Ptr a -> IO Int) -> (t1, t2, t3)
56createSFR_E10 s f = unsafePerformIO $ do 56createSFR_E10 s f = unsafePerformIO $ do
57 let sd = sizeOf (0::Double) 57 let sd = sizeOf (0::Double)
58 let si = sizeOf (0::Int) 58 let si = sizeOf (0::Int)
59 p <- mallocBytes (2*sd + si) 59 p <- mallocBytes (2*sd + si)
60 f p // check "createSFR_E10" [] 60 f p // check s []
61 val <- peekByteOff p 0 61 val <- peekByteOff p 0
62 err <- peekByteOff p sd 62 err <- peekByteOff p sd
63 expo <- peekByteOff p (2*sd) 63 expo <- peekByteOff p (2*sd)
diff --git a/lib/GSL/gsl-aux.c b/lib/GSL/gsl-aux.c
index fc95e6f..eec3e5f 100644
--- a/lib/GSL/gsl-aux.c
+++ b/lib/GSL/gsl-aux.c
@@ -64,6 +64,11 @@
64#define BAD_FILE 1003 64#define BAD_FILE 1003
65 65
66 66
67void no_abort_on_error() {
68 gsl_set_error_handler_off();
69}
70
71
67int toScalarR(int code, KRVEC(x), RVEC(r)) { 72int toScalarR(int code, KRVEC(x), RVEC(r)) {
68 REQUIRES(rn==1,BAD_SIZE); 73 REQUIRES(rn==1,BAD_SIZE);
69 DEBUGMSG("toScalarR"); 74 DEBUGMSG("toScalarR");
diff --git a/lib/GSL/gsl-aux.h b/lib/GSL/gsl-aux.h
index 2be07c1..8f8618a 100644
--- a/lib/GSL/gsl-aux.h
+++ b/lib/GSL/gsl-aux.h
@@ -10,6 +10,8 @@
10#define KCVEC(A) int A##n, const gsl_complex*A##p 10#define KCVEC(A) int A##n, const gsl_complex*A##p
11#define KCMAT(A) int A##r, int A##c, const gsl_complex* A##p 11#define KCMAT(A) int A##r, int A##c, const gsl_complex* A##p
12 12
13void no_abort_on_error();
14
13int toScalarR(int code, KRVEC(x), RVEC(r)); 15int toScalarR(int code, KRVEC(x), RVEC(r));
14/* norm2, absdif, maximum, posmax, etc. */ 16/* norm2, absdif, maximum, posmax, etc. */
15 17
@@ -58,10 +60,3 @@ int minimizeWithDeriv(double f(int, double*), void df(int, double*, double*),
58 KRVEC(xi), RMAT(sol)); 60 KRVEC(xi), RMAT(sol));
59 61
60int deriv(int code, double f(double, void*), double x, double h, double * result, double * abserr); 62int deriv(int code, double f(double, void*), double x, double h, double * result, double * abserr);
61
62double gsl_sf_erf(double);
63double gsl_sf_erf_Z(double);
64double gsl_sf_gamma(double);
65
66int gsl_sf_bessel_J0_e(double, double*); // hmmm...
67int gsl_sf_exp_e10_e(double, double*); // HMMMMM...