summaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorAlex Lang <me@alang.ca>2013-05-08 14:33:24 +0900
committerAlex Lang <me@alang.ca>2013-05-08 14:33:24 +0900
commit10ad01cecff6fef58fd9b8678b122172cc84a7f2 (patch)
treeceb2c926fbc82d53f24ce0a89d42f56f520f84bf
parentedb20f3c9993d9c9c349f6a001317592949154f1 (diff)
Implement uniRootJ for one-dimensional root-finding with derivatives
-rw-r--r--lib/Numeric/GSL/Root.hs37
-rw-r--r--lib/Numeric/GSL/gsl-aux.c79
2 files changed, 114 insertions, 2 deletions
diff --git a/lib/Numeric/GSL/Root.hs b/lib/Numeric/GSL/Root.hs
index 25ec5f5..6da15e5 100644
--- a/lib/Numeric/GSL/Root.hs
+++ b/lib/Numeric/GSL/Root.hs
@@ -46,6 +46,7 @@ main = do
46 46
47module Numeric.GSL.Root ( 47module Numeric.GSL.Root (
48 uniRoot, UniRootMethod(..), 48 uniRoot, UniRootMethod(..),
49 uniRootJ, UniRootMethodJ(..),
49 root, RootMethod(..), 50 root, RootMethod(..),
50 rootJ, RootMethodJ(..), 51 rootJ, RootMethodJ(..),
51) where 52) where
@@ -88,6 +89,38 @@ foreign import ccall safe "root"
88 c_root:: CInt -> FunPtr (Double -> Double) -> Double -> CInt -> Double -> Double -> TM 89 c_root:: CInt -> FunPtr (Double -> Double) -> Double -> CInt -> Double -> Double -> TM
89 90
90------------------------------------------------------------------------- 91-------------------------------------------------------------------------
92data UniRootMethodJ = UNewton
93 | Secant
94 | Steffenson
95 deriving (Enum, Eq, Show, Bounded)
96
97uniRootJ :: UniRootMethodJ
98 -> Double
99 -> Int
100 -> (Double -> Double)
101 -> (Double -> Double)
102 -> Double
103 -> (Double, Matrix Double)
104uniRootJ method epsrel maxit fun dfun x = uniRootJGen (fi (fromEnum method)) fun
105 dfun x epsrel maxit
106
107uniRootJGen 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
119foreign import ccall safe "rootj"
120 c_rootj :: CInt -> FunPtr (Double -> Double) -> FunPtr (Double -> Double)
121 -> Double -> CInt -> Double -> TM
122
123-------------------------------------------------------------------------
91 124
92data RootMethod = Hybrids 125data RootMethod = Hybrids
93 | Hybrid 126 | Hybrid
@@ -151,7 +184,7 @@ rootJGen m f jac xi epsabs maxit = unsafePerformIO $ do
151 jp <- mkVecMatfun (aux_vTom (checkdim2 n . fromLists . jac . toList)) 184 jp <- mkVecMatfun (aux_vTom (checkdim2 n . fromLists . jac . toList))
152 rawpath <- vec xiv $ \xiv' -> 185 rawpath <- vec xiv $ \xiv' ->
153 createMIO maxit (2*n+1) 186 createMIO maxit (2*n+1)
154 (c_rootj m fp jp epsabs (fi maxit) // xiv') 187 (c_multirootj m fp jp epsabs (fi maxit) // xiv')
155 "multiroot" 188 "multiroot"
156 let it = round (rawpath @@> (maxit-1,0)) 189 let it = round (rawpath @@> (maxit-1,0))
157 path = takeRows it rawpath 190 path = takeRows it rawpath
@@ -161,7 +194,7 @@ rootJGen m f jac xi epsabs maxit = unsafePerformIO $ do
161 return (take n $ drop 1 sol, path) 194 return (take n $ drop 1 sol, path)
162 195
163foreign import ccall safe "multirootj" 196foreign import ccall safe "multirootj"
164 c_rootj:: CInt -> FunPtr TVV -> FunPtr TVM -> Double -> CInt -> TVM 197 c_multirootj:: CInt -> FunPtr TVV -> FunPtr TVM -> Double -> CInt -> TVM
165 198
166------------------------------------------------------- 199-------------------------------------------------------
167 200
diff --git a/lib/Numeric/GSL/gsl-aux.c b/lib/Numeric/GSL/gsl-aux.c
index 27ffbdc..e727c91 100644
--- a/lib/Numeric/GSL/gsl-aux.c
+++ b/lib/Numeric/GSL/gsl-aux.c
@@ -1104,6 +1104,85 @@ int root(int method, double f(double),
1104 OK 1104 OK
1105} 1105}
1106 1106
1107typedef struct {
1108 double (*f)(double);
1109 double (*jf)(double);
1110} uniTfjf;
1111
1112double f_aux_uni(double x, void *pars) {
1113 uniTfjf * fjf = ((uniTfjf*) pars);
1114 return (fjf->f)(x);
1115}
1116
1117double jf_aux_uni(double x, void * pars) {
1118 uniTfjf * fjf = ((uniTfjf*) pars);
1119 return (fjf->jf)(x);
1120}
1121
1122void 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
1127int 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
1107typedef void TrawfunV(int, double*, int, double*); 1186typedef void TrawfunV(int, double*, int, double*);
1108 1187
1109int only_f_aux_multiroot(const gsl_vector*x, void *pars, gsl_vector*y) { 1188int only_f_aux_multiroot(const gsl_vector*x, void *pars, gsl_vector*y) {