chiark / gitweb /
curveopt: wip, demos more of the same bug
[moebius3.git] / findcurve.c
1 /*
2  */
3
4 #include <math.h>
5 #include <string.h>
6
7 #include <gsl/gsl_errno.h>
8 #include <gsl/gsl_sf.h>
9 #include <gsl/gsl_siman.h>
10 #include <gsl/gsl_randist.h>
11
12 #include "symbolic.c"
13
14 #define J_END_COL(i) \
15   for (j=0; j<PN; j++) gsl_matrix_set(J,i,j,J_COL[j]);
16
17 static inline _Bool IS_SMALL(double v) {
18   return v < 1E6;
19 }
20
21 static inline double sinc(double x) {
22   return gsl_sf_sinc(x / M_PI);
23 }
24
25 static int NP;
26 static double *INPUT; /* dyanmic array, on main's stack */
27 static double PREP[NPREP];
28
29 static void printcore(const double *X) {
30   int i, j;
31   DECLARE_F_G;
32   CALCULATE_F_G;
33   printf("[");
34   for (i=0; i<NP; i++)
35     for (j=0; j<3; j++)
36       printf(" %25.18g,", POINT(i)[j]);
37   printf(" ]\n");
38 }
39
40 static void prepare(double X[] /* startpoint */) {
41   /* fills in PREP and startpoint */
42   PREPARE;
43 }
44
45 //#define DEBUG
46
47 static double cb_Efunc(void *xp) {
48   const double *X = xp;
49   DECLARE_F_G;
50   CALCULATE_F_G;
51
52 #ifdef DEBUG
53   int i,j;
54   printf(" Efunc\n");
55   for (j=0; j<3; j++) {
56     for (i=0; i<NP; i++)
57       printf(" %7.4f", POINT(i)[j]);
58     printf("\n");
59   }
60   printf("             ");
61 #endif
62
63   double e = 0;
64   int P;
65   for (P=0; P<NP-3; P++) {
66     double P_cost;
67     CALCULATE_COST;
68 #ifdef DEBUG
69     printf(" %7.4f", P_cost);
70 #endif
71     e += P_cost;
72   }
73 #ifdef DEBUG
74   printf("\n");
75 #endif
76   return e;
77 }
78
79 static void cb_step(const gsl_rng *rng, void *xp, double step_size) {
80   double *x = xp;
81   int i;
82   double step[NX];
83   gsl_ran_dir_nd(rng, NX, step);
84   for (i=0; i<NX; i++)
85     x[i] += step_size * step[i];
86   //printf("\n cb_step %p %10.7f [", xp, step_size);
87   //for (i=0; i<N; i++) printf(" %10.7f,", step[i]);
88   //printf("]\n");
89 }
90
91 static double cb_metric(void *xp, void *yp) {
92   const double *x=xp, *y=yp;
93   int i;
94   double s;
95   for (i=0, s=0; i<NX; i++) {
96     double d = x[i] - y[i];
97     s += d*d;
98   }
99   return sqrt(s);
100 }
101
102 static void __attribute__((unused)) cb_print(void *xp) {
103   const double *x = xp;
104   printf("\n");
105   printcore(x);
106 }
107
108 static double scan1double(void) {
109   double v;
110   int r;
111
112   r = scanf("%lf",&v);
113   if (ferror(stdin)) { perror("stdin"); exit(-1); }
114   if (feof(stdin)) exit(0);
115   if (r!=1) { fputs("bad input\n",stderr); exit(-1); }
116   return v;
117 }
118
119 int main(int argc, const char *const *argv) {
120   double epsilon;
121   int i;
122
123   NP = atoi(argv[1]);
124   epsilon = atof(argv[2]);
125
126   gsl_rng *rng = gsl_rng_alloc(gsl_rng_ranlxd2);
127
128   double input[NINPUT]; INPUT = input;
129   double startpoint[NX];
130
131   for (;;) {
132     /* NINPUT + 1 doubles: startpoint, epsilon for residual */
133     for (i=0; i<NINPUT; i++)
134       INPUT[i] = scan1double();
135
136     gsl_rng_set(rng,0);
137
138     gsl_siman_params_t siman_params = {
139       .k = 1.0,
140       .t_initial = 0.5,
141       .mu_t = 1.001,
142       .t_min = epsilon * 1E-3,
143       .iters_fixed_T = 100,
144       .n_tries = 10,
145       .step_size = 0.05,
146     };
147
148     prepare(startpoint);
149
150     gsl_siman_solve(rng,
151                     startpoint,
152                     cb_Efunc, cb_step, cb_metric,
153                     0, // cb_print,
154                     0,0,0, sizeof(startpoint), siman_params);
155
156     printcore(startpoint);
157
158     printf("[]\n");
159     fflush(stdout);
160   }
161 }