From 78d8de381f0b1e99ad98364709cbf876689628b2 Mon Sep 17 00:00:00 2001
From: Jaron Kent-Dobias <jaron@kent-dobias.com>
Date: Fri, 29 Jun 2018 15:43:10 -0400
Subject: completely removed measurement during the simulation, opting to just
 save binary data points to files throughout

---
 lib/cluster_finite.c |   4 +-
 lib/cluster_finite.h |   2 +-
 lib/initial_finite.c |  22 ++++++--
 lib/initial_finite.h |   3 ++
 lib/measurement.c    | 145 ++++++++++++++++++++++++++++++++++++++++++++++++---
 lib/measurement.h    |  18 ++++++-
 6 files changed, 178 insertions(+), 16 deletions(-)

(limited to 'lib')

diff --git a/lib/cluster_finite.c b/lib/cluster_finite.c
index f11a3ea..71396e0 100644
--- a/lib/cluster_finite.c
+++ b/lib/cluster_finite.c
@@ -62,14 +62,14 @@ v_t flip_cluster_finite(state_finite_t *s, v_t v0, q_t rot_ind, gsl_rng *r) {
           s->M[rot_s_old]--;
           s->M[rot_s_new]++;
 
-          s->E += - s->H[rot_s_new] + s->H[rot_s_old];
         } else {
           q_t diff_old = (s_old + s->q - sn) % s->q;
           q_t diff_new = (s_new + s->q - sn) % s->q;
 
           prob = s->J_probs[diff_new * s->q + diff_old];
 
-          s->E += - s->J[diff_new] + s->J[diff_old];
+          s->B[diff_old]--;
+          s->B[diff_new]++;
         }
 
         if (gsl_rng_uniform(r) < prob) { // and with probability ps[e]...
diff --git a/lib/cluster_finite.h b/lib/cluster_finite.h
index ad45ed3..701c197 100644
--- a/lib/cluster_finite.h
+++ b/lib/cluster_finite.h
@@ -38,7 +38,7 @@ typedef struct {
   double *H_probs;
   q_t *spins;
   q_t *R;
-  double E;
+  v_t *B;
   v_t *M;
 } state_finite_t;
 
diff --git a/lib/initial_finite.c b/lib/initial_finite.c
index f286dcc..fb120f0 100644
--- a/lib/initial_finite.c
+++ b/lib/initial_finite.c
@@ -58,9 +58,10 @@ state_finite_t *initial_finite_prepare_ising(D_t D, L_t L, double T, double *H)
   s->spins = (q_t *)calloc(s->nv, sizeof(q_t));
   s->R = initialize_R(2);
 
-  s->E = - ((double)s->ne) * s->J[0] - ((double)s->nv) * s->H[0];
   s->M = (v_t *)calloc(2, sizeof(v_t));
   s->M[0] = s->nv; // everyone starts in state 0, remember?
+  s->B = (v_t *)calloc(2, sizeof(v_t));
+  s->B[0] = s->ne;
 
   return s;
 }
@@ -98,9 +99,10 @@ state_finite_t *initial_finite_prepare_potts(D_t D, L_t L, q_t q, double T, doub
   s->spins = (q_t *)calloc(s->nv, sizeof(q_t));
   s->R = initialize_R(q);
 
-  s->E = - ((double)s->ne) * s->J[0] - ((double)s->nv) * s->H[0];
   s->M = (v_t *)calloc(q, sizeof(v_t));
   s->M[0] = s->nv; // everyone starts in state 0, remember?
+  s->B = (v_t *)calloc(q, sizeof(v_t));
+  s->B[0] = s->ne; // everyone starts in state 0, remember?
 
   return s;
 }
@@ -142,9 +144,10 @@ state_finite_t *initial_finite_prepare_clock(D_t D, L_t L, q_t q, double T, doub
   s->spins = (q_t *)calloc(s->nv, sizeof(q_t));
   s->R = initialize_R(q);
 
-  s->E = - ((double)s->ne) * s->J[0] - ((double)s->nv) * s->H[0];
   s->M = (v_t *)calloc(q, sizeof(v_t));
   s->M[0] = s->nv; // everyone starts in state 0, remember?
+  s->B = (v_t *)calloc(q, sizeof(v_t));
+  s->B[0] = s->ne; // everyone starts in state 0, remember?
 
   return s;
 }
@@ -189,13 +192,23 @@ state_finite_t *initial_finite_prepare_dgm(D_t D, L_t L, q_t q, double T, double
   s->spins = (q_t *)calloc(s->nv, sizeof(q_t));
   s->R = initialize_R(q);
 
-  s->E = - ((double)s->ne) * s->J[0] - ((double)s->nv) * s->H[0];
   s->M = (v_t *)calloc(q, sizeof(v_t));
   s->M[0] = s->nv; // everyone starts in state 0, remember?
 
   return s;
 }
 
+double state_finite_energy(state_finite_t *s) {
+  double E = 0;
+
+  for (q_t i = 0; i < s->q; i++) {
+    E += s->J[i] * s->B[i];
+    E += s->H[i] * s->M[i];
+  }
+
+  return -E;
+}
+
 void state_finite_free(state_finite_t *s) {
   graph_free(s->g);
   free(s->J);
@@ -205,6 +218,7 @@ void state_finite_free(state_finite_t *s) {
   free(s->spins);
   free(s->R);
   free(s->M);
+  free(s->B);
   free(s->transformations);
   free(s);
 }
diff --git a/lib/initial_finite.h b/lib/initial_finite.h
index 65414cd..542f923 100644
--- a/lib/initial_finite.h
+++ b/lib/initial_finite.h
@@ -7,6 +7,8 @@
 #include "dihedral.h"
 #include "cluster_finite.h"
 
+static char *finite_model_t_strings[] = {"ISING", "POTTS", "CLOCK", "DGM"};
+
 typedef enum {
   ISING,
   POTTS,
@@ -21,4 +23,5 @@ state_finite_t *initial_finite_prepare_dgm(D_t D, L_t L, q_t q, double T, double
 
 void state_finite_free(state_finite_t *s);
 
+double state_finite_energy(state_finite_t *s);
 
diff --git a/lib/measurement.c b/lib/measurement.c
index ad824f6..b30cf6b 100644
--- a/lib/measurement.c
+++ b/lib/measurement.c
@@ -1,6 +1,15 @@
 
+#include "convex.h"
 #include "measurement.h"
 
+meas_t *meas_initialize(count_t W) {
+  meas_t *m = (meas_t *)calloc(1, sizeof(meas_t));
+  m->W = W;
+  m->xx = (double *)calloc(2 * W + 1, sizeof(double));
+
+  return m;
+}
+
 double add_to_avg(double mx, double x, count_t n) {
   return mx * (n / (n + 1.0)) + x / (n + 1.0);
 }
@@ -10,24 +19,42 @@ void meas_update(meas_t *m, double x) {
 
   m->x = add_to_avg(m->x, x, n);
   m->x2 = add_to_avg(m->x2, pow(x, 2), n);
+  m->x4 = add_to_avg(m->x4, pow(x, 4), n);
 
   m->m2 = add_to_avg(m->m2, pow(x - m->x, 2), n);
   m->m4 = add_to_avg(m->m4, pow(x - m->x, 4), n);
 
-  /*
-  if (n > 1) {
-    double s2 = n / (n - 1.) * (m->x2 - pow(m->x, 2));
-    m->dx = sqrt(s2 / n);
-    m->c = s2;
-    m->dc = sqrt((m->m4 - (n - 3.)/(n - 1.) * pow(m->m2, 2)) / n);
+  dll_t *tmp_window = m->x_window;
+  dll_t *pos_save;
+  count_t t = 0;
+
+  while (tmp_window != NULL) {
+    m->xx[t] = add_to_avg(m->xx[t], x * (tmp_window->x), m->n - t - 1);
+    t++;
+    if (t == 2 * (m->W)) {
+      pos_save = tmp_window;
+    }
+    tmp_window = tmp_window->next;
   }
-  */
+
+  if (t == 2 * (m->W) + 1) {
+    if (2 * (m->W) + 1 == 1) {
+      free(m->x_window);
+      m->x_window = NULL;
+    } else {
+      free(pos_save->next);
+      pos_save->next = NULL;
+    }
+  }
+
+  stack_push_d(&(m->x_window), x);
 
   (m->n)++;
 }
 
 double meas_dx(const meas_t *m) {
-  return sqrt(1. / (m->n - 1.) * (m->x2 - pow(m->x, 2)));
+  return 2 * get_tau(m) * Cxx(m, 0) / m->n;
+//  return sqrt(1. / (m->n - 1.) * (m->x2 - pow(m->x, 2)));
 }
 
 double meas_c(const meas_t *m) {
@@ -74,3 +101,105 @@ double rho(const autocorr_t *o, count_t i) {
   return (o->OO[i] - pow(o->O, 2)) / (o->O2 - pow(o->O, 2));
 }
 
+double Cxx(const meas_t *m, count_t t) {
+  return m->xx[t] - pow(m->x, 2);
+}
+
+double rho_m(const meas_t *m, count_t t) {
+  return Cxx(m, t) / Cxx(m, 0);
+}
+
+double get_tau(const meas_t *m) {
+  double *Gammas = (double *)malloc((m->W + 1) * sizeof(double));
+
+  Gammas[0] = 1 + rho_m(m, 0);
+  for (uint64_t i = 0; i < m->W; i++) {
+    Gammas[1 + i] = rho_m(m, 2 * i + 1) + rho_m(m, 2 * i + 2);
+  }
+
+  uint64_t n;
+  for (n = 0; n < m->W + 1; n++) {
+    if (Gammas[n] <= 0) {
+      break;
+    }
+  }
+
+  double *conv_Gamma = get_convex_minorant(n, Gammas);
+
+  double tau = - 0.5;
+
+  for (uint64_t i = 0; i < n + 1; i++) {
+    tau += conv_Gamma[i];
+  }
+
+  free(Gammas);
+
+  return tau;
+}
+
+void print_meas(const meas_t *m, const char *sym, FILE *outfile) {
+  fprintf(outfile, "%s-><|n->%" PRIcount ",x->%.15f,x^2->%.15f,x^4->%.15f,xx->{", sym, m->n, m->x, m->x2, m->x4);
+  for (count_t i = 0; i < 2 * (m->W) + 1; i++) {
+    fprintf(outfile, "%.15f", m->xx[i]);
+    if (i < 2 * (m->W)) {
+      fprintf(outfile, ",");
+    }
+  }
+  fprintf(outfile, "}|>");
+}
+
+void print_vec_meas(q_t q, const meas_t **m, const char *sym, FILE *outfile) {
+  fprintf(outfile, "%s-><|n->{", sym);
+  for (q_t i = 0; i < q; i++) {
+    fprintf(outfile, "%" PRIcount, m[i]->n);
+    if (i < q - 1) {
+      fprintf(outfile, ",");
+    }
+  }
+  fprintf(outfile, "},x->{");
+  for (q_t i = 0; i < q; i++) {
+    fprintf(outfile, "%.15f", m[i]->x);
+    if (i < q - 1) {
+      fprintf(outfile, ",");
+    }
+  }
+  fprintf(outfile, "},x^2->{");
+  for (q_t i = 0; i < q; i++) {
+    fprintf(outfile, "%.15f", m[i]->x2);
+    if (i < q - 1) {
+      fprintf(outfile, ",");
+    }
+  }
+  fprintf(outfile, "},x^4->{");
+  for (q_t i = 0; i < q; i++) {
+    fprintf(outfile, "%.15f", m[i]->x4);
+    if (i < q - 1) {
+      fprintf(outfile, ",");
+    }
+  }
+  fprintf(outfile, "},xx->{");
+  for (q_t i = 0; i < q; i++) {
+    fprintf(outfile, "{");
+    for (count_t j = 0; j < 2 * (m[i]->W) + 1; j++) {
+      fprintf(outfile, "%.15f", m[i]->xx[j]);
+      if (j < 2 * (m[i]->W)) {
+        fprintf(outfile, ",");
+      }
+    }
+    fprintf(outfile, "}");
+    if (i < q - 1) {
+      fprintf(outfile, ",");
+    }
+  }
+  fprintf(outfile, "}|>");
+}
+
+void free_meas(meas_t *m) {
+  free(m->xx);
+  while (m->x_window != NULL) {
+    stack_pop_d(&(m->x_window));
+  }
+  free(m);
+}
+
+
diff --git a/lib/measurement.h b/lib/measurement.h
index eaa260b..d9bd52e 100644
--- a/lib/measurement.h
+++ b/lib/measurement.h
@@ -3,16 +3,21 @@
 
 #include <math.h>
 #include <stdlib.h>
+#include <stdio.h>
 
 #include "types.h"
 #include "stack.h"
 
 typedef struct {
-  uint64_t n;
+  count_t n;
   double x;
   double x2;
+  double x4;
   double m2;
   double m4;
+  count_t W;
+  double *xx;
+  dll_t *x_window;
 } meas_t;
 
 typedef struct {
@@ -36,3 +41,14 @@ void update_autocorr(autocorr_t *OO, double O);
 
 double rho(const autocorr_t *o, uint64_t i);
 
+void print_meas(const meas_t *m, const char *sym, FILE *outfile);
+void print_vec_meas(q_t q, const meas_t **m, const char *sym, FILE *outfile);
+
+void free_meas(meas_t *m);
+
+meas_t *meas_initialize(count_t W);
+
+double get_tau(const meas_t *m);
+
+double Cxx(const meas_t *m, count_t t);
+
-- 
cgit v1.2.3-70-g09d2