HILA
Loading...
Searching...
No Matches
suN_gauge_csf.cpp
1#include "hila.h"
2#include "gauge/polyakov.h"
6#include "gauge/staples.h"
7#include "gauge/cs_flow.h"
11
12#ifndef NCOLOR
13#define NCOLOR 3
14#endif
15
16using ftype = float;
18
19// define a struct to hold the input parameters: this
20// makes it simpler to pass the values around
21struct parameters {
22 ftype beta; // inverse gauge coupling
23 int n_traj; // number of trajectories to generate
24 int n_therm; // number of thermalization trajectories (counts only accepted traj.)
25 int n_update; // number of heat-bath sweeps per "trajectory"
26 int n_overrelax; // number of overrelaxation sweeps between heat-batch sweeps
27 int gflow_freq; // number of trajectories between gflow measurements
28 ftype gflow_max_l; // flow scale at which gradient flow stops
29 ftype gflow_l_step; // flow scale interval between flow measurements
30 ftype gflow_a_accu; // desired absolute accuracy of gradient flow integration steps
31 ftype gflow_r_accu; // desired relative accuracy of gradient flow integration steps
32 int num_csflows; // number of repetitions of gradient flow with different cs-flow-evolved configs.
33 ftype csflow_t_step; // flow scale at which gradient flow stops
34 ftype csflow_a_accu; // desired absolute accuracy of gradient flow integration steps
35 ftype csflow_r_accu; // desired relative accuracy of gradient flow integration steps
36 int n_save; // number of trajectories between config. check point
37 std::string config_file;
38 ftype time_offset;
39};
40
41
42///////////////////////////////////////////////////////////////////////////////////
43// heat-bath functions
44
45/**
46 * @brief Wrapper update function
47 * @details Gauge Field update sweep with randomly chosen parities and directions
48 *
49 * @tparam group
50 * @param U GaugeField to update
51 * @param p Parameter struct
52 * @param relax If true evolves GaugeField with over relaxation if false then with heat bath
53 */
54template <typename group>
55void update(GaugeField<group> &U, const parameters &p, bool relax) {
56 for (int i = 0; i < 2 * NDIM; ++i) {
57 int tdp = hila::broadcast((int)(hila::random() * 2 * NDIM));
58 int tdir = tdp / 2;
59 int tpar = 1 + (tdp % 2);
60 // hila::out0 << " " << Parity(tpar) << " -- " << Direction(tdir);
61 update_parity_dir(U, p, Parity(tpar), Direction(tdir), relax);
62 }
63 // hila::out0 << "\n";
64}
65
66/**
67 * @brief Wrapper function to updated GaugeField per direction
68 * @details Computes first staplesum, then uses computed result to evolve GaugeField either with
69 * over relaxation or heat bath
70 *
71 * @tparam group
72 * @param U GaugeField to evolve
73 * @param p parameter struct
74 * @param par Parity
75 * @param d Direction to evolve
76 * @param relax If true evolves GaugeField with over relaxation if false then with heat bath
77 */
78template <typename group>
79void update_parity_dir(GaugeField<group> &U, const parameters &p, Parity par, Direction d,
80 bool relax) {
81
82 static hila::timer hb_timer("Heatbath");
83 static hila::timer or_timer("Overrelax");
84 static hila::timer staples_timer("Staplesum");
85
86 Field<group> staples;
87
88 staples_timer.start();
89
90 staplesum(U, staples, d, par);
91 staples_timer.stop();
92
93 if (relax) {
94
95 or_timer.start();
96
97 onsites(par) {
98#ifdef SUN_OVERRELAX_dFJ
99 suN_overrelax_dFJ(U[d][X], staples[X], p.beta);
100#else
101 suN_overrelax(U[d][X], staples[X]);
102#endif
103 }
104 or_timer.stop();
105
106 } else {
107
108 hb_timer.start();
109 onsites(par) {
110 suN_heatbath(U[d][X], staples[X], p.beta);
111 }
112 hb_timer.stop();
113 }
114}
115
116/**
117 * @brief Evolve gauge field
118 * @details Evolution happens by means of heat bath and over relaxation. For each heatbath update
119 * (p.n_update) we update on average also p.n_overrelax times with over relaxation.
120 *
121 * @tparam group
122 * @param U
123 * @param p
124 */
125template <typename group>
126void do_trajectory(GaugeField<group> &U, const parameters &p) {
127
128 for (int n = 0; n < p.n_update; n++) {
129 for (int i = 0; i <= p.n_overrelax; i++) {
130 bool relax = hila::broadcast((int)(hila::random() * (1 + p.n_overrelax)) != 0);
131 update(U, p, relax);
132 // hila::out0 << relax << "\n";
133 }
134 }
136}
137
138// heat-bath functions
139///////////////////////////////////////////////////////////////////////////////////
140// measurement functions
141
142template <typename group>
143void measure_stuff(const GaugeField<group> &U) {
144 // perform measurements on current gauge and momentum pair (U, E) and
145 // print results in formatted form to standard output
146 static bool first = true;
147 if (first) {
148 // print legend for measurement output
149 hila::out0 << "LMEAS: plaq P.real P.imag\n";
150 first = false;
151 }
152 auto plaq = measure_s_wplaq(U) / (lattice.volume() * NDIM * (NDIM - 1) / 2);
153 auto poly = measure_polyakov(U, e_t);
154 hila::out0 << string_format("MEAS % 0.6e % 0.6e % 0.6e", plaq, poly.real(), poly.imag())
155 << '\n';
156}
157
158// end measurement functions
159///////////////////////////////////////////////////////////////////////////////////
160// load/save config functions
161
162template <typename group>
163void checkpoint(const GaugeField<group> &U, int trajectory, const parameters &p) {
164 double t = hila::gettime();
165 // name of config with extra suffix
166 std::string config_file =
167 p.config_file + "_" + std::to_string(abs((trajectory + 1) / p.n_save) % 2);
168 // save config
169 U.config_write(config_file);
170 // write run_status file
171 if (hila::myrank() == 0) {
172 std::ofstream outf;
173 outf.open("run_status", std::ios::out | std::ios::trunc);
174 outf << "trajectory " << trajectory + 1 << '\n';
175 outf << "seed " << static_cast<uint64_t>(hila::random() * (1UL << 61)) << '\n';
176 outf << "time " << hila::gettime() << '\n';
177 // write config name to status file:
178 outf << "config name " << config_file << '\n';
179 outf.close();
180 }
181 std::stringstream msg;
182 msg << "Checkpointing, time " << hila::gettime() - t;
183 hila::timestamp(msg.str().c_str());
184}
185
186template <typename group>
187bool restore_checkpoint(GaugeField<group> &U, int &trajectory, parameters &p) {
188 uint64_t seed;
189 bool ok = true;
190 p.time_offset = 0;
191 hila::input status;
192 if (status.open("run_status", false, false)) {
193 hila::out0 << "RESTORING FROM CHECKPOINT:\n";
194 trajectory = status.get("trajectory");
195 seed = status.get("seed");
196 p.time_offset = status.get("time");
197 // get config name with suffix from status file:
198 std::string config_file = status.get("config name");
199 status.close();
200 hila::seed_random(seed);
201 U.config_read(config_file);
202 ok = true;
203 } else {
204 std::ifstream in;
205 in.open(p.config_file, std::ios::in | std::ios::binary);
206 if (in.is_open()) {
207 in.close();
208 hila::out0 << "READING initial config\n";
209 U.config_read(p.config_file);
210 ok = true;
211 } else {
212 ok = false;
213 }
214 }
215 return ok;
216}
217
218// end load/save config functions
219///////////////////////////////////////////////////////////////////////////////////
220
221
222int main(int argc, char **argv) {
223
224 // hila::initialize should be called as early as possible
225 hila::initialize(argc, argv);
226
227 // hila provides an input class hila::input, which is
228 // a convenient way to read in parameters from input files.
229 // parameters are presented as key - value pairs, as an example
230 // " lattice size 64, 64, 64, 64"
231 // is read below.
232 //
233 // Values are broadcast to all MPI nodes.
234 //
235 // .get() -method can read many different input types,
236 // see file "input.h" for documentation
237
238 hila::out0 << "SU(" << mygroup::size() << ") heat-bath with overrelaxation\n";
239
240 hila::out0 << "Using floating point epsilon: " << fp<ftype>::epsilon << "\n";
241
242 parameters p;
243
244 hila::input par("parameters");
245
246 CoordinateVector lsize;
247 // reads NDIM numbers
248 lsize = par.get("lattice size");
249 // inverse gauge coupling
250 p.beta = par.get("beta");
251 // number of trajectories
252 p.n_traj = par.get("number of trajectories");
253 // number of heat-bath (HB) sweeps per trajectory
254 p.n_update = par.get("updates in trajectory");
255 // number of overrelaxation sweeps petween HB sweeps
256 p.n_overrelax = par.get("overrelax steps");
257 // number of thermalization trajectories
258 p.n_therm = par.get("thermalization trajs");
259 // gradient flow frequency (number of traj. between w. flow measurement)
260 p.gflow_freq = par.get("gflow freq");
261 // gradient flow max. flow-distance
262 p.gflow_max_l = par.get("gflow max lambda");
263 // gradient flow flow-distance step size
264 p.gflow_l_step = par.get("gflow lambda step");
265 // gradient flow absolute accuracy (per integration step)
266 p.gflow_a_accu = par.get("gflow abs. accuracy");
267 // gradient flow relative accuracy (per integration step)
268 p.gflow_r_accu = par.get("gflow rel. accuracy");
269 // cs flow frequency (number of traj. between w. flow measurement)
270 p.num_csflows = par.get("number of csflows");
271 // cs flow flow-distance step size
272 p.csflow_t_step = par.get("csflow t step");
273 // cs flow absolute accuracy (per integration step)
274 p.csflow_a_accu = par.get("csflow abs. accuracy");
275 // cs flow relative accuracy (per integration step)
276 p.csflow_r_accu = par.get("csflow rel. accuracy");
277 // random seed = 0 -> get seed from time
278 long seed = par.get("random seed");
279 // save config and checkpoint
280 p.n_save = par.get("trajs/saved");
281 // measure surface properties and print "profile"
282 p.config_file = par.get("config name");
283
284 par.close(); // file is closed also when par goes out of scope
285
286 // set up the lattice
287 lattice.setup(lsize);
288
289 // We need random number here
290 hila::seed_random(seed);
291
292 // Alloc gauge field and momenta (E)
295
296 // use negative trajectory for thermal
297 int start_traj = -p.n_therm;
298
299 if (!restore_checkpoint(U, start_traj, p)) {
300 U = 1;
301 }
302
303
304 hila::timer update_timer("Updates");
305 hila::timer measure_timer("Measurements");
306 hila::timer gf_timer("Gradient Flow");
307 hila::timer csf_timer("CS Flow");
308
309 ftype t_step0 = 0;
310 for (int trajectory = start_traj; trajectory < p.n_traj; ++trajectory) {
311
312 ftype ttime = hila::gettime();
313
314 update_timer.start();
315
316 do_trajectory(U, p);
317
318 // put sync here in order to get approx gpu timing
319 hila::synchronize_threads();
320 update_timer.stop();
321
322 measure_timer.start();
323
324 hila::out0 << "Measure_start " << trajectory << '\n';
325
326 measure_stuff(U);
327
328 hila::out0 << "Measure_end " << trajectory << '\n';
329
330 measure_timer.stop();
331
332 if (trajectory >= 0) {
333
334 if (p.gflow_freq > 0 && trajectory % p.gflow_freq == 0) {
335
336
337 int gtrajectory = trajectory / p.gflow_freq;
338 if (p.gflow_l_step > 0) {
339
340
341 for (int ics = 0; ics < p.num_csflows; ++ics) {
342
343 foralldir(d) onsites(ALL) E[d][X].gaussian_random();
344
345 csf_timer.start();
346
347 do_cs_flow_adapt(U, E, p.csflow_t_step, p.csflow_a_accu, p.csflow_r_accu);
348
349 csf_timer.stop();
350
351
352
353 gf_timer.start();
354
355 int nflow_steps = (int)(p.gflow_max_l / p.gflow_l_step);
356
357 ftype gftime = hila::gettime();
358 hila::out0 << "Gflow_start " << gtrajectory << '\n';
359
361
362 ftype t_step = t_step0;
363 measure_gradient_flow_stuff(V, (ftype)0.0, t_step);
364 t_step = do_gradient_flow_adapt(V, (ftype)0.0, p.gflow_l_step, p.gflow_a_accu,
365 p.gflow_r_accu, t_step);
366 measure_gradient_flow_stuff(V, p.gflow_l_step, t_step);
367 t_step0 = t_step;
368
369 for (int i = 1; i < nflow_steps; ++i) {
370
371 t_step =
372 do_gradient_flow_adapt(V, i * p.gflow_l_step, (i + 1) * p.gflow_l_step,
373 p.gflow_a_accu, p.gflow_r_accu, t_step);
374
375 measure_gradient_flow_stuff(V, (i + 1) * p.gflow_l_step, t_step);
376 }
377
378 hila::out0 << "Gflow_end " << gtrajectory << " time " << std::setprecision(3)
379 << hila::gettime() - gftime << '\n';
380
381 gf_timer.stop();
382 }
383
384 }
385 }
386 }
387
388 if (p.n_save > 0 && (trajectory + 1) % p.n_save == 0) {
389 checkpoint(U, trajectory, p);
390 }
391 }
392
394}
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
void config_write(const std::string &filename) const
config_write writes the gauge field to file, with additional "verifying" header
Definition gaugefield.h:155
static constexpr int size()
Returns size of Vector or square Matrix.
Definition matrix.h:248
Class for SU(N) matrix.
Definition sun_matrix.h:110
hila::input - Class for parsing runtime parameter files.
Definition input.h:52
void close()
Closes input parameter file.
Definition input.cpp:79
returntype get(const std::string &key)
Get next value in stack of read in input string from parameters file.
Definition input.h:269
bool open(const std::string &fname, bool use_cmdline=true, bool exit_on_error=true)
Open file that parameters are read from.
Definition input.cpp:22
void setup(const CoordinateVector &siz)
General lattice setup.
Definition lattice.cpp:33
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
double random()
Real valued uniform random number generator.
Definition hila_gpu.cpp:121
int myrank()
rank of this node
Definition com_mpi.cpp:234
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
T broadcast(T &var, int rank=0)
Broadcast the value of var to all MPI ranks from rank (default=0).
Definition com_mpi.h:152
double gettime()
Definition timing.cpp:163
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:73
void update_parity_dir(GaugeField< group > &U, const parameters &p, Parity par, Direction d, bool relax)
Wrapper function to updated GaugeField per direction.
void do_trajectory(GaugeField< group > &U, const parameters &p)
Evolve gauge field.
void measure_stuff(const GaugeField< group > &U, const parameters &p)
Measures Polyakov lines and Wilson action.
Definition suN_gauge.cpp:41
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