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