1#ifndef POLYAKOV_SURFACE_H_
2#define POLYAKOV_SURFACE_H_
8void spectraldensity_surface(std::vector<float> &surf, std::vector<double> &npow,
9 std::vector<int> &hits) {
12 static bool first =
true;
14 static fftw_plan fftwplan;
16 int area = lattice.size(e_x) * lattice.size(e_y);
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);
29 for (
int i = 0; i < area; i++) {
33 fftw_execute(fftwplan);
35 int pow_size = npow.size();
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);
43 int k = x * x + y * y;
45 npow[k] += buf[i].
squarenorm() / (area * area);
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);
74 for (
int plane = lattice.size(e_t) - 2; plane >= 0; plane--) {
82#pragma hila safe_access(polyakov)
84 if (X.coordinate(e_t) == plane) {
85 polyakov[X] = Ut[X] * polyakov[X + e_t];
90 onsites(
ALL)
if (X.coordinate(e_t) == 0) {
91 polyakov_field[X] = trace(polyakov[X]);
97void smear_polyakov_field(
Field<T> &polyakov_field,
int nsmear,
float smear_coeff,
double twist) {
101 for (
int i = 0; i < nsmear; i++) {
102 onsites(
ALL)
if (X.coordinate(e_t) == 0) {
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)
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)
113 smear_coeff * (polyakov_field[X + e_z] +
114 expi(2 * M_PI * twist / NCOLOR) * polyakov_field[X - e_z]);
116 pl2[X] += smear_coeff * (polyakov_field[X + e_z] + polyakov_field[X - e_z]);
119 onsites(
ALL)
if (X.coordinate(e_t) == 0) {
120 polyakov_field[X] = pl2[X] / (1 + 6 * smear_coeff);
136void measure_polyakov_profile(
Field<T> &polyakov_field, std::vector<T> &surface_origin_profile,
137 std::vector<T> &surface_average_profile) {
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];
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; });
162template <
typename group>
163void measure_polyakov_surface(
GaugeField<group> &U,
const parameters &p,
int traj) {
165 static bool first =
true;
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",
173 print_formatted_numbers(
174 dummy,
"PRO " + std::to_string(sl) +
" Average polyakov smeared",
true);
185 for (
int s = 0; s < 20; s++) {
187 U[e_t][
ALL] = (U[e_t][X] + 0.5 * staples[X]) / (1 + 6 * 0.5);
190 measure_polyakov_field(U[e_t], polyakov_field);
196 measure_polyakov_field(U[e_t], polyakov_field);
201 std::vector<MyType> surface_origin_profile, surface_average_profile;
205 for (
int sl = 0; sl < p.n_smear.size(); sl++) {
208 int smear = p.n_smear.at(sl);
209 smear_polyakov_field(polyakov_field, smear - prev_smear, p.smear_coeff, p.twist_coeff);
213 if (p.z_smear.at(sl) > 0) {
215 for (
int j = 0; j < p.z_smear.at(sl); j++) {
216 double twist = p.twist_coeff;
219 onsites(
ALL)
if (X.coordinate(e_t) == 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]);
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]);
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]);
236 onsites(
ALL)
if (X.coordinate(e_t) == 0) {
237 polyakov_field_z[X] = sub_polyakov_field[X] / (1 + 2 * p.smear_coeff);
242 measure_polyakov_profile(polyakov_field_z, surface_origin_profile, surface_average_profile);
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();
279 hila::synchronize_threads();
282 int area = lattice.size(e_x) * lattice.size(e_y);
298 std::vector<float> surf_interpolated;
301 surf_interpolated.resize(area);
307 std::vector<MyType> polyakov_3D_volume;
308 std::vector<MyType> line(lattice.size(e_z));
314 for (
int y = 0; y < lattice.size(e_y); y++) {
317 polyakov_3D_volume = polyakov_field_z.
get_slice({-1, y, -1, 0});
319 for (
int x = 0; x < lattice.size(e_x); x++) {
323 for (
int z = 0; z < lattice.size(e_z); z++) {
324 line[z] = polyakov_3D_volume[x + lattice.size(e_x) * (z)];
339 for (
int i = 0; i < line.size(); i++) {
340 if (min > line[i].
abs()) {
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();
352 double interpolated_min =
353 minloc - (x_3 - x_1) / ((2.0 * (x_3 + x_1 - 2.0 * x_2)));
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);
370 surf_interpolated[x + y * lattice.size(e_x)] = interpolated_min;
395 constexpr int pow_size = 80;
396 std::vector<double> npow(pow_size);
397 std::vector<int> hits(pow_size);
399 spectraldensity_surface(surf_interpolated, npow, hits);
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);
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);
T squarenorm() const
Compute square norm of Complex number.
The field class implements the standard methods for accessing Fields. Hilapp replaces the parity acce...
std::vector< T > get_slice(const CoordinateVector &c, bool broadcast=false) const
and get a slice (subvolume)
T abs(const Complex< T > &a)
Return absolute value of Complex number.
constexpr Parity ALL
bit pattern: 011
int myrank()
rank of this node
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.
int z_ind(int z)
Helper function to get valid z-coordinate index.