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
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
|
#include "measurements.hpp"
#include <iostream>
#include <cstdio>
class badcycleException: public std::exception
{
virtual const char* what() const throw()
{
return "Could not find a valid cycle on the broken system.";
}
} badcycleex;
void update_distribution_file(std::string id, const std::vector<uint64_t>& data, std::string model_string) {
std::string filename = model_string + id + ".dat";
std::ifstream file(filename);
std::vector<uint64_t> data_old(data.size(), 0);
if (file.is_open()) {
for (unsigned i = 0; i < data.size(); i++) {
uint64_t num;
file >> num;
data_old[i] = num;
}
file.close();
}
std::ofstream file_out(filename);
for (unsigned i = 0; i < data.size(); i++) {
file_out <<std::fixed<< data_old[i] + data[i] << " ";
}
file_out.close();
}
template <class T>
void update_field_file(std::string id, const std::vector<T>& data, std::string model_string) {
std::string filename = model_string + id + ".dat";
std::ifstream file(filename);
std::vector<T> data_old(data.size(), 0);
if (file.is_open()) {
for (unsigned j = 0; j < data.size(); j++) {
file >> data_old[j];
}
file.close();
}
std::ofstream file_out(filename);
for (unsigned j = 0; j < data.size(); j++) {
file_out << data_old[j] + data[j] << " ";
}
file_out.close();
}
fftw_complex* data_transform(unsigned Mx, unsigned My, fftw_plan forward_plan, double *fftw_forward_in, fftw_complex *fftw_forward_out) {
fftw_execute(forward_plan);
fftw_complex* output = (fftw_complex*)malloc(Mx * (My / 2 + 1) * sizeof(fftw_complex));
for (unsigned i = 0; i < Mx * (My / 2 + 1); i++) {
output[i][0] = fftw_forward_out[i][0];
output[i][1] = fftw_forward_out[i][1];
}
return output;
}
template <class T>
void correlation(unsigned Mx, unsigned My, std::vector<std::vector<T>>& data, const fftw_complex* tx1, const fftw_complex* tx2, fftw_plan reverse_plan, fftw_complex *fftw_reverse_in, double *fftw_reverse_out) {
for (unsigned i = 0; i < Mx * (My / 2 + 1); i++) {
fftw_reverse_in[i][0] = tx1[i][0] * tx2[i][0] + tx1[i][1] * tx2[i][1];
fftw_reverse_in[i][1] = tx1[i][0] * tx2[i][1] - tx1[i][1] * tx2[i][0];
}
fftw_execute(reverse_plan);
for (unsigned j = 0; j < data.size(); j++) {
for (unsigned i = 0; i < (Mx / 2 + 1) * (My / 2 + 1); i++) {
data[j][i] += (T)pow(fftw_reverse_out[Mx * (i / (Mx / 2 + 1)) + i % (Mx / 2 + 1)] / (Mx * My), j + 1);
}
}
}
void autocorrelation2(double Lx, double Ly, unsigned Mx, unsigned My, std::vector<uint64_t>& data, const std::list<graph::coordinate>& pos, std::array<unsigned, 2> count) {
for (std::list<graph::coordinate>::const_iterator it1 = pos.begin(); it1 != pos.end(); it1++) {
for (std::list<graph::coordinate>::const_iterator it2 = it1; it2 != pos.end(); it2++) {
double Δx_tmp = fabs(it1->x - it2->x);
double Δx = Δx_tmp < Lx / 2 ? Δx_tmp : Lx - Δx_tmp;
double Δy_tmp = fabs(it1->y - it2->y);
double Δy = Δy_tmp < Ly / 2 ? Δy_tmp : Ly - Δy_tmp;
if (count[0] % 2 == 0) {
data[(unsigned)(Mx * (Δx / Lx)) + (Mx / 2) * (unsigned)(My * (Δy / Ly))]++;
} else {
data[(unsigned)(Mx * (Δy / Ly)) + (Mx / 2) * (unsigned)(My * (Δx / Lx))]++;
}
}
}
}
void correlation2(double Lx, double Ly, unsigned Mx, unsigned My, std::vector<uint64_t>& data, const std::list<graph::coordinate>& pos1, const std::list<graph::coordinate>& pos2) {
for (std::list<graph::coordinate>::const_iterator it1 = pos1.begin(); it1 != pos1.end(); it1++) {
for (std::list<graph::coordinate>::const_iterator it2 = pos2.begin(); it2 != pos2.end(); it2++) {
double Δx_tmp = fabs(it1->x - it2->x);
double Δx = Δx_tmp < Lx / 2 ? Δx_tmp : Lx - Δx_tmp;
double Δy_tmp = fabs(it1->y - it2->y);
double Δy = Δy_tmp < Ly / 2 ? Δy_tmp : Ly - Δy_tmp;
data[floor((Mx * Δx) / Lx) + (Mx / 2) * floor((My * Δy) / Ly)]++;
}
}
}
template <class T>
void autocorrelation(unsigned Mx, unsigned My, std::vector<std::vector<T>>& out_data, fftw_plan forward_plan, double *fftw_forward_in, fftw_complex *fftw_forward_out, fftw_plan reverse_plan, fftw_complex *fftw_reverse_in, double *fftw_reverse_out) {
fftw_execute(forward_plan);
for (unsigned i = 0; i < My * (Mx / 2 + 1); i++) {
fftw_reverse_in[i][0] = pow(fftw_forward_out[i][0], 2) + pow(fftw_forward_out[i][1], 2);
fftw_reverse_in[i][1] = 0.0;
}
fftw_execute(reverse_plan);
for (unsigned j = 0; j < out_data.size(); j++) {
for (unsigned i = 0; i < (Mx / 2 + 1) * (My / 2 + 1); i++) {
out_data[j][i] += (T)pow(fftw_reverse_out[Mx * (i / (Mx / 2 + 1)) + i % (Mx / 2 + 1)] / (Mx * My), j + 1);
}
}
}
unsigned edge_r_to_ind(graph::coordinate r, double Lx, double Ly, unsigned Mx, unsigned My) {
return floor((Mx * r.x) / Lx) + Mx * floor((My * r.y) / Ly);
}
ma::ma(unsigned n, double a, double beta, double weight) :
sn(2 * n),
ss(2 * n),
sm(2 * n),
sl(2 * n),
sb(n + 1),
sd(3 * n),
sa(3 * n),
sA(3 * n),
si(10000),
sI(10000),
cn(pow((unsigned)sqrt(n), 2)),
cs(pow((unsigned)sqrt(n), 2)),
cm(pow((unsigned)sqrt(n), 2)),
cl(pow((unsigned)sqrt(n), 2)),
cb(pow((unsigned)sqrt(n), 2)),
ca(pow((unsigned)sqrt(n), 2)),
cA(pow((unsigned)sqrt(n), 2)),
cp(pow((unsigned)sqrt(n), 2)),
cq(pow((unsigned)sqrt(n), 2))
{
if (beta != 0.0) {
model_string = "fracture_" + std::to_string(n) + "_" + std::to_string(a) + "_" + std::to_string(beta) + "_" + std::to_string(weight) + "_";
} else {
model_string = "fracture_" + std::to_string(n) + "_" + std::to_string(a) + "_INF_" + std::to_string(weight) + "_";
}
}
ma::ma(unsigned Lx, unsigned Ly, double beta, double weight) :
sn(Lx * Ly / 2),
ss(Lx * Ly / 2),
sm(Lx * Ly / 2),
sl(Lx * Ly / 2),
sb(Lx * Ly / 2 + 1),
sd(Lx * Ly),
sa(Lx * Ly),
sA(Lx * Ly),
si(10000),
sI(10000),
cn(Lx * Ly),
cs(Lx * Ly),
cm(Lx * Ly),
cl(Lx * Ly),
cb(Lx * Ly),
ca(Lx * Ly),
cA(Lx * Ly),
cp(Lx * Ly),
cq(Lx * Ly)
{
if (beta != 0.0) {
model_string = "fracture_" + std::to_string(Lx) + "_" + std::to_string(Ly) + "_" + std::to_string(beta) + "_" + std::to_string(weight) + "_";
} else {
model_string = "fracture_" + std::to_string(Lx) + "_" + std::to_string(Ly) + "_INF_" + std::to_string(weight) + "_";
}
}
ma::~ma() {
update_distribution_file("sn", sn, model_string);
update_distribution_file("ss", ss, model_string);
update_distribution_file("sm", sm, model_string);
update_distribution_file("sl", sl, model_string);
update_distribution_file("sb", sb, model_string);
update_distribution_file("sd", sd, model_string);
update_distribution_file("sa", sa, model_string);
update_distribution_file("sA", sA, model_string);
update_distribution_file("si", si, model_string);
update_distribution_file("sI", sI, model_string);
update_field_file("cl", cl, model_string);
update_field_file("cm", cm, model_string);
update_field_file("cs", cs, model_string);
update_field_file("cn", cn, model_string);
update_field_file("cb", cb, model_string);
update_field_file("ca", ca, model_string);
update_field_file("cA", cA, model_string);
update_field_file("cp", cp, model_string);
update_field_file("cq", cq, model_string);
}
void ma::pre_fracture(const network&) {
lv = std::numeric_limits<long double>::lowest();
avalanches = {};
num = 0;
}
void ma::bond_broken(const network& net, const current_info& cur, unsigned i) {
long double c = net.thresholds[i] - logl(cur.currents[i]);
if (c > lv) {
lv = c;
avalanches.push_back({i});
} else {
avalanches.back().push_back(i);
}
for (unsigned j = 0; j < cur.currents.size(); j++) {
if (logl(cur.currents[j]) >= -100 && logl(cur.currents[j]) < 0 && (!net.backbone[j] || j == i)) {
si[(unsigned)(10000 * (logl(cur.currents[j]) + 100) / 100)]++;
} else if (logl(cur.currents[j]) >= -100 && logl(cur.currents[j]) < 0 && net.backbone[j] && j != i) {
sI[(unsigned)(10000 * (logl(cur.currents[j]) + 100) / 100)]++;
}
}
num++;
}
void ma::post_fracture(network &n) {
auto crack = find_minimal_crack(n, avalanches.back().back());
std::list<graph::coordinate> cl_cs;
sl[crack.second.size() - 1]++;
for (unsigned e : crack.second) {
cl_cs.push_back(n.G.dual_edges[e].r);
}
autocorrelation2(n.G.L.x, n.G.L.y, 2 * sqrt(cl.size()), 2 * sqrt(cl.size()), cl, cl_cs, crack.first);
std::vector<std::list<graph::coordinate>> components(n.G.dual_vertices.size());
for (unsigned i = 0; i < n.G.dual_vertices.size(); i++) {
components[n.C.findroot(i)].push_back(n.G.dual_vertices[i].r);
}
unsigned crack_component = n.C.findroot(n.G.dual_edges[avalanches.back().back()].v[0]);
for (unsigned i = 0; i < n.G.dual_vertices.size(); i++) {
if (components[i].size() > 0) {
if (i != crack_component) {
sm[components[i].size() - 1]++;
autocorrelation2(n.G.L.x, n.G.L.y, 2 * sqrt(cm.size()), 2 * sqrt(cm.size()), cm, components[i], crack.first);
autocorrelation2(n.G.L.x, n.G.L.y, 2 * sqrt(cp.size()), 2 * sqrt(cp.size()), cp, components[i], {0, 1});
} else {
ss[components[i].size() - 1]++;
autocorrelation2(n.G.L.x, n.G.L.y, 2 * sqrt(cs.size()), 2 * sqrt(cs.size()), cs, components[i], crack.first);
}
sn[components[i].size() - 1]++;
autocorrelation2(n.G.L.x, n.G.L.y, 2 * sqrt(cn.size()), 2 * sqrt(cn.size()), cn, components[i], crack.first);
}
}
std::vector<bool> vertex_in(n.G.vertices.size());
for (unsigned i = 0; i < n.G.edges.size(); i++) {
if (!n.backbone[i]) {
vertex_in[n.G.edges[i].v[0]] = true;
vertex_in[n.G.edges[i].v[1]] = true;
}
}
unsigned bb_size = 0;
std::list<graph::coordinate> cb_co;
for (unsigned i = 0; i < n.G.vertices.size(); i++) {
if (vertex_in[i]) {
bb_size++;
cb_co.push_back(n.G.vertices[i].r);
}
}
sb[bb_size]++;
autocorrelation2(n.G.L.x, n.G.L.y, 2 * sqrt(cb.size()), 2 * sqrt(cb.size()), cb, cb_co, crack.first);
auto av_it = avalanches.rbegin();
av_it++;
while (av_it != avalanches.rend()) {
sa[(*av_it).size() - 1]++;
std::list<graph::coordinate> ca_co;
for (unsigned e : (*av_it)) {
ca_co.push_back(n.G.edges[e].r);
}
autocorrelation2(n.G.L.x, n.G.L.y, 2 * sqrt(ca.size()), 2 * sqrt(ca.size()), ca, ca_co, crack.first);
autocorrelation2(n.G.L.x, n.G.L.y, 2 * sqrt(cq.size()), 2 * sqrt(cq.size()), cq, ca_co, {0,1});
av_it++;
}
sA[avalanches.back().size() - 1]++;
std::list<graph::coordinate> cA_co;
for (unsigned e : avalanches.back()) {
cA_co.push_back(n.G.edges[e].r);
}
autocorrelation2(n.G.L.x, n.G.L.y, 2 * sqrt(cA.size()), 2 * sqrt(cA.size()), cA, cA_co, crack.first);
sd[num - 1]++;
}
|