chiark / gitweb /
helixish: attempt at the whole thing
[moebius3.git] / findcurve.c
index 169116a826df4969cd9cfdab8c1ae9ad0c037829..c038182e685ca954002d0e0473338077c01f5af4 100644 (file)
@@ -5,11 +5,12 @@
 
 #include <gsl/gsl_errno.h>
 #include <gsl/gsl_sf.h>
-#include <gsl/gsl_multiroots.h>
+#include <gsl/gsl_siman.h>
+#include <gsl/gsl_randist.h>
 
 #include "symbolic.c"
 
-#define X(i) gsl_vector_get(x,i)
+#define X(i) x[i]
 #define J_END_COL(i) \
   for (j=0; j<N; j++) gsl_matrix_set(J,i,j,J_COL[j]);
 
@@ -18,33 +19,65 @@ static inline _Bool IS_SMALL(double v) {
 }
 
 static inline double sinc(double x) {
-  return gsl_sf_sinc(x / M_PI);  
+  return gsl_sf_sinc(x / M_PI);
 }
 
-static int cb_f(const gsl_vector *x, void *params, gsl_vector *f) {
-  double *target = params;
-  double F[6];
+static double target[N];
 
+static double cb_Efunc(void *xp) {
+  const double *x = xp;
+  double F[N], e;
   int i;
   X_EXTRACT;
   F_POPULATE;
-  for (i=0; i<N; i++) gsl_vector_set(f,i, F[i] - target[i]);
-  return 0;
+  for (i=0, e=0; i<N; i++) {
+    double d = F[i] - target[i];
+    e += d*d;
+  }
+  //printf("\n cb_Efunc %p %10.7f [", xp, e);
+  //for (i=0; i<N; i++) printf(" %10.7f,", x[i]);
+  //printf("]\n");
+  return e;
+}
+
+static void cb_step(const gsl_rng *rng, void *xp, double step_size) {
+  double *x = xp;
+  int i;
+  double step[N];
+  gsl_ran_dir_nd(rng, N,step);
+  for (i=0; i<N; i++)
+    x[i] += step_size * step[i];
+  //printf("\n cb_step %p %10.7f [", xp, step_size);
+  //for (i=0; i<N; i++) printf(" %10.7f,", step[i]);
+  //printf("]\n");
+}
+
+static double cb_metric(void *xp, void *yp) {
+  const double *x=xp, *y=yp;
+  int i;
+  double s;
+  for (i=0, s=0; i<N; i++) {
+    double d = x[i] - y[i];
+    s += d*d;
+  }
+  return sqrt(s);
 }
 
-static int cb_df(const gsl_vector *x, void *params, gsl_matrix *J) {
-  double J_COL[N];
-  int j;
+static void printcore(const double *x) {
+  int i;
+  double F[N];
   X_EXTRACT;
-  J_POPULATE;
-  return 0;
+  F_POPULATE;
+  printf("[");
+  for (i=0; i<6; i++) printf(" %.18g,", x[i]);
+  for (i=0; i<6; i++) printf(" %.18g,", F[i]);
+  printf(" ]\n");
 }
 
-static int cb_fdf(const gsl_vector *x, void *params,
-                 gsl_vector *f, gsl_matrix *J) {
-  cb_f(x,params,f);
-  cb_df(x,params,J);
-  return 0;
+static void __attribute__((unused)) cb_print(void *xp) {
+  const double *x = xp;
+  printf("\n");
+  printcore(x);
 }
 
 static double scan1double(void) {
@@ -59,43 +92,42 @@ static double scan1double(void) {
 }
 
 int main(void) {
-  double target[N], epsilon;
+  double epsilon;
   int i;
 
-  gsl_multiroot_function_fdf function;
-  function.f = cb_f;
-  function.df = cb_df;
-  function.fdf = cb_fdf;
-  function.n = N;
-  function.params = target;
-
-  gsl_vector *startpoint = gsl_vector_alloc(N);
+  double startpoint[N];
 
-  gsl_multiroot_fdfsolver *solver
-    = gsl_multiroot_fdfsolver_alloc(gsl_multiroot_fdfsolver_hybridsj, N);
+  gsl_rng *rng = gsl_rng_alloc(gsl_rng_ranlxd2);
 
   for (;;) {
     /* 2N+1 doubles: target, initial guess, epsilon for residual */
     for (i=0; i<N; i++)
       target[i] = scan1double();
     for (i=0; i<N; i++)
-      gsl_vector_set(startpoint,i,scan1double());
+      startpoint[i] = scan1double();
     epsilon = scan1double();
 
-    gsl_multiroot_fdfsolver_set(solver, &function, startpoint);
+    gsl_rng_set(rng,0);
+
+    gsl_siman_params_t siman_params = {
+      .k = 1.0,
+      .t_initial = 0.5,
+      .mu_t = 1.001,
+      .t_min = epsilon * 1E-3,
+      .iters_fixed_T = 100,
+      .n_tries = 10,
+      .step_size = 0.05,
+    };
 
-    int r;
-    for (;;) {
-      r = gsl_multiroot_fdfsolver_iterate(solver);
-      if (r) break;
+    gsl_siman_solve(rng,
+                   startpoint,
+                   cb_Efunc, cb_step, cb_metric,
+                   0, //cb_print,
+                   0,0,0, sizeof(startpoint), siman_params);
 
-      r = gsl_multiroot_test_residual(solver->f, epsilon);
-      if (r != GSL_CONTINUE) break;
-    }
+    printcore(startpoint);
 
-    if (r) {
-      fprintf(stderr,"ERROR %s\n",gsl_strerror(r));
-      exit(-1);
-    }
+    printf("[]\n");
+    fflush(stdout);
   }
 }