summaryrefslogtreecommitdiff
path: root/lib/Numeric/GSL/gsl-aux.c
diff options
context:
space:
mode:
Diffstat (limited to 'lib/Numeric/GSL/gsl-aux.c')
-rw-r--r--lib/Numeric/GSL/gsl-aux.c65
1 files changed, 61 insertions, 4 deletions
diff --git a/lib/Numeric/GSL/gsl-aux.c b/lib/Numeric/GSL/gsl-aux.c
index 61b6a54..27ffbdc 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>
@@ -1047,9 +1048,65 @@ int minimizeD(int method, double f(int, double*), int df(int, double*, int, doub
1047 1048
1048//--------------------------------------------------------------- 1049//---------------------------------------------------------------
1049 1050
1051double only_f_aux_root(double x, void *pars) {
1052 double (*f)(double) = (double (*)(double)) pars;
1053 return f(x);
1054}
1055
1056int 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
1050typedef void TrawfunV(int, double*, int, double*); 1107typedef void TrawfunV(int, double*, int, double*);
1051 1108
1052int only_f_aux_root(const gsl_vector*x, void *pars, gsl_vector*y) { 1109int only_f_aux_multiroot(const gsl_vector*x, void *pars, gsl_vector*y) {
1053 TrawfunV * f = (TrawfunV*) pars; 1110 TrawfunV * f = (TrawfunV*) pars;
1054 double* p = (double*)calloc(x->size,sizeof(double)); 1111 double* p = (double*)calloc(x->size,sizeof(double));
1055 double* q = (double*)calloc(y->size,sizeof(double)); 1112 double* q = (double*)calloc(y->size,sizeof(double));
@@ -1066,14 +1123,14 @@ int only_f_aux_root(const gsl_vector*x, void *pars, gsl_vector*y) {
1066 return 0; //hmmm 1123 return 0; //hmmm
1067} 1124}
1068 1125
1069int root(int method, void f(int, double*, int, double*), 1126int multiroot(int method, void f(int, double*, int, double*),
1070 double epsabs, int maxit, 1127 double epsabs, int maxit,
1071 KRVEC(xi), RMAT(sol)) { 1128 KRVEC(xi), RMAT(sol)) {
1072 REQUIRES(solr == maxit && solc == 1+2*xin,BAD_SIZE); 1129 REQUIRES(solr == maxit && solc == 1+2*xin,BAD_SIZE);
1073 DEBUGMSG("root_only_f"); 1130 DEBUGMSG("root_only_f");
1074 gsl_multiroot_function my_func; 1131 gsl_multiroot_function my_func;
1075 // extract function from pars 1132 // extract function from pars
1076 my_func.f = only_f_aux_root; 1133 my_func.f = only_f_aux_multiroot;
1077 my_func.n = xin; 1134 my_func.n = xin;
1078 my_func.params = f; 1135 my_func.params = f;
1079 size_t iter = 0; 1136 size_t iter = 0;
@@ -1175,7 +1232,7 @@ int fjf_aux(const gsl_vector * x, void * pars, gsl_vector * f, gsl_matrix * g) {
1175 return 0; 1232 return 0;
1176} 1233}
1177 1234
1178int rootj(int method, int f(int, double*, int, double*), 1235int multirootj(int method, int f(int, double*, int, double*),
1179 int jac(int, double*, int, int, double*), 1236 int jac(int, double*, int, int, double*),
1180 double epsabs, int maxit, 1237 double epsabs, int maxit,
1181 KRVEC(xi), RMAT(sol)) { 1238 KRVEC(xi), RMAT(sol)) {