diff options
Diffstat (limited to 'lib/Numeric/GSL')
-rw-r--r-- | lib/Numeric/GSL/Internal.hs | 3 | ||||
-rw-r--r-- | lib/Numeric/GSL/Root.hs | 81 | ||||
-rw-r--r-- | lib/Numeric/GSL/gsl-aux.c | 154 |
3 files changed, 220 insertions, 18 deletions
diff --git a/lib/Numeric/GSL/Internal.hs b/lib/Numeric/GSL/Internal.hs index 84417ce..69a9750 100644 --- a/lib/Numeric/GSL/Internal.hs +++ b/lib/Numeric/GSL/Internal.hs | |||
@@ -36,6 +36,9 @@ foreign import ccall safe "wrapper" | |||
36 | foreign import ccall safe "wrapper" | 36 | foreign import ccall safe "wrapper" |
37 | mkDoubleVecVecfun :: (Double -> TVV) -> IO (FunPtr (Double -> TVV)) | 37 | mkDoubleVecVecfun :: (Double -> TVV) -> IO (FunPtr (Double -> TVV)) |
38 | 38 | ||
39 | foreign import ccall safe "wrapper" | ||
40 | mkDoublefun :: (Double -> Double) -> IO (FunPtr (Double -> Double)) | ||
41 | |||
39 | aux_vTov :: (Vector Double -> Vector Double) -> TVV | 42 | aux_vTov :: (Vector Double -> Vector Double) -> TVV |
40 | aux_vTov f n p nr r = g where | 43 | aux_vTov f n p nr r = g where |
41 | v = f x | 44 | v = f x |
diff --git a/lib/Numeric/GSL/Root.hs b/lib/Numeric/GSL/Root.hs index cd2982a..6da15e5 100644 --- a/lib/Numeric/GSL/Root.hs +++ b/lib/Numeric/GSL/Root.hs | |||
@@ -45,6 +45,8 @@ main = do | |||
45 | ----------------------------------------------------------------------------- | 45 | ----------------------------------------------------------------------------- |
46 | 46 | ||
47 | module Numeric.GSL.Root ( | 47 | module Numeric.GSL.Root ( |
48 | uniRoot, UniRootMethod(..), | ||
49 | uniRootJ, UniRootMethodJ(..), | ||
48 | root, RootMethod(..), | 50 | root, RootMethod(..), |
49 | rootJ, RootMethodJ(..), | 51 | rootJ, RootMethodJ(..), |
50 | ) where | 52 | ) where |
@@ -58,6 +60,68 @@ import System.IO.Unsafe(unsafePerformIO) | |||
58 | 60 | ||
59 | ------------------------------------------------------------------------- | 61 | ------------------------------------------------------------------------- |
60 | 62 | ||
63 | data UniRootMethod = Bisection | ||
64 | | FalsePos | ||
65 | | Brent | ||
66 | deriving (Enum, Eq, Show, Bounded) | ||
67 | |||
68 | uniRoot :: UniRootMethod | ||
69 | -> Double | ||
70 | -> Int | ||
71 | -> (Double -> Double) | ||
72 | -> Double | ||
73 | -> Double | ||
74 | -> (Double, Matrix Double) | ||
75 | uniRoot method epsrel maxit fun xl xu = uniRootGen (fi (fromEnum method)) fun xl xu epsrel maxit | ||
76 | |||
77 | uniRootGen m f xl xu epsrel maxit = unsafePerformIO $ do | ||
78 | fp <- mkDoublefun f | ||
79 | rawpath <- createMIO maxit 4 | ||
80 | (c_root m fp epsrel (fi maxit) xl xu) | ||
81 | "root" | ||
82 | let it = round (rawpath @@> (maxit-1,0)) | ||
83 | path = takeRows it rawpath | ||
84 | [sol] = toLists $ dropRows (it-1) path | ||
85 | freeHaskellFunPtr fp | ||
86 | return (sol !! 1, path) | ||
87 | |||
88 | foreign import ccall safe "root" | ||
89 | c_root:: CInt -> FunPtr (Double -> Double) -> Double -> CInt -> Double -> Double -> TM | ||
90 | |||
91 | ------------------------------------------------------------------------- | ||
92 | data UniRootMethodJ = UNewton | ||
93 | | Secant | ||
94 | | Steffenson | ||
95 | deriving (Enum, Eq, Show, Bounded) | ||
96 | |||
97 | uniRootJ :: UniRootMethodJ | ||
98 | -> Double | ||
99 | -> Int | ||
100 | -> (Double -> Double) | ||
101 | -> (Double -> Double) | ||
102 | -> Double | ||
103 | -> (Double, Matrix Double) | ||
104 | uniRootJ method epsrel maxit fun dfun x = uniRootJGen (fi (fromEnum method)) fun | ||
105 | dfun x epsrel maxit | ||
106 | |||
107 | uniRootJGen m f df x epsrel maxit = unsafePerformIO $ do | ||
108 | fp <- mkDoublefun f | ||
109 | dfp <- mkDoublefun df | ||
110 | rawpath <- createMIO maxit 2 | ||
111 | (c_rootj m fp dfp epsrel (fi maxit) x) | ||
112 | "rootj" | ||
113 | let it = round (rawpath @@> (maxit-1,0)) | ||
114 | path = takeRows it rawpath | ||
115 | [sol] = toLists $ dropRows (it-1) path | ||
116 | freeHaskellFunPtr fp | ||
117 | return (sol !! 1, path) | ||
118 | |||
119 | foreign import ccall safe "rootj" | ||
120 | c_rootj :: CInt -> FunPtr (Double -> Double) -> FunPtr (Double -> Double) | ||
121 | -> Double -> CInt -> Double -> TM | ||
122 | |||
123 | ------------------------------------------------------------------------- | ||
124 | |||
61 | data RootMethod = Hybrids | 125 | data RootMethod = Hybrids |
62 | | Hybrid | 126 | | Hybrid |
63 | | DNewton | 127 | | DNewton |
@@ -82,8 +146,8 @@ rootGen m f xi epsabs maxit = unsafePerformIO $ do | |||
82 | fp <- mkVecVecfun (aux_vTov (checkdim1 n . fromList . f . toList)) | 146 | fp <- mkVecVecfun (aux_vTov (checkdim1 n . fromList . f . toList)) |
83 | rawpath <- vec xiv $ \xiv' -> | 147 | rawpath <- vec xiv $ \xiv' -> |
84 | createMIO maxit (2*n+1) | 148 | createMIO maxit (2*n+1) |
85 | (c_root m fp epsabs (fi maxit) // xiv') | 149 | (c_multiroot m fp epsabs (fi maxit) // xiv') |
86 | "root" | 150 | "multiroot" |
87 | let it = round (rawpath @@> (maxit-1,0)) | 151 | let it = round (rawpath @@> (maxit-1,0)) |
88 | path = takeRows it rawpath | 152 | path = takeRows it rawpath |
89 | [sol] = toLists $ dropRows (it-1) path | 153 | [sol] = toLists $ dropRows (it-1) path |
@@ -91,8 +155,8 @@ rootGen m f xi epsabs maxit = unsafePerformIO $ do | |||
91 | return (take n $ drop 1 sol, path) | 155 | return (take n $ drop 1 sol, path) |
92 | 156 | ||
93 | 157 | ||
94 | foreign import ccall safe "root" | 158 | foreign import ccall safe "multiroot" |
95 | c_root:: CInt -> FunPtr TVV -> Double -> CInt -> TVM | 159 | c_multiroot:: CInt -> FunPtr TVV -> Double -> CInt -> TVM |
96 | 160 | ||
97 | ------------------------------------------------------------------------- | 161 | ------------------------------------------------------------------------- |
98 | 162 | ||
@@ -120,8 +184,8 @@ rootJGen m f jac xi epsabs maxit = unsafePerformIO $ do | |||
120 | jp <- mkVecMatfun (aux_vTom (checkdim2 n . fromLists . jac . toList)) | 184 | jp <- mkVecMatfun (aux_vTom (checkdim2 n . fromLists . jac . toList)) |
121 | rawpath <- vec xiv $ \xiv' -> | 185 | rawpath <- vec xiv $ \xiv' -> |
122 | createMIO maxit (2*n+1) | 186 | createMIO maxit (2*n+1) |
123 | (c_rootj m fp jp epsabs (fi maxit) // xiv') | 187 | (c_multirootj m fp jp epsabs (fi maxit) // xiv') |
124 | "root" | 188 | "multiroot" |
125 | let it = round (rawpath @@> (maxit-1,0)) | 189 | let it = round (rawpath @@> (maxit-1,0)) |
126 | path = takeRows it rawpath | 190 | path = takeRows it rawpath |
127 | [sol] = toLists $ dropRows (it-1) path | 191 | [sol] = toLists $ dropRows (it-1) path |
@@ -129,9 +193,8 @@ rootJGen m f jac xi epsabs maxit = unsafePerformIO $ do | |||
129 | freeHaskellFunPtr jp | 193 | freeHaskellFunPtr jp |
130 | return (take n $ drop 1 sol, path) | 194 | return (take n $ drop 1 sol, path) |
131 | 195 | ||
132 | 196 | foreign import ccall safe "multirootj" | |
133 | foreign import ccall safe "rootj" | 197 | c_multirootj:: CInt -> FunPtr TVV -> FunPtr TVM -> Double -> CInt -> TVM |
134 | c_rootj:: CInt -> FunPtr TVV -> FunPtr TVM -> Double -> CInt -> TVM | ||
135 | 198 | ||
136 | ------------------------------------------------------- | 199 | ------------------------------------------------------- |
137 | 200 | ||
diff --git a/lib/Numeric/GSL/gsl-aux.c b/lib/Numeric/GSL/gsl-aux.c index fc14ff5..e727c91 100644 --- a/lib/Numeric/GSL/gsl-aux.c +++ b/lib/Numeric/GSL/gsl-aux.c | |||
@@ -32,6 +32,7 @@ | |||
32 | #include <gsl/gsl_complex_math.h> | 32 | #include <gsl/gsl_complex_math.h> |
33 | #include <gsl/gsl_rng.h> | 33 | #include <gsl/gsl_rng.h> |
34 | #include <gsl/gsl_randist.h> | 34 | #include <gsl/gsl_randist.h> |
35 | #include <gsl/gsl_roots.h> | ||
35 | #include <gsl/gsl_multifit_nlin.h> | 36 | #include <gsl/gsl_multifit_nlin.h> |
36 | #include <string.h> | 37 | #include <string.h> |
37 | #include <stdio.h> | 38 | #include <stdio.h> |
@@ -937,7 +938,7 @@ int minimize(int method, double f(int, double*), double tolsize, int maxit, | |||
937 | iter++; | 938 | iter++; |
938 | if (status) break; | 939 | if (status) break; |
939 | status = gsl_multimin_test_size (size, tolsize); | 940 | status = gsl_multimin_test_size (size, tolsize); |
940 | } while (status == GSL_CONTINUE && iter <= maxit); | 941 | } while (status == GSL_CONTINUE && iter < maxit); |
941 | int i,j; | 942 | int i,j; |
942 | for (i=iter; i<solr; i++) { | 943 | for (i=iter; i<solr; i++) { |
943 | solp[i*solc+0] = iter; | 944 | solp[i*solc+0] = iter; |
@@ -1033,7 +1034,7 @@ int minimizeD(int method, double f(int, double*), int df(int, double*, int, doub | |||
1033 | iter++; | 1034 | iter++; |
1034 | if (status) break; | 1035 | if (status) break; |
1035 | status = gsl_multimin_test_gradient (s->gradient, tolgrad); | 1036 | status = gsl_multimin_test_gradient (s->gradient, tolgrad); |
1036 | } while (status == GSL_CONTINUE && iter <= maxit); | 1037 | } while (status == GSL_CONTINUE && iter < maxit); |
1037 | int i,j; | 1038 | int i,j; |
1038 | for (i=iter; i<solr; i++) { | 1039 | for (i=iter; i<solr; i++) { |
1039 | solp[i*solc+0] = iter; | 1040 | solp[i*solc+0] = iter; |
@@ -1047,9 +1048,144 @@ int minimizeD(int method, double f(int, double*), int df(int, double*, int, doub | |||
1047 | 1048 | ||
1048 | //--------------------------------------------------------------- | 1049 | //--------------------------------------------------------------- |
1049 | 1050 | ||
1051 | double only_f_aux_root(double x, void *pars) { | ||
1052 | double (*f)(double) = (double (*)(double)) pars; | ||
1053 | return f(x); | ||
1054 | } | ||
1055 | |||
1056 | int root(int method, double f(double), | ||
1057 | double epsrel, int maxit, | ||
1058 | double xl, double xu, RMAT(sol)) { | ||
1059 | REQUIRES(solr == maxit && solc == 4,BAD_SIZE); | ||
1060 | DEBUGMSG("root_only_f"); | ||
1061 | gsl_function my_func; | ||
1062 | // extract function from pars | ||
1063 | my_func.function = only_f_aux_root; | ||
1064 | my_func.params = f; | ||
1065 | size_t iter = 0; | ||
1066 | int status; | ||
1067 | const gsl_root_fsolver_type *T; | ||
1068 | gsl_root_fsolver *s; | ||
1069 | // Starting point | ||
1070 | switch(method) { | ||
1071 | case 0 : {T = gsl_root_fsolver_bisection; printf("7\n"); break; } | ||
1072 | case 1 : {T = gsl_root_fsolver_falsepos; break; } | ||
1073 | case 2 : {T = gsl_root_fsolver_brent; break; } | ||
1074 | default: ERROR(BAD_CODE); | ||
1075 | } | ||
1076 | s = gsl_root_fsolver_alloc (T); | ||
1077 | gsl_root_fsolver_set (s, &my_func, xl, xu); | ||
1078 | do { | ||
1079 | double best, current_lo, current_hi; | ||
1080 | status = gsl_root_fsolver_iterate (s); | ||
1081 | best = gsl_root_fsolver_root (s); | ||
1082 | current_lo = gsl_root_fsolver_x_lower (s); | ||
1083 | current_hi = gsl_root_fsolver_x_upper (s); | ||
1084 | solp[iter*solc] = iter + 1; | ||
1085 | solp[iter*solc+1] = best; | ||
1086 | solp[iter*solc+2] = current_lo; | ||
1087 | solp[iter*solc+3] = current_hi; | ||
1088 | iter++; | ||
1089 | if (status) /* check if solver is stuck */ | ||
1090 | break; | ||
1091 | |||
1092 | status = | ||
1093 | gsl_root_test_interval (current_lo, current_hi, 0, epsrel); | ||
1094 | } | ||
1095 | while (status == GSL_CONTINUE && iter < maxit); | ||
1096 | int i; | ||
1097 | for (i=iter; i<solr; i++) { | ||
1098 | solp[i*solc+0] = iter; | ||
1099 | solp[i*solc+1]=0.; | ||
1100 | solp[i*solc+2]=0.; | ||
1101 | solp[i*solc+3]=0.; | ||
1102 | } | ||
1103 | gsl_root_fsolver_free(s); | ||
1104 | OK | ||
1105 | } | ||
1106 | |||
1107 | typedef struct { | ||
1108 | double (*f)(double); | ||
1109 | double (*jf)(double); | ||
1110 | } uniTfjf; | ||
1111 | |||
1112 | double f_aux_uni(double x, void *pars) { | ||
1113 | uniTfjf * fjf = ((uniTfjf*) pars); | ||
1114 | return (fjf->f)(x); | ||
1115 | } | ||
1116 | |||
1117 | double jf_aux_uni(double x, void * pars) { | ||
1118 | uniTfjf * fjf = ((uniTfjf*) pars); | ||
1119 | return (fjf->jf)(x); | ||
1120 | } | ||
1121 | |||
1122 | void fjf_aux_uni(double x, void * pars, double * f, double * g) { | ||
1123 | *f = f_aux_uni(x,pars); | ||
1124 | *g = jf_aux_uni(x,pars); | ||
1125 | } | ||
1126 | |||
1127 | int rootj(int method, double f(double), | ||
1128 | double df(double), | ||
1129 | double epsrel, int maxit, | ||
1130 | double x, RMAT(sol)) { | ||
1131 | REQUIRES(solr == maxit && solc == 2,BAD_SIZE); | ||
1132 | DEBUGMSG("root_fjf"); | ||
1133 | gsl_function_fdf my_func; | ||
1134 | // extract function from pars | ||
1135 | my_func.f = f_aux_uni; | ||
1136 | my_func.df = jf_aux_uni; | ||
1137 | my_func.fdf = fjf_aux_uni; | ||
1138 | uniTfjf stfjf; | ||
1139 | stfjf.f = f; | ||
1140 | stfjf.jf = df; | ||
1141 | my_func.params = &stfjf; | ||
1142 | size_t iter = 0; | ||
1143 | int status; | ||
1144 | const gsl_root_fdfsolver_type *T; | ||
1145 | gsl_root_fdfsolver *s; | ||
1146 | // Starting point | ||
1147 | switch(method) { | ||
1148 | case 0 : {T = gsl_root_fdfsolver_newton;; break; } | ||
1149 | case 1 : {T = gsl_root_fdfsolver_secant; break; } | ||
1150 | case 2 : {T = gsl_root_fdfsolver_steffenson; break; } | ||
1151 | default: ERROR(BAD_CODE); | ||
1152 | } | ||
1153 | s = gsl_root_fdfsolver_alloc (T); | ||
1154 | |||
1155 | gsl_root_fdfsolver_set (s, &my_func, x); | ||
1156 | |||
1157 | do { | ||
1158 | double x0; | ||
1159 | status = gsl_root_fdfsolver_iterate (s); | ||
1160 | x0 = x; | ||
1161 | x = gsl_root_fdfsolver_root(s); | ||
1162 | solp[iter*solc+0] = iter+1; | ||
1163 | solp[iter*solc+1] = x; | ||
1164 | |||
1165 | iter++; | ||
1166 | if (status) /* check if solver is stuck */ | ||
1167 | break; | ||
1168 | |||
1169 | status = | ||
1170 | gsl_root_test_delta (x, x0, 0, epsrel); | ||
1171 | } | ||
1172 | while (status == GSL_CONTINUE && iter < maxit); | ||
1173 | |||
1174 | int i; | ||
1175 | for (i=iter; i<solr; i++) { | ||
1176 | solp[i*solc+0] = iter; | ||
1177 | solp[i*solc+1]=0.; | ||
1178 | } | ||
1179 | gsl_root_fdfsolver_free(s); | ||
1180 | OK | ||
1181 | } | ||
1182 | |||
1183 | |||
1184 | //--------------------------------------------------------------- | ||
1185 | |||
1050 | typedef void TrawfunV(int, double*, int, double*); | 1186 | typedef void TrawfunV(int, double*, int, double*); |
1051 | 1187 | ||
1052 | int only_f_aux_root(const gsl_vector*x, void *pars, gsl_vector*y) { | 1188 | int only_f_aux_multiroot(const gsl_vector*x, void *pars, gsl_vector*y) { |
1053 | TrawfunV * f = (TrawfunV*) pars; | 1189 | TrawfunV * f = (TrawfunV*) pars; |
1054 | double* p = (double*)calloc(x->size,sizeof(double)); | 1190 | double* p = (double*)calloc(x->size,sizeof(double)); |
1055 | double* q = (double*)calloc(y->size,sizeof(double)); | 1191 | double* q = (double*)calloc(y->size,sizeof(double)); |
@@ -1066,14 +1202,14 @@ int only_f_aux_root(const gsl_vector*x, void *pars, gsl_vector*y) { | |||
1066 | return 0; //hmmm | 1202 | return 0; //hmmm |
1067 | } | 1203 | } |
1068 | 1204 | ||
1069 | int root(int method, void f(int, double*, int, double*), | 1205 | int multiroot(int method, void f(int, double*, int, double*), |
1070 | double epsabs, int maxit, | 1206 | double epsabs, int maxit, |
1071 | KRVEC(xi), RMAT(sol)) { | 1207 | KRVEC(xi), RMAT(sol)) { |
1072 | REQUIRES(solr == maxit && solc == 1+2*xin,BAD_SIZE); | 1208 | REQUIRES(solr == maxit && solc == 1+2*xin,BAD_SIZE); |
1073 | DEBUGMSG("root_only_f"); | 1209 | DEBUGMSG("root_only_f"); |
1074 | gsl_multiroot_function my_func; | 1210 | gsl_multiroot_function my_func; |
1075 | // extract function from pars | 1211 | // extract function from pars |
1076 | my_func.f = only_f_aux_root; | 1212 | my_func.f = only_f_aux_multiroot; |
1077 | my_func.n = xin; | 1213 | my_func.n = xin; |
1078 | my_func.params = f; | 1214 | my_func.params = f; |
1079 | size_t iter = 0; | 1215 | size_t iter = 0; |
@@ -1112,7 +1248,7 @@ int root(int method, void f(int, double*, int, double*), | |||
1112 | status = | 1248 | status = |
1113 | gsl_multiroot_test_residual (s->f, epsabs); | 1249 | gsl_multiroot_test_residual (s->f, epsabs); |
1114 | } | 1250 | } |
1115 | while (status == GSL_CONTINUE && iter <= maxit); | 1251 | while (status == GSL_CONTINUE && iter < maxit); |
1116 | 1252 | ||
1117 | int i,j; | 1253 | int i,j; |
1118 | for (i=iter; i<solr; i++) { | 1254 | for (i=iter; i<solr; i++) { |
@@ -1175,7 +1311,7 @@ int fjf_aux(const gsl_vector * x, void * pars, gsl_vector * f, gsl_matrix * g) { | |||
1175 | return 0; | 1311 | return 0; |
1176 | } | 1312 | } |
1177 | 1313 | ||
1178 | int rootj(int method, int f(int, double*, int, double*), | 1314 | int multirootj(int method, int f(int, double*, int, double*), |
1179 | int jac(int, double*, int, int, double*), | 1315 | int jac(int, double*, int, int, double*), |
1180 | double epsabs, int maxit, | 1316 | double epsabs, int maxit, |
1181 | KRVEC(xi), RMAT(sol)) { | 1317 | KRVEC(xi), RMAT(sol)) { |
@@ -1228,7 +1364,7 @@ int rootj(int method, int f(int, double*, int, double*), | |||
1228 | status = | 1364 | status = |
1229 | gsl_multiroot_test_residual (s->f, epsabs); | 1365 | gsl_multiroot_test_residual (s->f, epsabs); |
1230 | } | 1366 | } |
1231 | while (status == GSL_CONTINUE && iter <= maxit); | 1367 | while (status == GSL_CONTINUE && iter < maxit); |
1232 | 1368 | ||
1233 | int i,j; | 1369 | int i,j; |
1234 | for (i=iter; i<solr; i++) { | 1370 | for (i=iter; i<solr; i++) { |
@@ -1293,7 +1429,7 @@ int nlfit(int method, int f(int, double*, int, double*), | |||
1293 | 1429 | ||
1294 | status = gsl_multifit_test_delta (s->dx, s->x, epsabs, epsrel); | 1430 | status = gsl_multifit_test_delta (s->dx, s->x, epsabs, epsrel); |
1295 | } | 1431 | } |
1296 | while (status == GSL_CONTINUE && iter <= maxit); | 1432 | while (status == GSL_CONTINUE && iter < maxit); |
1297 | 1433 | ||
1298 | int i,j; | 1434 | int i,j; |
1299 | for (i=iter; i<solr; i++) { | 1435 | for (i=iter; i<solr; i++) { |