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 U, typename T>
136void measure_polyakov_profile(Field<U> &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()] += abs(polyakov_field[X]);
143 if (X.x() == 0 && X.y() == 0)
144 p_origin[X.z()] += abs(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 using MyType = Complex<float>;
166 Field<MyType> polyakov_field;
167 // 1. First we measure Polyakov field
168 if (0) {
169 // this section does local sums of poly lines
170 Field<group> staples, Ut;
171
172 Ut = U[e_t];
173
174 for (int s = 0; s < 20; s++) {
175 staplesum(U, staples, e_t);
176 U[e_t][ALL] = (U[e_t][X] + 0.5 * staples[X]) / (1 + 6 * 0.5);
177 }
178
179 measure_polyakov_field(U[e_t], polyakov_field);
180
181 U[e_t] = Ut;
182
183 } else {
184 // here standard non-link integrated field
185 measure_polyakov_field(U[e_t], polyakov_field);
186 }
187
188 hila::out0 << std::setprecision(5);
189
190 std::vector<double> surface_origin_profile, surface_average_profile;
191
192 // 2. Smearing subprocess
193 int prev_smear = 0;
194 for (int sl = 0; sl < p.n_smear.size(); sl++) {
195
196 // 2.1 Smear polyakov field
197 int smear = p.n_smear.at(sl);
198 smear_polyakov_field(polyakov_field, smear - prev_smear, p.smear_coeff, p.twist_coeff);
199 prev_smear = smear;
200
201 Field<MyType> polyakov_field_z = polyakov_field;
202 if (p.z_smear.at(sl) > 0) {
203 Field<MyType> sub_polyakov_field;
204 for (int j = 0; j < p.z_smear.at(sl); j++) {
205 double twist = p.twist_coeff;
206 twist /= NCOLOR;
207 // hila::out0 << twist << '\n';
208 onsites(ALL) if (X.coordinate(e_t) == 0) {
209
210 if (X.z() == 0)
211 sub_polyakov_field[X] =
212 polyakov_field_z[X] +
213 p.smear_coeff * (expi(-2 * M_PI * twist) * polyakov_field_z[X + e_z] +
214 polyakov_field_z[X - e_z]);
215 else if (X.z() == 1)
216 sub_polyakov_field[X] =
217 polyakov_field_z[X] +
218 p.smear_coeff * (polyakov_field_z[X + e_z] +
219 expi(2 * M_PI * twist) * polyakov_field_z[X - e_z]);
220 else
221 sub_polyakov_field[X] =
222 polyakov_field_z[X] +
223 p.smear_coeff * (polyakov_field_z[X + e_z] + polyakov_field_z[X - e_z]);
224 }
225 onsites(ALL) if (X.coordinate(e_t) == 0) {
226 polyakov_field_z[X] = sub_polyakov_field[X] / (1 + 2 * p.smear_coeff);
227 }
228 }
229 }
230
231 measure_polyakov_profile<MyType,double>(polyakov_field_z, surface_origin_profile, surface_average_profile);
232
233 // double inverse_surface_are = 1.0 / (lattice.size(e_x) *
234 // lattice.size(e_y)); for (int i = 0; i < surface_average_profile.size();
235 // i++) {
236 // surface_average_profile[i] *= inverse_surface_are;
237 // hila::out0 << "PRO " << sl << " " << i << " (" <<
238 // surface_average_profile[i].re <<
239 // ","
240 // << surface_average_profile[i].im << ") (" <<
241 // surface_origin_profile[i].re
242 // << "," << surface_origin_profile[i].im << ")\n";
243 // }
244 // for (int i = 0; i < surface_average_profile.size(); i++) {
245 // surface_average_profile[i] *= inverse_surface_are;
246 // hila::out0 << "Polyakov_smeared: " << surface_average_profile[i].re
247 // << " " << surface_average_profile[i].im << ", ";
248 // }
249 // hila::out0 << "PRO " << sl << "\n";
250
251 print_formatted_numbers(surface_origin_profile,"PRO " +
252 std::to_string(sl) + " Origin polyakov smeared", false);
253 print_formatted_numbers(surface_average_profile,"PRO " +
254 std::to_string(sl) + " Average polyakov smeared", false);
255 float min = 1e8;
256 int minloc_global;
257 for (int i = 0; i < surface_average_profile.size(); i++) {
258 if (min > surface_average_profile[i]) {
259 min = surface_average_profile[i];
260 minloc_global = i;
261 }
262 }
263 // hila::out0 << "PRO " << sl << " Min: " <<
264 // surface_average_profile[minloc].abs() << " " << minloc <<
265 // "\nPRO " << sl << " Max: " <<
266 // surface_average_profile[maxloc].abs()<< " " << maxloc <<
267 // std::endl;
268 hila::synchronize_threads();
269 // find the surface between minloc and maxloc
270 // float surface_level = max * 0.5; // assume min is really 0
271 int area = lattice.size(e_x) * lattice.size(e_y);
272
273 // hila::out0 << "Surface_level" << sl << ' ' << surface_level << '\n';
274
275 // int startloc;
276 // if (maxloc > minloc)
277 // startloc = (maxloc + minloc) / 2;
278 // else
279 // startloc =
280 // ((maxloc + minloc + lattice.size(e_z)) / 2) % lattice.size(e_z);
281
282 // starting positio for the other surface
283 // startloc2 = z_ind(startloc + lattice.size(e_z) / 2);
284
285 // hila::out0 << "Start location: " << startloc << std::endl;
286
287 std::vector<float> surf_interpolated; //, surf_discrete;
288 // Only allocate on first rank
289 if (hila::myrank() == 0) {
290 surf_interpolated.resize(area);
291 // surf_discrete.resize(area);
292 }
293
294 hila::out0 << std::setprecision(6);
295
296 std::vector<MyType> polyakov_3D_volume;
297 std::vector<MyType> line(lattice.size(e_z));
298
299 // get full xyz-volume t=0 slice to main node
300 // polyakov_3D_volume = polyakov_field_z.get_slice({-1, -1, -1, 0});
301 // hila::out0 << traj << " " << p.n_trajectories << std::endl;
302
303 for (int y = 0; y < lattice.size(e_y); y++) {
304 // get now full xz-plane polyakov line to main node
305 // reduces MPI calls compared with doing line-by-line
306 polyakov_3D_volume = polyakov_field_z.get_slice({-1, y, -1, 0});
307 if (hila::myrank() == 0) {
308 for (int x = 0; x < lattice.size(e_x); x++) {
309 // line = polyakov_field_z.get_slice({x, y, -1, 0});
310
311 // copy ploop data to line - x runs fastest
312 for (int z = 0; z < lattice.size(e_z); z++) {
313 line[z] = polyakov_3D_volume[x + lattice.size(e_x) * (z)];
314 }
315
316 // start search of the surface from the center between min and max
317 // int z = startloc;
318
319 // while (line[z_ind(z)].abs() > surface_level && startloc - z <
320 // lattice.size(e_z) * 0.4)
321 // z--;
322
323 // while (line[z_ind(z + 1)].abs() <= surface_level &&
324 // z - startloc < lattice.size(e_z) * 0.4)
325 // z++;
326 // min = 1e8;
327 // int minloc;
328 // for (int i = 0; i < line.size(); i++) {
329 // if (min > line[i].abs()) {
330 // min = line[i].abs();
331 // minloc = i;
332 // }
333 // }
334 // FIND LOCAL MIN BASED OFF OF GLOBAL MINIMA
335
336 // THINK HERE
337 int minloc;
338 if (line[z_ind(minloc_global)].abs() < line[z_ind(minloc_global+1)].abs() && line[z_ind(minloc_global)].abs() < line[z_ind(minloc_global-1)].abs()) {
339 minloc = minloc_global;
340 } else {
341 for (int i = 1; i < line.size()/2; i++) {
342 if (line[z_ind(minloc_global+i)].abs() < line[z_ind(minloc_global+i+1)].abs() && line[z_ind(minloc_global+i)].abs() < line[z_ind(minloc_global+i-1)].abs()) {
343 minloc = minloc_global+i;
344 break;
345 }
346 if (line[z_ind(minloc_global-i)].abs() < line[z_ind(minloc_global-i+1)].abs() && line[z_ind(minloc_global-i)].abs() < line[z_ind(minloc_global-i-1)].abs()) {
347 minloc = minloc_global-i;
348 break;
349 }
350 }
351 }
352 int z = minloc;
353 double x_1 = line[z_ind(minloc-1)].abs();
354 double x_2 = line[z_ind(minloc)].abs();
355 double x_3 = line[z_ind(minloc+1)].abs();
356 // hila::out0 << "Check loc: " << x_1 << " " << x_2 << " " << x_3 <<
357 // "\n";
358 //double interpolated_min = minloc - (1.0/2.0)*((x_2 - x_3)-(x_2 - x_1))/(-1*(x_2 - x_3)-(x_2 - x_1));
359 double interpolated_min =
360 minloc - (x_3 - x_1) / ((2.0 * (x_3 + x_1 - 2.0 * x_2)));
361
362// if (interpolated_min > minloc_global + 0.5 * lattice.size(e_z)) {
363// interpolated_min -= lattice.size(e_z);
364// } else if (interpolated_min < minloc_global - 0.5 * lattice.size(e_z)) {
365// interpolated_min += lattice.size(e_z);
366// }
367
368 surf_interpolated[x + y * lattice.size(e_x)] = interpolated_min;
369
370 // and locate the other surface - start from Lz/2 offset
371
372 // z = startloc2;
373
374 // while (line[z_ind(z)] <= surface_level &&
375 // startloc2 - z < lattice.size(e_z) * 0.4)
376 // z--;
377
378 // while (line[z_ind(z + 1)] > surface_level &&
379 // z - startloc2 < lattice.size(e_z) * 0.4)
380 // z++;
381
382 // //do linear interpolation
383 // surf[x + y * lattice.size(e_x)] = z;
384 // surf2[x + y * lattice.size(e_x)] =
385 // z +
386 // (surface_level - line[z_ind(z)]) / (line[z_ind(z + 1)] -
387 // line[z_ind(z)]);
388 }
389 }
390 }
391
392 if (hila::myrank() == 0) {
393 constexpr int pow_size = 80;
394 std::vector<double> npow(pow_size);
395 std::vector<int> hits(pow_size);
396 // wrap_surface(surf_interpolated);
397 spectraldensity_surface(surf_interpolated, npow, hits);
398 // spectraldensity_surface(surf2, npow, hits);
399
400 if (traj == p.n_trajectories - 1) {
401 write_fourier(npow, hits, pow_size, p.out_folder + "fourier_profile_" + std::to_string(smear) + "_" + std::to_string(z_smear),
402 APPEND_FILE::TRUE, CLOSE_FILE::TRUE);
403 write_surface(surf_interpolated, p.out_folder + "surface_smooth_" + std::to_string(smear) + "_" + std::to_string(z_smear),
404 APPEND_FILE::TRUE, CLOSE_FILE::TRUE);
405 } else {
406 write_fourier(npow, hits, pow_size, p.out_folder + "fourier_profile_" + std::to_string(smear) + "_" + std::to_string(z_smear),
407 APPEND_FILE::TRUE, CLOSE_FILE::FALSE);
408 write_surface(surf_interpolated, p.out_folder + "surface_smooth_" + std::to_string(smear) + "_" + std::to_string(z_smear),
409 APPEND_FILE::TRUE, CLOSE_FILE::FALSE);
410 }
411 }
412 // if (traj == p.n_trajectories - 1 && sl == 1) {
413 // write_surface(surf_discrete,"surface_discrete");
414 // write_surface(surf_interpolated,"surface_smooth");
415 // }
416 }
417}
418
419#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