diff options
| author | Jaron Kent-Dobias <jaron@kent-dobias.com> | 2018-06-29 15:43:10 -0400 | 
|---|---|---|
| committer | Jaron Kent-Dobias <jaron@kent-dobias.com> | 2018-06-29 15:43:10 -0400 | 
| commit | 78d8de381f0b1e99ad98364709cbf876689628b2 (patch) | |
| tree | 73605132c64124b6ad3a057a1ffec8d40d356832 | |
| parent | 3eb67e3bca774eb0441db60158e1968ad901273b (diff) | |
| download | c++-78d8de381f0b1e99ad98364709cbf876689628b2.tar.gz c++-78d8de381f0b1e99ad98364709cbf876689628b2.tar.bz2 c++-78d8de381f0b1e99ad98364709cbf876689628b2.zip  | |
completely removed measurement during the simulation, opting to just save binary data points to files throughout
| -rw-r--r-- | lib/cluster_finite.c | 4 | ||||
| -rw-r--r-- | lib/cluster_finite.h | 2 | ||||
| -rw-r--r-- | lib/initial_finite.c | 22 | ||||
| -rw-r--r-- | lib/initial_finite.h | 3 | ||||
| -rw-r--r-- | lib/measurement.c | 145 | ||||
| -rw-r--r-- | lib/measurement.h | 18 | ||||
| -rw-r--r-- | src/wolff_finite.c | 398 | 
7 files changed, 242 insertions, 350 deletions
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); + diff --git a/src/wolff_finite.c b/src/wolff_finite.c index 47fcc88..e41c326 100644 --- a/src/wolff_finite.c +++ b/src/wolff_finite.c @@ -1,93 +1,60 @@ +#include <time.h>  #include <getopt.h>  #include <initial_finite.h>  int main(int argc, char *argv[]) { -  L_t L = 128;    count_t N = (count_t)1e7; -  count_t min_runs = 10; -  count_t n = 3; + +  finite_model_t model = ISING; +    q_t q = 2;    D_t D = 2; +  L_t L = 128;    double T = 2.26918531421;    double *J = (double *)calloc(MAX_Q, sizeof(double));    J[0] = 1.0;    double *H = (double *)calloc(MAX_Q, sizeof(double)); -  double eps = 0; -  bool silent = false; -  bool snapshots = false; -  bool snapshot = false; -  bool record_autocorrelation = false; -  bool record_distribution = false; -  count_t W = 10; -  count_t ac_skip = 1; -  finite_model_t model = ISING; +  bool silent = false;    int opt;    q_t J_ind = 0;    q_t H_ind = 0; -  while ((opt = getopt(argc, argv, "N:n:D:L:q:T:J:H:m:e:IpsSPak:W:drt:")) != -1) { +  while ((opt = getopt(argc, argv, "N:t:q:D:L:T:J:H:s")) != -1) {      switch (opt) { -    case 'N': +    case 'N': // number of steps        N = (count_t)atof(optarg);        break; -    case 'n': -      n = (count_t)atof(optarg); +    case 't': // type of simulation +      model = (finite_model_t)atoi(optarg); +      break; +    case 'q': // number of states, if relevant +      q = atoi(optarg);        break; -    case 'D': +    case 'D': // dimension        D = atoi(optarg);        break; -    case 'L': +    case 'L': // linear size        L = atoi(optarg);        break; -    case 'q': -      q = atoi(optarg); -      break; -    case 'T': +    case 'T': // temperature         T = atof(optarg);        break; -    case 'J': +    case 'J': // couplings, if relevant. nth call couples states i and i + n        J[J_ind] = atof(optarg);        J_ind++;        break; -    case 'H': +    case 'H': // external field. nth call couples to state n        H[H_ind] = atof(optarg);        H_ind++;        break; -    case 'm': -      min_runs = atoi(optarg); -      break; -    case 'e': -      eps = atof(optarg); -      break; -    case 's': +    case 's': // don't print anything during simulation. speeds up slightly        silent = true;        break; -    case 'S': -      snapshots = true; -      break; -    case 'P': -      snapshot = true; -      break; -    case 'a': -      record_autocorrelation = true; -      break; -    case 'k': -      ac_skip = (count_t)atof(optarg); -      break; -    case 'W': -      W = (count_t)atof(optarg); -      break; -    case 'd': -      record_distribution = true; -      break; -    case 't': -      model = (finite_model_t)atoi(optarg); -      break;      default:        exit(EXIT_FAILURE);      } @@ -95,9 +62,6 @@ int main(int argc, char *argv[]) {    state_finite_t *s; -  gsl_rng *r = gsl_rng_alloc(gsl_rng_mt19937); -  gsl_rng_set(r, rand_seed()); -    switch (model) {      case ISING:        s = initial_finite_prepare_ising(D, L, T, H);  @@ -113,318 +77,84 @@ int main(int argc, char *argv[]) {        break;      default:        printf("Not a valid model!\n"); -      return 1; +      free(J); +      free(H); +      exit(EXIT_FAILURE);    }    free(J);    free(H); +  // initialize random number generator +  gsl_rng *r = gsl_rng_alloc(gsl_rng_mt19937); +  gsl_rng_set(r, rand_seed()); -  double diff = 1e31; -  count_t n_runs = 0; -  count_t n_steps = 0; +  unsigned long timestamp = (unsigned long)time(NULL); -  meas_t *E, *clust, **M, **sE, ***sM; +  char *filename_M = (char *)malloc(256 * sizeof(char)); +  char *filename_B = (char *)malloc(256 * sizeof(char)); +  char *filename_S = (char *)malloc(256 * sizeof(char)); -  M = (meas_t **)malloc(q * sizeof(meas_t *)); +  sprintf(filename_M, "wolff_%s_%" PRIq "_%" PRID "_%" PRIL "_%.8f", finite_model_t_strings[model], q, D, L, T);    for (q_t i = 0; i < q; i++) { -    M[i] = (meas_t *)calloc(1, sizeof(meas_t)); +    sprintf(filename_M + strlen(filename_M), "_%.8f", s->H[i]);    } -  E = calloc(1, sizeof(meas_t)); -  clust = calloc(1, sizeof(meas_t)); +  strcpy(filename_B, filename_M); +  strcpy(filename_S, filename_M); -  sE = (meas_t **)malloc(q * sizeof(meas_t *)); -  sM = (meas_t ***)malloc(q * sizeof(meas_t **)); +  sprintf(filename_M + strlen(filename_M), "_%lu_M.dat", timestamp); +  sprintf(filename_B + strlen(filename_B), "_%lu_B.dat", timestamp); +  sprintf(filename_S + strlen(filename_S), "_%lu_S.dat", timestamp); -  for (q_t i = 0; i < q; i++) { -    sE[i] = (meas_t *)calloc(1, sizeof(meas_t)); -    sM[i] = (meas_t **)malloc(q * sizeof(meas_t *)); -    for (q_t j = 0; j < q; j++) { -      sM[i][j] = (meas_t *)calloc(1, sizeof(meas_t)); -    } -  } - -  count_t *freqs = (count_t *)calloc(q, sizeof(count_t)); -  q_t cur_M = 0; +  FILE *outfile_M = fopen(filename_M, "wb"); +  FILE *outfile_B = fopen(filename_B, "wb"); +  FILE *outfile_S = fopen(filename_S, "wb"); -  autocorr_t *autocorr; -  if (record_autocorrelation) { -    autocorr = (autocorr_t *)calloc(1, sizeof(autocorr_t)); -    autocorr->W = 2 * W + 1; -    autocorr->OO = (double *)calloc(2 * W + 1, sizeof(double)); -  } +  free(filename_M); +  free(filename_B); +  free(filename_S); -  count_t *mag_dist; -  if (record_distribution) { -    mag_dist = (count_t *)calloc(s->nv + 1, sizeof(count_t)); -  } +  v_t cluster_size = 0;    if (!silent) printf("\n"); -  while (((diff > eps || diff != diff) && n_runs < N) || n_runs < min_runs) { -    if (!silent) printf("\033[F\033[JWOLFF: sweep %" PRIu64 -           ", dH/H = %.4f, dM/M = %.4f, dC/C = %.4f, dX/X = %.4f, cps: %.1f\n", -           n_runs, fabs(meas_dx(E) / E->x), meas_dx(M[0]) / M[0]->x, meas_dc(E) / meas_c(E), meas_dc(M[0]) / meas_c(M[0]), s->nv / clust->x); +  for (count_t steps = 0; steps < N; steps++) { +    if (!silent) printf("\033[F\033[JWOLFF: sweep %" PRIu64 " / %" PRIu64 ": E = %.2f, B_0 = %" PRIv ", M_0 = %" PRIv ", S = %" PRIv "\n", steps, N, state_finite_energy(s), s->B[0], s->M[0], cluster_size); -    count_t n_flips = 0; - -    while (n_flips / s->nv < n) { -      v_t v0 = gsl_rng_uniform_int(r, s->nv); -      R_t step; +    v_t v0 = gsl_rng_uniform_int(r, s->nv); +    R_t step; -      bool changed = false; -      while (!changed) { -        step = gsl_rng_uniform_int(r, s->n_transformations); -        if (symmetric_act(s->transformations + q * step, s->spins[v0]) != s->spins[v0]) { -          changed = true; -        } +    bool changed = false; +    while (!changed) { +      step = gsl_rng_uniform_int(r, s->n_transformations); +      if (symmetric_act(s->transformations + q * step, s->spins[v0]) != s->spins[v0]) { +        changed = true;        } - -      v_t tmp_flips = flip_cluster_finite(s, v0, step, r); -      n_flips += tmp_flips; - -      if (n_runs > 0) { -        n_steps++; -        meas_update(clust, tmp_flips); - -        if (record_autocorrelation && n_steps % ac_skip == 0) { -          update_autocorr(autocorr, s->E); -        } - -      } -      } -    for (q_t i = 0; i < q; i++) { -      meas_update(M[i], s->M[i]); -    } -    meas_update(E, s->E); +    cluster_size = flip_cluster_finite(s, v0, step, r); -    q_t n_at_max = 0; -    q_t max_M_i = 0; -    v_t max_M = 0; +    // v_t is never going to be bigger than 32 bits, but since it's specified +    // as a fast time many machines will actually have it be 64 bits. we cast +    // it down here to halve space. -    for (q_t i = 0; i < q; i++) { -      if (s->M[i] > max_M) { -        n_at_max = 1; -        max_M_i = i; -        max_M = s->M[i]; -      } else if (s->M[i] == max_M) { -        n_at_max++; -      } +    for (q_t i = 0; i < q - 1; i++) { +      fwrite(&(s->M[i]), sizeof(uint32_t), 1, outfile_M); // if we know the occupation of the first q - 1 states, we know the occupation of the last +      fwrite(&(s->B[i]), sizeof(uint32_t), 1, outfile_B); // if we know the occupation of the first q - 1 states, we know the occupation of the last      } -    if (record_distribution) { -      mag_dist[s->M[0]]++; -    } - -    if (n_at_max == 1) { -      for (q_t i = 0; i < q; i++) { -        meas_update(sM[max_M_i][i], s->M[i]); -      } -      meas_update(sE[max_M_i], s->E); -      freqs[max_M_i]++; -    } +    fwrite(&cluster_size, sizeof(uint32_t), 1, outfile_S); -    diff = fabs(meas_dx(clust) / clust->x); - -    n_runs++;    }    if (!silent) {      printf("\033[F\033[J");    } -  printf("WOLFF: sweep %" PRIu64 -         ", dH/H = %.4f, dM/M = %.4f, dC/C = %.4f, dX/X = %.4f, cps: %.1f\n", -         n_runs, fabs(meas_dx(E) / E->x), meas_dx(M[0]) / M[0]->x, meas_dc(E) / meas_c(E), meas_dc(M[0]) / meas_c(M[0]), s->nv / clust->x); - -  if (snapshots) { -    FILE *snapfile = fopen("snapshots.m", "a"); -    fprintf(snapfile, "\n"); -  } - -  if (snapshot) { -    q_t *R_inv = symmetric_invert(q, s->R); -    FILE *snapfile = fopen("snapshot.m", "a"); -    fprintf(snapfile, "{{"); -    for (L_t i = 0; i < L; i++) { -      fprintf(snapfile, "{"); -      for (L_t j = 0; j < L; j++) { -        fprintf(snapfile, "%" PRIq, symmetric_act(R_inv, s->spins[L * i + j])); -        if (j != L - 1) { -          fprintf(snapfile, ","); -        } -      } -      fprintf(snapfile, "}"); -      if (i != L - 1) { -        fprintf(snapfile, ","); -      } -    } -    fprintf(snapfile, "}}\n"); -    fclose(snapfile); -  } - -  double tau = 0; -  int tau_failed = 0; - -  if (record_autocorrelation) { -    double *Gammas = (double *)malloc((W + 1) * sizeof(double)); +  printf("WOLFF: sweep %" PRIu64 " / %" PRIu64 ": E = %.2f, B_0 = %" PRIv ", M_0 = %" PRIv ", S = %" PRIv "\n", N, N, state_finite_energy(s), s->B[0], s->M[0], cluster_size); -    Gammas[0] = 1 + rho(autocorr, 0); -    for (uint64_t i = 0; i < W; i++) { -      Gammas[1 + i] = rho(autocorr, 2 * i + 1) + rho(autocorr, 2 * i + 2); -    }  +  fclose(outfile_M); +  fclose(outfile_B); +  fclose(outfile_S); -    uint64_t n; -    for (n = 0; n < W + 1; n++) { -      if (Gammas[n] <= 0) { -        break; -      } -    } - -    if (n == W + 1) { -      printf("WARNING: correlation function never hit the noise floor.\n"); -      tau_failed = 1; -    } - -    if (n < 2) { -      printf("WARNING: correlation function only has one nonnegative term.\n"); -      tau_failed = 2; -    } - -    double *conv_Gamma = get_convex_minorant(n, Gammas); - -    double ttau = - 0.5; - -    for (uint64_t i = 0; i < n + 1; i++) { -      ttau += conv_Gamma[i]; -    } -     -    tau = ttau * ac_skip * clust->x / s->nv; -     -    free(Gammas); -    free(autocorr->OO); -    while (autocorr->Op != NULL) { -      stack_pop_d(&(autocorr->Op)); -    } -    free(autocorr); -  } - -  if (tau_failed) { -    //tau = 0; -  } - -  { -  FILE *outfile = fopen("out.m", "a"); - -  fprintf(outfile, "<|N->%" PRIcount ",n->%" PRIcount ",D->%" PRID ",L->%" PRIL ",q->%" PRIq ",T->%.15f,J->{", N, n, D, L, q, T); -  for (q_t i = 0; i < q; i++) { -    fprintf(outfile, "%.15f", s->J[i]); -    if (i != q-1) { -      fprintf(outfile, ","); -    } -  } -  fprintf(outfile, "},H->{"); -  for (q_t i = 0; i < q; i++) { -    fprintf(outfile, "%.15f", s->H[i]); -    if (i != q-1) { -      fprintf(outfile, ","); -    } -  } -  fprintf(outfile, "},E->%.15f,\\[Delta]E->%.15f,C->%.15f,\\[Delta]C->%.15f,M->{", E->x / s->nv, meas_dx(E) / s->nv, meas_c(E) / s->nv, meas_dc(E) / s->nv); -  for (q_t i = 0; i < q; i++) { -    fprintf(outfile, "%.15f", M[i]->x / s->nv); -    if (i != q-1) { -      fprintf(outfile, ","); -    } -  } -  fprintf(outfile, "},\\[Delta]M->{"); -  for (q_t i = 0; i < q; i++) { -    fprintf(outfile, "%.15f", meas_dx(M[i]) / s->nv); -    if (i != q-1) { -      fprintf(outfile, ","); -    } -  } -  fprintf(outfile, "},\\[Chi]->{"); -  for (q_t i = 0; i < q; i++) { -    fprintf(outfile, "%.15f", meas_c(M[i]) / s->nv); -    if (i != q-1) { -      fprintf(outfile, ","); -    } -  } -  fprintf(outfile, "},\\[Delta]\\[Chi]->{"); -  for (q_t i = 0; i < q; i++) { -    fprintf(outfile, "%.15f", meas_dc(M[i]) / s->nv); -    if (i != q-1) { -      fprintf(outfile, ","); -    } -  } -  for (q_t i = 0; i < q; i++) { -    fprintf(outfile, "},Subscript[E,%" PRIq "]->%.15f,Subscript[\\[Delta]E,%" PRIq "]->%.15f,Subscript[C,%" PRIq "]->%.15f,Subscript[\\[Delta]C,%" PRIq "]->%.15f,Subscript[M,%" PRIq "]->{", i, sE[i]->x / s->nv, i, meas_dx(sE[i]) / s->nv, i, meas_c(sE[i]) / s->nv, i, meas_dc(sE[i]) / s->nv, i); -    for (q_t j = 0; j < q; j++) { -      fprintf(outfile, "%.15f", sM[i][j]->x / s->nv); -      if (j != q-1) { -        fprintf(outfile, ","); -      } -    } -    fprintf(outfile, "},Subscript[\\[Delta]M,%" PRIq "]->{", i); -    for (q_t j = 0; j < q; j++) { -      fprintf(outfile, "%.15f", meas_dx(sM[i][j]) / s->nv); -      if (j != q-1) { -        fprintf(outfile, ","); -      } -    } -    fprintf(outfile, "},Subscript[\\[Chi],%" PRIq "]->{", i); -    for (q_t j = 0; j < q; j++) { -      fprintf(outfile, "%.15f", meas_c(sM[i][j]) / s->nv); -      if (j != q-1) { -        fprintf(outfile, ","); -      } -    } -    fprintf(outfile, "},Subscript[\\[Delta]\\[Chi],%" PRIq "]->{", i); -    for (q_t j = 0; j < q; j++) { -      fprintf(outfile, "%.15f", meas_dc(sM[i][j]) / s->nv); -      if (j != q-1) { -        fprintf(outfile, ","); -      } -    } -  } -  fprintf(outfile,"}"); -  for (q_t i = 0; i < q; i++) { -    fprintf(outfile, ",Subscript[f,%" PRIq "]->%.15f,Subscript[\\[Delta]f,%" PRIq "]->%.15f", i, (double)freqs[i] / (double)n_runs, i, sqrt(freqs[i]) / (double)n_runs); -  } -  fprintf(outfile, ",Subscript[n,\"clust\"]->%.15f,Subscript[\\[Delta]n,\"clust\"]->%.15f,Subscript[m,\"clust\"]->%.15f,Subscript[\\[Delta]m,\"clust\"]->%.15f,\\[Tau]->%.15f,\\[Tau]s->%d", clust->x / s->nv, meas_dx(clust) / s->nv, meas_c(clust) / s->nv, meas_dc(clust) / s->nv,tau,tau_failed); -  if (record_distribution) { -    fprintf(outfile, ",S->{"); -    for (v_t i = 0; i < s->nv + 1; i++) { -      fprintf(outfile, "%" PRIcount, mag_dist[i]); -      if (i != s->nv) { -        fprintf(outfile, ","); -      } -    } -    fprintf(outfile, "}"); -    free(mag_dist); -  } -  fprintf(outfile, "|>\n"); - -  fclose(outfile); -  } - -  free(E); -  free(clust); -  for (q_t i = 0; i < q; i++) { -    free(M[i]); -    for (q_t j = 0; j < q; j++) { -      free(sM[i][j]); -    } -    free(sM[i]); -  } -  free(M); -  free(sM); -  for (q_t i = 0; i < q; i++) { -    free(sE[i]); -  } -  free(freqs); -  free(sE);    state_finite_free(s);    gsl_rng_free(r);  | 
