Loading [MathJax]/extensions/TeX/AMSsymbols.js
HILA
All Classes Namespaces Files Functions Variables Typedefs Enumerations Friends Macros Pages
twist_specific_methods.h
1#ifndef TWIST_SPECIFIC_METHODS_H_
2#define TWIST_SPECIFIC_METHODS_H_
3
4#include "hila.h"
5
6/**
7 * @brief Sum the staples of link matrices to direction dir by inducing twist at
8 * z=0 t=0 on x-y plane
9 *
10 * @tparam T
11 * @param U GaugeField to compute staples for
12 * @param staples Filed to compute staplesum into at each lattice point
13 * @param d1 Direction to compute staplesum for
14 * @param par Parity to compute staplesum for
15 * @param twist_coeff Integer to rotate phase with
16 */
17template <typename T>
18void staplesum_twist(const GaugeField<T> &U, Field<T> &staples, Direction d1,
19 int twist_coeff, Parity par = ALL) {
20
21 Field<T> lower;
22 // GaugeField<double> twist = 0;
23 // onsites(ALL) {
24 // if (X.z() == 0 && X.t() == 0) {
25 // twist[e_z][X] = twist_coeff;
26 // twist[e_t][X] = -twist_coeff;
27 // }
28 // }
29
30 bool first = true;
31 staples = 0;
32
33 foralldir(d2) if (d2 != d1) {
34
35 // anticipate that these are needed
36 // not really necessary, but may be faster
37 U[d2].start_gather(d1, ALL);
38 U[d1].start_gather(d2, par);
39
40 double twist;
41 if (d1 == e_z && d2 == e_t)
42 twist = twist_coeff;
43 else if (d1 == e_t && d2 == e_z)
44 twist = -twist_coeff;
45 else
46 twist = 0;
47
48 twist /= NCOLOR;
49
50 // calculate first lower 'U' of the staple sum
51 // do it on opp parity
52 onsites(opp_parity(par)) {
53 lower[X] = U[d2][X].dagger() * U[d1][X] * U[d2][X + d1];
54 if (X.z() == 0 && X.t() == 0)
55 lower[X] *= expi(-2 * M_PI * twist);
56 }
57
58 // calculate then the upper 'n', and add the lower
59 // lower could also be added on a separate loop
60 onsites(par) {
61 auto upper = U[d2][X] * U[d1][X + d2] * U[d2][X + d1].dagger();
62 if (X.z() == 0 && X.t() == 0)
63 upper *= expi(2 * M_PI * twist);
64
65 staples[X] += upper + lower[X - d2];
66 }
67 }
68}
69
70/**
71 * @brief Computes Wilson action
72 * @details \f{align}{ S &= \beta\sum_{\textbf{dir}_1 < \textbf{dir}_2}\sum_{X}
73 * \frac{1}{N} \Re\mathrm{Tr}\left[ 1- U_{\textbf{dir}_1 \textbf{dir}_2}(X)
74 * \right] \f} Where \f$\beta = 2N/g^2\f$
75 *
76 * @return double
77 */
78template <typename T>
79std::vector<double> measure_plaq_with_z(GaugeField<T> U, int twist_coeff) {
80 ReductionVector<double> plaq_vec(lattice.size(e_z) + 1);
81 plaq_vec = 0;
82 plaq_vec.allreduce(false);
83 GaugeField<double> twist = 0;
84
85 onsites(ALL) {
86 if (X.z() == 0 && X.t() == 0) {
87 twist[e_z][X] = -twist_coeff;
88 }
89 }
90
91 foralldir(dir1) foralldir(dir2) if (dir1 < dir2) {
92
93 onsites(ALL) {
94 double p;
95 p = 1.0 - real(expi(2 * M_PI * (twist[dir1][X] / NCOLOR)) *
96 trace(U[dir1][X] * U[dir2][X + dir1] *
97 U[dir1][X + dir2].dagger() * U[dir2][X].dagger())) /
98 T::size();
99 plaq_vec[lattice.size(e_z)] += p;
100 plaq_vec[X.z()] += p;
101 }
102 }
103 // plaq_vec[lattice.size(e_z)] =
104 // plaq.value() / (lattice.volume() * NDIM * (NDIM - 1) / 2);
105 for (int i = 0; i < plaq_vec.size(); i++) {
106 plaq_vec[i] /= (lattice.volume() * NDIM * (NDIM - 1)) / 2;
107 }
108
109 return plaq_vec.vector();
110}
111
112
113/**
114 * @brief Measure Polyakov lines to direction dir bining based on z index
115 * @details Naive implementation, includes extra communication
116 * @tparam T GaugeField Group
117 * @param U GaugeField to measure
118 * @param dir Direction
119 * @return Complex<double>
120 */
121template <typename T>
122std::vector<Complex<double>>
123measure_polyakov_with_z(const GaugeField<T> &U,
124 Direction dir = Direction(NDIM - 1)) {
125
126 Field<T> polyakov = U[dir];
127
128 // mult links so that polyakov[X.dir == 0] contains the polyakov loop
129 for (int plane = lattice.size(dir) - 2; plane >= 0; plane--) {
130
131 // safe_access(polyakov) pragma allows the expression below, otherwise
132 // hilapp would reject it because X and X+dir can refer to the same
133 // site on different "iterations" of the loop. However, here this
134 // is restricted on single dir-plane so it works but we must tell it to
135 // hilapp.
136
137#pragma hila safe_access(polyakov)
138 onsites(ALL) {
139 if (X.coordinate(dir) == plane) {
140 polyakov[X] = U[dir][X] * polyakov[X + dir];
141 }
142 }
143 }
144
145 Complex<double> ploop = 0;
146 ReductionVector<Complex<double>> ploop_z(lattice.size(e_z) + 1);
147 ploop_z = 0;
148 onsites(ALL) if (X.coordinate(dir) == 0) {
149 Complex<double> p = trace(polyakov[X]);
150 ploop_z[lattice.size(e_z)] += p;
151 ploop_z[X.z()] += p;
152 }
153 for (int i = 0; i < ploop_z.size(); i++) {
154 ploop_z[i] /= (lattice.volume() / lattice.size(dir));
155 }
156 // return average polyakov
157 return ploop_z.vector();
158 // return ploop;
159}
160
161/**
162 * @brief Measure Polyakov lines to direction dir bining based on z index
163 * @details Naive implementation, includes extra communication
164 * @tparam T GaugeField Group
165 * @param U GaugeField to measure
166 * @param dir Direction
167 * @return Complex<double>
168 */
169template <typename T>
170std::vector<double>
171measure_polyakov_with_z_abs(const GaugeField<T> &U,
172 Direction dir = Direction(NDIM - 1)) {
173
174 Field<T> polyakov = U[dir];
175
176 // mult links so that polyakov[X.dir == 0] contains the polyakov loop
177 for (int plane = lattice.size(dir) - 2; plane >= 0; plane--) {
178
179 // safe_access(polyakov) pragma allows the expression below, otherwise
180 // hilapp would reject it because X and X+dir can refer to the same
181 // site on different "iterations" of the loop. However, here this
182 // is restricted on single dir-plane so it works but we must tell it to
183 // hilapp.
184
185#pragma hila safe_access(polyakov)
186 onsites(ALL) {
187 if (X.coordinate(dir) == plane) {
188 polyakov[X] = U[dir][X] * polyakov[X + dir];
189 }
190 }
191 }
192
193 double ploop = 0;
194 ReductionVector<double> ploop_z(lattice.size(e_z) + 1);
195 ploop_z = 0;
196 onsites(ALL) if (X.coordinate(dir) == 0) {
197 double p = abs(trace(polyakov[X]));
198 ploop_z[lattice.size(e_z)] += p;
199 ploop_z[X.z()] += p;
200 }
201 for (int i = 0; i < ploop_z.size(); i++) {
202 ploop_z[i] /= (lattice.size(e_x) * lattice.size(e_y));
203 }
204 // return average polyakov
205 return ploop_z.vector();
206 // return ploop;
207}
208
209#endif
Array< n, m, hila::arithmetic_type< T > > real(const Array< n, m, T > &arg)
Return real part of Array.
Definition array.h:688
Complex definition.
Definition cmplx.h:56
The field class implements the standard methods for accessing Fields. Hilapp replaces the parity acce...
Definition field.h:61
Gauge field class.
Definition gaugefield.h:22
Complex< T > dagger(const Complex< T > &val)
Return dagger of Complex number.
Definition cmplx.h:1223
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
Parity
Parity enum with values EVEN, ODD, ALL; refers to parity of the site. Parity of site (x,...
#define foralldir(d)
Macro to loop over (all) Direction(s)
Definition coordinates.h:78
Direction
Enumerator for direction that assigns integer to direction to be interpreted as unit vector.
Definition coordinates.h:34
constexpr Parity ALL
bit pattern: 011