HILA
Loading...
Searching...
No Matches
suN_gauge.cpp
Go to the documentation of this file.
1/**
2 * @file suN_gauge.cpp
3 * @brief Application to simulate \f$ SU(N) \f$ Gauge field.
4 * @details Simple application which Generates \f$ SU(N) \f$ GaugeField using \ref staplesum, \ref
5 * suN_overrelax and \ref suN_heatbath. Each evolution, the application measures the Wilson action
6 * using GaugeField::measure_plaq and Polyakov lines using \ref measure_polyakov.
7 *
8 */
9#include "hila.h"
10#include "gauge/staples.h"
11#include "gauge/polyakov.h"
12#include "gauge/stout_smear.h"
13#include "gauge/sun_heatbath.h"
14#include "gauge/sun_overrelax.h"
15
16#include <fftw3.h>
17
18// local includes
19#include "parameters.h"
20#include "checkpoint.h"
21
22/**
23 * @brief Helper function to get valid z-coordinate index
24 *
25 * @param z
26 * @return int
27 */
28int z_ind(int z) {
29 return (z + lattice.size(e_z)) % lattice.size(e_z);
30}
31
32/**
33 * @brief Measures Polyakov lines and Wilson action
34 *
35 * @tparam group
36 * @param U GaugeField to measure
37 * @param p Parameter struct
38 */
39template <typename group>
40void measure_stuff(const GaugeField<group> &U, const parameters &p) {
41
42 static bool first = true;
43
44 if (first) {
45 hila::out0 << "Legend:";
46 hila::out0 << " plaq P.real P.imag\n";
47
48 first = false;
49 }
50
51 auto poly = measure_polyakov(U);
52
53 auto plaq = U.measure_plaq() / (lattice.volume() * NDIM * (NDIM - 1) / 2);
54
55 hila::out0 << "MEAS " << std::setprecision(8);
56
57 // write the -(polyakov potential) first, this is used as a weight factor in aa
58
59 hila::out0 << plaq << ' ' << poly << '\n';
60}
61
62/**
63 * @brief Wrapper update function
64 * @details Updates Gauge Field one direction at a time first EVEN then ODD parity
65 *
66 * @tparam group
67 * @param U GaugeField to update
68 * @param p Parameter struct
69 * @param relax If true evolves GaugeField with over relaxation if false then with heat bath
70 */
71template <typename group>
72void update(GaugeField<group> &U, const parameters &p, bool relax) {
73
74 foralldir(d) {
75 for (Parity par : {EVEN, ODD}) {
76
77 update_parity_dir(U, p, par, d, relax);
78 }
79 }
80}
81
82/**
83 * @brief Wrapper function to updated GaugeField per direction
84 * @details Computes first staplesum, then uses computed result to evolve GaugeField either with
85 * over relaxation or heat bath
86 *
87 * @tparam group
88 * @param U GaugeField to evolve
89 * @param p parameter struct
90 * @param par Parity
91 * @param d Direction to evolve
92 * @param relax If true evolves GaugeField with over relaxation if false then with heat bath
93 */
94template <typename group>
95void update_parity_dir(GaugeField<group> &U, const parameters &p, Parity par, Direction d,
96 bool relax) {
97
98 static hila::timer hb_timer("Heatbath");
99 static hila::timer or_timer("Overrelax");
100 static hila::timer staples_timer("Staplesum");
101
102 Field<group> staples;
103
104 staples_timer.start();
105
106 staplesum(U, staples, d, par);
107 staples_timer.stop();
108
109 if (relax) {
110
111 or_timer.start();
112
113 onsites(par) {
114#ifdef SUN_OVERRELAX_dFJ
115 suN_overrelax_dFJ(U[d][X], staples[X], p.beta);
116#else
117 suN_overrelax(U[d][X], staples[X]);
118#endif
119 }
120 or_timer.stop();
121
122 } else {
123
124 hb_timer.start();
125 onsites(par) {
126 suN_heatbath(U[d][X], staples[X], p.beta);
127 }
128 hb_timer.stop();
129 }
130}
131
132/**
133 * @brief Evolve gauge field
134 * @details Evolution happens by means of heat bath and over relaxation. For each heatbath update
135 * (p.n_update) we update p.n_overrelax times with over relaxation.
136 *
137 * @tparam group
138 * @param U
139 * @param p
140 */
141template <typename group>
142void do_trajectory(GaugeField<group> &U, const parameters &p) {
143
144 for (int n = 0; n < p.n_update; n++) {
145 for (int i = 0; i < p.n_overrelax; i++) {
146 update(U, p, true);
147 }
148 update(U, p, false);
149 }
151}
152
153
154int main(int argc, char **argv) {
155
156 // hila::initialize should be called as early as possible
157 hila::initialize(argc, argv);
158
159 // hila provides an input class hila::input, which is
160 // a convenient way to read in parameters from input files.
161 // parameters are presented as key - value pairs, as an example
162 // " lattice size 64, 64, 64, 64"
163 // is read below.
164 //
165 // Values are broadcast to all MPI nodes.
166 //
167 // .get() -method can read many different input types,
168 // see file "input.h" for documentation
169
170 parameters p;
171
172 hila::out0 << "SU(" << mygroup::size() << ") heat bath + overrelax update\n";
173
174 hila::input par("parameters");
175
176 CoordinateVector lsize;
177 lsize = par.get("lattice size"); // reads NDIM numbers
178
179 p.beta = par.get("beta");
180 // deltab sets system to different beta on different sides, by beta*(1 +- deltab)
181 // use for initial config generation only
182 p.deltab = par.get("delta beta fraction");
183 // trajectory length in steps
184 p.n_overrelax = par.get("overrelax steps");
185 p.n_update = par.get("updates in trajectory");
186 p.n_trajectories = par.get("trajectories");
187 p.n_thermal = par.get("thermalization");
188
189 // random seed = 0 -> get seed from time
190 uint64_t seed = par.get("random seed");
191 // save config and checkpoint
192 p.n_save = par.get("traj/saved");
193 // measure surface properties and print "profile"
194 p.config_file = par.get("config name");
195
196 par.close(); // file is closed also when par goes out of scope
197
198 // setting up the lattice is convenient to do after reading
199 // the parameter
200 lattice.setup(lsize);
201
202 // Alloc gauge field
204
205 U = 1;
206
207 // use negative trajectory for thermal
208 int start_traj = -p.n_thermal;
209
210 hila::timer update_timer("Updates");
211 hila::timer measure_timer("Measurements");
212
213 restore_checkpoint(U, start_traj, p);
214
215 // We need random number here
216 if (!hila::is_rng_seeded())
217 hila::seed_random(seed);
218
219 for (int trajectory = start_traj; trajectory < p.n_trajectories; trajectory++) {
220
221 double ttime = hila::gettime();
222
223 update_timer.start();
224
225 double acc = 0;
226
227 do_trajectory(U, p);
228
229 // put sync here in order to get approx gpu timing
230 hila::synchronize_threads();
231 update_timer.stop();
232
233 // trajectory is negative during thermalization
234 if (trajectory >= 0) {
235 measure_timer.start();
236
237 hila::out0 << "Measure_start " << trajectory << '\n';
238
239 measure_stuff(U, p);
240
241 hila::out0 << "Measure_end " << trajectory << std::endl;
242
243 measure_timer.stop();
244 }
245
246 if (p.n_save > 0 && (trajectory + 1) % p.n_save == 0) {
247 checkpoint(U, trajectory, p);
248 }
249 }
250
252}
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
void reunitarize_gauge()
Reunitarize Gauge Field consisting of matrices.
Definition gaugefield.h:94
double measure_plaq() const
Computes Wilson action.
Definition gaugefield.h:108
static constexpr int size()
Returns size of Vector or square Matrix.
Definition matrix.h:242
hila::input - Class for parsing runtime parameter files.
Definition input.h:52
void setup(const CoordinateVector &siz)
General lattice setup.
Definition lattice.cpp:33
Parity
Parity enum with values EVEN, ODD, ALL; refers to parity of the site. Parity of site (x,...
constexpr Parity EVEN
bit pattern: 001
#define foralldir(d)
Macro to loop over (all) Direction(s)
Definition coordinates.h:78
constexpr Parity ODD
bit pattern: 010
Direction
Enumerator for direction that assigns integer to direction to be interpreted as unit vector.
Definition coordinates.h:34
std::ostream out0
This writes output only from main process (node 0)
void initialize(int argc, char **argv)
Read in command line arguments. Initialise default stream and MPI communication.
void seed_random(uint64_t seed, bool device_rng=true)
Seed random generators with 64-bit unsigned value. On MPI shuffles the seed so that different MPI ran...
Definition random.cpp:86
double gettime()
Definition timing.cpp:163
bool is_rng_seeded()
Check if RNG is seeded already.
Definition random.cpp:239
void finishrun()
Normal, controlled exit - all nodes must call this. Prints timing information and information about c...
Complex< double > measure_polyakov(const GaugeField< T > &U, Direction dir=Direction(NDIM - 1))
Measure Polyakov lines to direction dir.
Definition polyakov.h:17
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
void update(GaugeField< group > &U, const parameters &p, bool relax)
Wrapper update function.
Definition suN_gauge.cpp:72
void update_parity_dir(GaugeField< group > &U, const parameters &p, Parity par, Direction d, bool relax)
Wrapper function to updated GaugeField per direction.
Definition suN_gauge.cpp:95
void do_trajectory(GaugeField< group > &U, const parameters &p)
Evolve gauge field.
int z_ind(int z)
Helper function to get valid z-coordinate index.
Definition suN_gauge.cpp:28
void measure_stuff(const GaugeField< group > &U, const parameters &p)
Measures Polyakov lines and Wilson action.
Definition suN_gauge.cpp:40
double suN_heatbath(SU< N, T > &U, const SU< N, T > &staple, double beta)
heatbath
void suN_overrelax(SU< N, T > &U, const SU< N, T > &staple)
Overrelaxation using SU(2) subgroups