summaryrefslogtreecommitdiff
path: root/lib/cluster.h
blob: 2de17e585424c329ebc1cedf55104edcaeff3263 (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72

#pragma once

#include <assert.h>
#include <fftw3.h>
#include <float.h>
#include <gsl/gsl_randist.h>
#include <gsl/gsl_rng.h>
#include <inttypes.h>
#include <math.h>
#include <stdbool.h>
#include <string.h>
#include <sys/types.h>

#include "types.h"
#include "rand.h"
#include "stack.h"
#include "convex.h"
#include "graph.h"
#include "tree.h"
#include "measurement.h"
#include "orthogonal.h"
#include "dihedral.h"
#include "dihinf.h"

typedef struct {
  graph_t *g;
  q_t *spins;
  double T;
  double *J;
  double *H;
  double *J_probs;
  double *H_probs;
  dihedral_t *R;
  double E;
  v_t *M;
  q_t q;
} ising_state_t;

typedef struct {
  graph_t *g;
  h_t *spins;
  double T;
  double (*J)(h_t);
  double (*H)(double *, h_t);
  double *H_info;
  dihinf_t *R;
  double E;
  h_t M;
} dgm_state_t;

typedef struct {
  graph_t *g;
  double *spins;
  double T;
  double (*J)(double);
  double (*H)(q_t, double *, double *);
  double *H_info;
  double *R;
  double E;
  double *M;
  q_t n;
} vector_state_t;

v_t flip_cluster(ising_state_t *s, v_t v0, q_t s1, gsl_rng *r);

v_t flip_cluster_vector(vector_state_t *s, v_t v0, double *rot, gsl_rng *r);

v_t flip_cluster_dgm(dgm_state_t *s, v_t v0, h_t rot, gsl_rng *r);

graph_t *graph_add_ext(const graph_t *g);