HILA
Loading...
Searching...
No Matches
polyakov_surface.h
1#ifndef POLYAKOV_SURFACE_H_
2#define POLYAKOV_SURFACE_H_
3#include "hila.h"
4#include "utility.h"
5#include <stdio.h>
6#include <stdlib.h>
7
8void spectraldensity_surface(std::vector<float> &surf, std::vector<double> &npow,
9 std::vector<int> &hits) {
10
11 // do fft for the surface
12 static bool first = true;
13 static Complex<double> *buf;
14 static fftw_plan fftwplan;
15
16 int area = lattice.size(e_x) * lattice.size(e_y);
17
18 if (first) {
19 first = false;
20
21 buf = (Complex<double> *)fftw_malloc(sizeof(Complex<double>) * area);
22
23 // note: we had x as the "fast" dimension, but fftw wants the 2nd dim to be
24 // the "fast" one. thus, first y, then x.
25 fftwplan = fftw_plan_dft_2d(lattice.size(e_y), lattice.size(e_x), (fftw_complex *)buf,
26 (fftw_complex *)buf, FFTW_FORWARD, FFTW_ESTIMATE);
27 }
28
29 for (int i = 0; i < area; i++) {
30 buf[i] = surf[i];
31 }
32
33 fftw_execute(fftwplan);
34
35 int pow_size = npow.size();
36
37 for (int i = 0; i < area; i++) {
38 int x = i % lattice.size(e_x);
39 int y = i / lattice.size(e_x);
40 x = (x <= lattice.size(e_x) / 2) ? x : (lattice.size(e_x) - x);
41 y = (y <= lattice.size(e_y) / 2) ? y : (lattice.size(e_y) - y);
42
43 int k = x * x + y * y;
44 if (k < pow_size) {
45 npow[k] += buf[i].squarenorm() / (area * area);
46 hits[k]++;
47 }
48 }
49}
50
51void wrap_surface(std::vector<float> &surface) {
52 float sum = std::accumulate(surface.begin(), surface.end(), 0.0f);
53 float average = sum / surface.size();
54 int half_lattice_z = lattice.size(e_z) / 2;
55 auto [min_it, max_it] = std::minmax_element(surface.begin(), surface.end());
56 float min_value = *min_it;
57 float max_value = *max_it;
58 if (abs(max_value - min_value) > lattice.size(e_z) * 0.9) {
59 for (auto &value : surface) {
60 if (average < half_lattice_z && value > half_lattice_z) {
61 value -= lattice.size(e_z);
62 } else if (average > half_lattice_z && value < half_lattice_z) {
63 value += lattice.size(e_z);
64 }
65 }
66 }
67}
68
69template <typename T>
70void measure_polyakov_field(const Field<T> &Ut, Field<Complex<float>> &polyakov_field) {
71 Field<T> polyakov = Ut;
72
73 // mult links so that polyakov[X.dir == 0] contains the polyakov loop
74 for (int plane = lattice.size(e_t) - 2; plane >= 0; plane--) {
75
76 // safe_access(polyakov) pragma allows the expression below, otherwise
77 // hilapp would reject it because X and X+dir can refer to the same
78 // site on different "iterations" of the loop. However, here this
79 // is restricted on single dir-plane so it works but we must tell it to
80 // hilapp.
81
82#pragma hila safe_access(polyakov)
83 onsites(ALL) {
84 if (X.coordinate(e_t) == plane) {
85 polyakov[X] = Ut[X] * polyakov[X + e_t];
86 }
87 }
88 }
89
90 onsites(ALL) if (X.coordinate(e_t) == 0) {
91 polyakov_field[X] = trace(polyakov[X]);
92 }
93}
94
95////////////////////////////////////////////////////////////////////
96template <typename T>
97void smear_polyakov_field(Field<T> &polyakov_field, int nsmear, float smear_coeff, double twist) {
98 if (nsmear > 0) {
99 Field<T> pl2 = 0;
100
101 for (int i = 0; i < nsmear; i++) {
102 onsites(ALL) if (X.coordinate(e_t) == 0) {
103
104 pl2[X] = polyakov_field[X] +
105 smear_coeff * (polyakov_field[X + e_x] + polyakov_field[X - e_x] +
106 polyakov_field[X + e_y] + polyakov_field[X - e_y]);
107 if (X.coordinate(e_z) == 0)
108 pl2[X] +=
109 smear_coeff * (expi(-2 * M_PI * twist / NCOLOR) * polyakov_field[X + e_z] +
110 polyakov_field[X - e_z]);
111 else if (X.coordinate(e_z) == 1)
112 pl2[X] +=
113 smear_coeff * (polyakov_field[X + e_z] +
114 expi(2 * M_PI * twist / NCOLOR) * polyakov_field[X - e_z]);
115 else
116 pl2[X] += smear_coeff * (polyakov_field[X + e_z] + polyakov_field[X - e_z]);
117 }
118
119 onsites(ALL) if (X.coordinate(e_t) == 0) {
120 polyakov_field[X] = pl2[X] / (1 + 6 * smear_coeff);
121 }
122 }
123 }
124}
125
126/**
127 * @brief Polyakov average and value over z-index
128 *
129 * @tparam T
130 * @param polyakov_field
131 * @param surface_origin_profile z-index of Polyakov field at \f$x=0\f$ and
132 * \f$y=0\f$
133 * @param surface_average_profile z-index Polyakov field sum
134 */
135template <typename T>
136void measure_polyakov_profile(Field<T> &polyakov_field, std::vector<T> &surface_origin_profile,
137 std::vector<T> &surface_average_profile) {
138 ReductionVector<T> p_surface_average(lattice.size(e_z)), p_origin(lattice.size(e_z));
139 p_surface_average.allreduce(false);
140 p_origin.allreduce(false);
141 onsites(ALL) if (X.coordinate(e_t) == 0) {
142 p_surface_average[X.z()] += polyakov_field[X];
143 if (X.x() == 0 && X.y() == 0)
144 p_origin[X.z()] += polyakov_field[X];
145 }
146 double inverse_surface_are = 1.0 / (lattice.size(e_x) * lattice.size(e_y));
147 surface_origin_profile = p_origin.vector();
148 surface_average_profile = p_surface_average.vector();
149 std::transform(surface_average_profile.begin(), surface_average_profile.end(),
150 surface_average_profile.begin(),
151 [inverse_surface_are](T value) { return value * inverse_surface_are; });
152}
153
154/**
155 * @brief Measrue surface tension through Polyakov loop fluctuations
156 *
157 * @tparam group gauge field group
158 * @param U gauge field
159 * @param p parameter list
160 * @param traj current update step
161 */
162template <typename group>
163void measure_polyakov_surface(GaugeField<group> &U, const parameters &p, int traj) {
164 // Keep track of printing label
165 static bool first = true;
166 using MyType = Complex<float>;
167 Field<MyType> polyakov_field;
168 if (first) {
169 std::vector<MyType> dummy(lattice.size(e_z));
170 for (int sl = 0; sl < p.n_smear.size(); sl++) {
171 print_formatted_numbers(dummy, "PRO " + std::to_string(sl) + " Origin polyakov smeared",
172 true);
173 print_formatted_numbers(
174 dummy, "PRO " + std::to_string(sl) + " Average polyakov smeared", true);
175 first = false;
176 }
177 }
178 // 1. First we measure Polyakov field
179 if (0) {
180 // this section does local sums of poly lines
181 Field<group> staples, Ut;
182
183 Ut = U[e_t];
184
185 for (int s = 0; s < 20; s++) {
186 staplesum(U, staples, e_t);
187 U[e_t][ALL] = (U[e_t][X] + 0.5 * staples[X]) / (1 + 6 * 0.5);
188 }
189
190 measure_polyakov_field(U[e_t], polyakov_field);
191
192 U[e_t] = Ut;
193
194 } else {
195 // here standard non-link integrated field
196 measure_polyakov_field(U[e_t], polyakov_field);
197 }
198
199 hila::out0 << std::setprecision(5);
200
201 std::vector<MyType> surface_origin_profile, surface_average_profile;
202
203 // 2. Smearing subprocess
204 int prev_smear = 0;
205 for (int sl = 0; sl < p.n_smear.size(); sl++) {
206
207 // 2.1 Smear polyakov field
208 int smear = p.n_smear.at(sl);
209 smear_polyakov_field(polyakov_field, smear - prev_smear, p.smear_coeff, p.twist_coeff);
210 prev_smear = smear;
211
212 Field<MyType> polyakov_field_z = polyakov_field;
213 if (p.z_smear.at(sl) > 0) {
214 Field<MyType> sub_polyakov_field;
215 for (int j = 0; j < p.z_smear.at(sl); j++) {
216 double twist = p.twist_coeff;
217 twist /= NCOLOR;
218 // hila::out0 << twist << '\n';
219 onsites(ALL) if (X.coordinate(e_t) == 0) {
220
221 if (X.z() == 0)
222 sub_polyakov_field[X] =
223 polyakov_field_z[X] +
224 p.smear_coeff * (expi(-2 * M_PI * twist) * polyakov_field_z[X + e_z] +
225 polyakov_field_z[X - e_z]);
226 else if (X.z() == 1)
227 sub_polyakov_field[X] =
228 polyakov_field_z[X] +
229 p.smear_coeff * (polyakov_field_z[X + e_z] +
230 expi(2 * M_PI * twist) * polyakov_field_z[X - e_z]);
231 else
232 sub_polyakov_field[X] =
233 polyakov_field_z[X] +
234 p.smear_coeff * (polyakov_field_z[X + e_z] + polyakov_field_z[X - e_z]);
235 }
236 onsites(ALL) if (X.coordinate(e_t) == 0) {
237 polyakov_field_z[X] = sub_polyakov_field[X] / (1 + 2 * p.smear_coeff);
238 }
239 }
240 }
241
242 measure_polyakov_profile(polyakov_field_z, surface_origin_profile, surface_average_profile);
243
244 // double inverse_surface_are = 1.0 / (lattice.size(e_x) *
245 // lattice.size(e_y)); for (int i = 0; i < surface_average_profile.size();
246 // i++) {
247 // surface_average_profile[i] *= inverse_surface_are;
248 // hila::out0 << "PRO " << sl << " " << i << " (" <<
249 // surface_average_profile[i].re <<
250 // ","
251 // << surface_average_profile[i].im << ") (" <<
252 // surface_origin_profile[i].re
253 // << "," << surface_origin_profile[i].im << ")\n";
254 // }
255 // for (int i = 0; i < surface_average_profile.size(); i++) {
256 // surface_average_profile[i] *= inverse_surface_are;
257 // hila::out0 << "Polyakov_smeared: " << surface_average_profile[i].re
258 // << " " << surface_average_profile[i].im << ", ";
259 // }
260 // hila::out0 << "PRO " << sl << "\n";
261
262 // print_formatted_numbers(surface_origin_profile,"PRO " +
263 // std::to_string(sl) + " Origin polyakov smeared", false);
264 // print_formatted_numbers(surface_average_profile,"PRO " +
265 // std::to_string(sl) + " Average polyakov smeared", false);
266 float min = 1e8;
267 int minloc_global;
268 for (int i = 0; i < surface_average_profile.size(); i++) {
269 if (min > surface_average_profile[i].abs()) {
270 min = surface_average_profile[i].abs();
271 minloc_global = i;
272 }
273 }
274 // hila::out0 << "PRO " << sl << " Min: " <<
275 // surface_average_profile[minloc].abs() << " " << minloc <<
276 // "\nPRO " << sl << " Max: " <<
277 // surface_average_profile[maxloc].abs()<< " " << maxloc <<
278 // std::endl;
279 hila::synchronize_threads();
280 // find the surface between minloc and maxloc
281 // float surface_level = max * 0.5; // assume min is really 0
282 int area = lattice.size(e_x) * lattice.size(e_y);
283
284 // hila::out0 << "Surface_level" << sl << ' ' << surface_level << '\n';
285
286 // int startloc;
287 // if (maxloc > minloc)
288 // startloc = (maxloc + minloc) / 2;
289 // else
290 // startloc =
291 // ((maxloc + minloc + lattice.size(e_z)) / 2) % lattice.size(e_z);
292
293 // starting positio for the other surface
294 // startloc2 = z_ind(startloc + lattice.size(e_z) / 2);
295
296 // hila::out0 << "Start location: " << startloc << std::endl;
297
298 std::vector<float> surf_interpolated; //, surf_discrete;
299 // Only allocate on first rank
300 if (hila::myrank() == 0) {
301 surf_interpolated.resize(area);
302 // surf_discrete.resize(area);
303 }
304
305 hila::out0 << std::setprecision(6);
306
307 std::vector<MyType> polyakov_3D_volume;
308 std::vector<MyType> line(lattice.size(e_z));
309
310 // get full xyz-volume t=0 slice to main node
311 // polyakov_3D_volume = polyakov_field_z.get_slice({-1, -1, -1, 0});
312 // hila::out0 << traj << " " << p.n_trajectories << std::endl;
313
314 for (int y = 0; y < lattice.size(e_y); y++) {
315 // get now full xz-plane polyakov line to main node
316 // reduces MPI calls compared with doing line-by-line
317 polyakov_3D_volume = polyakov_field_z.get_slice({-1, y, -1, 0});
318 if (hila::myrank() == 0) {
319 for (int x = 0; x < lattice.size(e_x); x++) {
320 // line = polyakov_field_z.get_slice({x, y, -1, 0});
321
322 // copy ploop data to line - x runs fastest
323 for (int z = 0; z < lattice.size(e_z); z++) {
324 line[z] = polyakov_3D_volume[x + lattice.size(e_x) * (z)];
325 }
326
327 // start search of the surface from the center between min and max
328 // int z = startloc;
329
330 // while (line[z_ind(z)].abs() > surface_level && startloc - z <
331 // lattice.size(e_z) * 0.4)
332 // z--;
333
334 // while (line[z_ind(z + 1)].abs() <= surface_level &&
335 // z - startloc < lattice.size(e_z) * 0.4)
336 // z++;
337 min = 1e8;
338 int minloc;
339 for (int i = 0; i < line.size(); i++) {
340 if (min > line[i].abs()) {
341 min = line[i].abs();
342 minloc = i;
343 }
344 }
345 int z = minloc;
346 auto x_1 = line[z_ind(minloc - 1)].abs();
347 auto x_2 = line[z_ind(minloc)].abs();
348 auto x_3 = line[z_ind(minloc + 1)].abs();
349 // hila::out0 << "Check loc: " << x_1 << " " << x_2 << " " << x_3 <<
350 // "\n"; double interpolated_min = minloc - (1.0/2.0)*((x_2 -
351 // x_3)-(x_2 - x_1))/(-1*(x_2 - x_3)-(x_2 - x_1));
352 double interpolated_min =
353 minloc - (x_3 - x_1) / ((2.0 * (x_3 + x_1 - 2.0 * x_2)));
354
355 if (interpolated_min > minloc_global + 0.5 * lattice.size(e_z)) {
356 interpolated_min -= lattice.size(e_z);
357 } else if (interpolated_min < minloc_global - 0.5 * lattice.size(e_z)) {
358 interpolated_min += lattice.size(e_z);
359 }
360 // hila::out0 << line[z_ind(minloc-1)].abs() << " " <<
361 // line[z_ind(minloc)].abs()
362 // << " " <<line[z_ind(minloc+1)].abs() << " " <<
363 // interpolated_min
364 // <<" " << minloc<< std::endl; do linear interpolation
365 // surf_discrete[x + y * lattice.size(e_x)] = z;
366 // surf_interpolated[x + y * lattice.size(e_x)] =
367 // z +
368 // (surface_level - line[z_ind(z)].abs()) / (line[z_ind(z +
369 // 1)].abs() - line[z_ind(z)].abs());
370 surf_interpolated[x + y * lattice.size(e_x)] = interpolated_min;
371
372 // and locate the other surface - start from Lz/2 offset
373
374 // z = startloc2;
375
376 // while (line[z_ind(z)] <= surface_level &&
377 // startloc2 - z < lattice.size(e_z) * 0.4)
378 // z--;
379
380 // while (line[z_ind(z + 1)] > surface_level &&
381 // z - startloc2 < lattice.size(e_z) * 0.4)
382 // z++;
383
384 // //do linear interpolation
385 // surf[x + y * lattice.size(e_x)] = z;
386 // surf2[x + y * lattice.size(e_x)] =
387 // z +
388 // (surface_level - line[z_ind(z)]) / (line[z_ind(z + 1)] -
389 // line[z_ind(z)]);
390 }
391 }
392 }
393
394 if (hila::myrank() == 0) {
395 constexpr int pow_size = 80;
396 std::vector<double> npow(pow_size);
397 std::vector<int> hits(pow_size);
398 // wrap_surface(surf_interpolated);
399 spectraldensity_surface(surf_interpolated, npow, hits);
400 // spectraldensity_surface(surf2, npow, hits);
401
402 if (traj == p.n_trajectories - 1) {
403 write_fourier(npow, hits, pow_size, "fourier_profile_" + std::to_string(smear),
404 APPEND_FILE::TRUE, CLOSE_FILE::TRUE);
405 write_surface(surf_interpolated, "surface_smooth_" + std::to_string(smear),
406 APPEND_FILE::TRUE, CLOSE_FILE::TRUE);
407 } else {
408 write_fourier(npow, hits, pow_size, "fourier_profile_" + std::to_string(smear),
409 APPEND_FILE::TRUE, CLOSE_FILE::FALSE);
410 write_surface(surf_interpolated, "surface_smooth_" + std::to_string(smear),
411 APPEND_FILE::TRUE, CLOSE_FILE::FALSE);
412 }
413 }
414 // if (traj == p.n_trajectories - 1 && sl == 1) {
415 // write_surface(surf_discrete,"surface_discrete");
416 // write_surface(surf_interpolated,"surface_smooth");
417 // }
418 }
419}
420
421#endif
Complex definition.
Definition cmplx.h:56
T squarenorm() const
Compute square norm of Complex number.
Definition cmplx.h:248
The field class implements the standard methods for accessing Fields. Hilapp replaces the parity acce...
Definition field.h:61
std::vector< T > get_slice(const CoordinateVector &c, bool broadcast=false) const
and get a slice (subvolume)
Definition field_comm.h:703
Gauge field class.
Definition gaugefield.h:22
Complex< T > expi(T a)
Definition cmplx.h:1413
T abs(const Complex< T > &a)
Return absolute value of Complex number.
Definition cmplx.h:1187
constexpr Parity ALL
bit pattern: 011
int myrank()
rank of this node
Definition com_mpi.cpp:234
std::ostream out0
This writes output only from main process (node 0)
void staplesum(const GaugeField< T > &U, Field< T > &staples, Direction d1, Parity par=ALL)
Sum the staples of link matrices to direction dir.
Definition staples.h:28
int z_ind(int z)
Helper function to get valid z-coordinate index.
Definition suN_gauge.cpp:29