HILA
Loading...
Searching...
No Matches
suN_gauge_twist.cpp
1#include "gauge/staples.h"
2#include "gauge/polyakov.h"
3#include "gauge/stout_smear.h"
6#include "hila.h"
7#include "multicanonical.h"
8
9//#include "gauge/polyakov.h"
10
11#include <fftw3.h>
12#include <numeric>
13#include <stdio.h>
14#include <stdlib.h>
15
16// local includes
17#include "checkpoint.h"
18#include "parameters.h"
19#include "polyakov_surface.h"
20#include "twist_specific_methods.h"
21#include "utility.h"
22// /**
23// * @brief Helper function to get valid z-coordinate index
24// *
25// * @param z
26// * @return int
27// */
28// int z_ind(int z) { return (z + lattice.size(e_z)) % lattice.size(e_z); }
29
30/**
31 * @brief Measures Polyakov lines and Wilson action
32 *
33 * @tparam group
34 * @param U GaugeField to measure
35 * @param p Parameter struct
36 */
37template <typename group>
38void measure_plaq(const GaugeField<group> &U, const parameters &p) {
39
40
41 auto plaq_z = measure_plaq_with_z(
42 U, p.twist_coeff);
43
44 print_formatted_numbers(plaq_z, "plaquette z", false, true);
45 hila::out0 << "plaquette: " << plaq_z.back() << '\n';
46}
47
48/**
49 * @brief Measures Polyakov lines and Wilson action with multicanonical weights
50 *
51 * @tparam group
52 * @param U GaugeField to measure
53 * @param p Parameter struct
54 */
55template <typename group>
56void measure_plaq_multicanonical(const GaugeField<group> &U,
57 const parameters &p) {
58
59 auto plaq = measure_plaq_with_z(
60 U, p.twist_coeff); /// (lattice.volume() * NDIM * (NDIM - 1) / 2);
61
62 print_formatted_numbers(plaq, "plaquette", false, true);
63 hila::out0 << "muca_plaquette: " << hila::muca::weight(plaq.back()) << '\n';
64}
65
66template <typename group>
67void measure_poly(const GaugeField<group> &U, const parameters &p) {
68 auto poly_z = measure_polyakov_with_z(U);
69 auto poly_abs = measure_polyakov_with_z_abs(U);
70
71 print_formatted_numbers(poly_z, "polyakov z", false, true);
72 print_formatted_numbers(poly_abs, "polyakov abs z", false, true);
73 hila::out0 << "polyakov: " << poly_z.back() << '\n';
74 hila::out0 << "polyakov abs: " << poly_abs.back() << '\n';
75}
76
77template <typename group>
78void measure_poly_multicanonical(const GaugeField<group> &U,
79 const parameters &p) {
80 auto poly = measure_polyakov_with_z(U);
81 auto poly_abs = measure_polyakov_with_z_abs(U);
82
83 print_formatted_numbers(poly, "polyakov z", false, true);
84 print_formatted_numbers(poly_abs, "polyakov abs z", false, true);
85 hila::out0 << "polyakov: " << poly.back() << '\n';
86 hila::out0 << "polyakov abs: " << poly_abs.back() << '\n';
87
88 hila::out0 << "muca_polyakov: " << hila::muca::weight(poly_abs.back())
89 << '\n';
90}
91
92/**
93 * @brief Wrapper update function
94 * @details Updates Gauge Field one direction at a time first EVEN then ODD
95 * parity
96 *
97 * @tparam group
98 * @param U GaugeField to update
99 * @param p Parameter struct
100 * @param relax If true evolves GaugeField with over relaxation if false then
101 * with heat bath
102 */
103template <typename group>
104void update(GaugeField<group> &U, const parameters &p, bool relax) {
105 foralldir(d) {
106 for (Parity par : {EVEN, ODD}) {
107
108 update_parity_dir(U, p, par, d, relax);
109 }
110 }
111}
112
113template <typename group>
114void update_multicanonical(GaugeField<group> &U, const parameters &p,
115 bool relax) {
116 foralldir(d) {
117 for (Parity par : {EVEN, ODD}) {
118 if (p.muca_poly) {
119 if (d == e_t) {
120 auto U_old = U;
121 update_parity_dir(U, p, par, d, relax);
122 auto OP = measure_polyakov_with_z_abs(U);
123 auto OP_old = measure_polyakov_with_z_abs(U_old);
124 // hila::out0 << "OP old: " << OP_old.back() << " OP new: " <<
125 // OP.back() << std::endl;
126 if (!hila::muca::accept_reject(OP_old.back(), OP.back())) {
127 U = U_old;
128 }
129 } else {
130 update_parity_dir(U, p, par, d, relax);
131 }
132 }
133 if (p.muca_action) {
134 auto U_old = U;
135 update_parity_dir(U, p, par, d, relax);
136 auto OP = measure_plaq_with_z(U, p.twist_coeff);
137 auto OP_old = measure_plaq_with_z(U_old, p.twist_coeff);
138 // hila::out0 << "OP old: " << OP_old.back() << " OP new: " << OP.back()
139 // << std::endl;
140 if (!hila::muca::accept_reject(OP_old.back(), OP.back())) {
141 U = U_old;
142 }
143 }
144 }
145 }
146}
147
148/**
149 * @brief Wrapper function to updated GaugeField per direction
150 * @details Computes first staplesum, then uses computed result to evolve
151 * GaugeField either with over relaxation or heat bath
152 *
153 * @tparam group
154 * @param U GaugeField to evolve
155 * @param p parameter struct
156 * @param par Parity
157 * @param d Direction to evolve
158 * @param relax If true evolves GaugeField with over relaxation if false then
159 * with heat bath
160 */
161template <typename group>
162void update_parity_dir(GaugeField<group> &U, const parameters &p, Parity par,
163 Direction d, bool relax) {
164
165 static hila::timer hb_timer("Heatbath");
166 static hila::timer or_timer("Overrelax");
167 static hila::timer staples_timer("Staplesum");
168
169 Field<group> staples;
170
171 staples_timer.start();
172
173 // staplesum_twist(U, staples, d, p.twist_coeff,par);
174 staplesum_twist(U, staples, d, p.twist_coeff, par);
175
176 staples_timer.stop();
177
178 if (relax) {
179
180 or_timer.start();
181 onsites(par) { suN_overrelax(U[d][X], staples[X]); }
182 or_timer.stop();
183
184 } else {
185
186 hb_timer.start();
187 onsites(par) { suN_heatbath(U[d][X], staples[X], p.beta); }
188 hb_timer.stop();
189 }
190}
191
192/**
193 * @brief Evolve gauge field
194 * @details Evolution happens by means of heat bath and over relaxation. For
195 * each heatbath update (p.n_update) we update p.n_overrelax times with over
196 * relaxation.
197 *
198 * @tparam group
199 * @param U
200 * @param p
201 */
202template <typename group>
203void do_trajectory(GaugeField<group> &U, const parameters &p) {
204 auto U_old = U;
205 for (int n = 0; n < p.n_update; n++) {
206 for (int i = 0; i < p.n_overrelax; i++) {
207 update(U, p, true);
208 }
209
210 update(U, p, false);
211 }
213}
214
215/**
216 * @brief Evolve gauge field with multicanonical update
217 * @details Evolution happens by means of heat bath and over relaxation. For
218 * each heatbath update (p.n_update) we update p.n_overrelax times with over
219 * relaxation.
220 *
221 * @tparam group
222 * @param U
223 * @param p
224 */
225template <typename group>
226void do_trajectory_multicanonical(GaugeField<group> &U, const parameters &p) {
227 for (int n = 0; n < p.n_update; n++) {
228 for (int i = 0; i < p.n_overrelax; i++) {
229 update_multicanonical(U, p, true);
230 }
231
232 update_multicanonical(U, p, false);
233 }
235}
236
237/**
238 * @brief Create weight function with multicanonical method
239 *
240 * @param U Gauge field
241 * @param p Parameters struct
242 */
243void iterate_weights_multicanonical(GaugeField<mygroup> U,
244 const parameters &p) {
245
246 hila::out0 << "here also\n";
247 bool iterate_status = true;
248 while (iterate_status) {
249 for (int i = 0; i < p.muca_updates; i++) {
250 do_trajectory_multicanonical(U, p);
251 }
252 if (p.muca_action) {
253 auto OP = measure_plaq_with_z(U, p.twist_coeff);
254 auto P = measure_polyakov_with_z(U);
255 hila::out0 << "Order parameter: " << OP.back() << std::endl;
256 hila::out0 << "polyakov muca: " << P.back() << std::endl;
257 iterate_status = hila::muca::iterate_weights(OP.back());
258 } else {
259 auto OP = measure_polyakov_with_z_abs(U);
260 auto Plaq = measure_plaq_with_z(U, p.twist_coeff);
261 auto Poly = measure_polyakov_with_z(U);
262 hila::out0 << "Order parameter: " << abs(OP.back()) << std::endl;
263 hila::out0 << "plaquette muca: "<< Plaq.back() << std::endl;
264 hila::out0 << "polyakov muca: "<< Poly.back() << std::endl;
265
266 iterate_status = hila::muca::iterate_weights(OP.back());
267 }
268 }
269 hila::muca::write_weight_function(hila::muca::generate_outfile_name());
270}
271
272int main(int argc, char **argv) {
273 parameters p;
274 CoordinateVector lsize;
275 hila::cmdline.add_flag("-beta", "Temparature");
276 hila::cmdline.add_flag("-twist", "Twist coefficient");
277 hila::cmdline.add_flag("-config", "Config filename");
278 hila::cmdline.add_flag("-ntraj", "Number of trajectories");
279 hila::cmdline.add_flag("-lsize", "Lattice size");
280 hila::cmdline.add_flag("-out-folder", "Output folder");
281 hila::initialize(argc, argv);
282
283 hila::out0 << "SU(" << mygroup::size() << ") heat bath + overrelax update\n";
284
285 hila::input par("parameters");
286
287 lsize = par.get("lattice size"); // reads NDIM numbers
288
289 p.beta = par.get("beta");
290
291 // deltab sets system to different beta on different sides, by beta*(1 +-
292 // deltab) use for initial config generation only
293 p.deltab = par.get("delta beta fraction");
294 // trajectory length in steps
295 p.n_overrelax = par.get("overrelax steps");
296 p.n_update = par.get("updates in trajectory");
297 p.n_trajectories = par.get("trajectories");
298 p.n_thermal = par.get("thermalization");
299
300 // random seed = 0 -> get seed from time
301 long seed = par.get("random seed");
302 // save config and checkpoint
303 p.n_save = par.get("traj/saved");
304 // measure surface properties and print "profile"
305 p.config_file = par.get("config name");
306 p.twist_coeff = par.get("twist coeff");
307
308 if (par.get_item("updates/profile meas", {"off", "%i"}) == 1) {
309 p.n_profile = par.get();
310 } else {
311 p.n_profile = 0;
312 }
313 if (p.n_profile) {
314 p.n_smear = par.get("smearing steps");
315 p.smear_coeff = par.get("smear coefficient");
316 p.z_smear = par.get("z smearing steps");
317 p.n_surface = par.get("traj/surface");
318 p.n_dump_polyakov = par.get("traj/polyakov dump");
319
320 if (p.n_smear.size() != p.z_smear.size()) {
321 hila::out0 << "Error in input file: number of values in 'smearing "
322 "steps' != 'z "
323 "smearing steps'\n";
325 }
326
327 } else {
328 p.n_dump_polyakov = 0;
329 }
330 p.muca_action = par.get("muca_action");
331 p.muca_poly = par.get("muca_poly");
332 p.muca_updates = par.get("muca_updates");
333 p.measure_surface = par.get("measure_surface");
334 std::string initial_state = par.get("initial_state");
335
336 par.close(); // file is closed also when par goes out of scope
337
338 if (hila::cmdline.flag_present("-beta")) {
339 p.beta = hila::cmdline.get_double("-beta");
340 hila::out0 << "beta " << p.beta << std::endl;
341 }
342 if (hila::cmdline.flag_present("-twist")) {
343 p.twist_coeff = hila::cmdline.get_int("-twist");
344 hila::out0 << "twist coeff " << p.twist_coeff << std::endl;
345 }
346 if (hila::cmdline.flag_present("-config")) {
347 p.config_file = hila::cmdline.get_string("-config");
348 }
349 if (hila::cmdline.flag_present("-ntraj")) {
350 p.n_trajectories = hila::cmdline.get_int("-ntraj");
351 hila::out0 << "Trajectories " << p.twist_coeff << std::endl;
352 }
353 if (hila::cmdline.flag_present("-lsize")) {
354 for (int i = 0; i < NDIM; i++) {
355 lsize[i] = hila::cmdline.get_int("-lsize", i);
356 }
357 hila::out0 << "Lattice size ";
358 for (int i = 0; i < NDIM; i++) {
359 hila::out0 << lsize[i] << " ";
360 }
361 hila::out0 << std::endl;
362 }
363 if (hila::cmdline.flag_present("-out-folder")) {
364 p.out_folder = hila::cmdline.get_string("-out-folder");
365 } else {
366 p.out_folder = "./";
367 }
368 // setting up the lattice is convenient to do after reading
369 // the parameter
370 lattice.setup(lsize);
371 // Alloc gauge field
372
373 // U = initial_state;
374 // use negative trajectory for thermal
375 int start_traj = -p.n_thermal;
376
377 hila::timer update_timer("Updates");
378 hila::timer measure_timer("Measurements");
379 hila::timer muca_timer("Muca");
381
382 // We need random number here
383 if (!hila::is_rng_seeded())
384 hila::seed_random(seed);
385
386 if (initial_state == "cold") {
387 foralldir(d) onsites(ALL) { U[d][X].gaussian_random(); }
388 } else if (initial_state == "hot") {
389 U = 1;
390 } else if (initial_state == "both") {
391 foralldir(d) onsites(EVEN) { U[d][X].gaussian_random(); }
392 foralldir(d) onsites(ODD) { U[d][X] = 1; }
393 }
394
395 // restore_checkpoint(U, start_traj, p);
396
397 // muca_timer.start();
398
399 if (p.muca_action || p.muca_poly) {
400 hila::muca::initialise("muca_parameters");
401 std::ofstream MFile;
402 MFile.open("muca_measurements", std::ios_base::app);
403 iterate_weights_multicanonical(U, p);
404 }
405 // muca_timer.stop();
406 hila::out0 << "MEASURE start\n";
407 void (*do_trajectory_ptr)(GaugeField<mygroup> &, const parameters &);
408 void (*measure_plaquette_ptr)(const GaugeField<mygroup> &,
409 const parameters &);
410 void (*measure_polyakov_ptr)(const GaugeField<mygroup> &, const parameters &);
411 if (p.muca_action) {
412 do_trajectory_ptr = do_trajectory_multicanonical<mygroup>;
413 measure_plaquette_ptr = measure_plaq_multicanonical<mygroup>;
414 measure_polyakov_ptr = measure_poly<mygroup>;
415 } else if (p.muca_poly) {
416 do_trajectory_ptr = do_trajectory_multicanonical<mygroup>;
417 measure_plaquette_ptr = measure_plaq<mygroup>;
418 measure_polyakov_ptr = measure_poly_multicanonical<mygroup>;
419 } else {
420 do_trajectory_ptr = do_trajectory<mygroup>;
421 measure_plaquette_ptr = measure_plaq<mygroup>;
422 measure_polyakov_ptr = measure_poly<mygroup>;
423 }
424
425 for (int trajectory = start_traj; trajectory < p.n_trajectories;
426 trajectory++) {
427 double ttime = hila::gettime();
428
429 update_timer.start();
430
431 double acc = 0;
432 do_trajectory_ptr(U, p);
433
434 // put sync here in order to get approx gpu timing
435 hila::synchronize_threads();
436 update_timer.stop();
437
438 // trajectory is negative during thermalization
439 if (trajectory >= 0) {
440 measure_timer.start();
441
442 measure_plaquette_ptr(U, p);
443
444 measure_polyakov_ptr(U, p);
445
446 if (p.measure_surface)
447 measure_polyakov_surface(U, p, trajectory);
448
449 // hila::out0 << "Measure_end " << trajectory << std::endl;
450
451 measure_timer.stop();
452 }
453
454 if (p.n_save > 0 && (trajectory + 1) % p.n_save == 0) {
455 checkpoint(U, trajectory, p);
456 }
457 }
458
459 hila::out0 << "MEASURE end\n";
460 hila::out0 << expi(4.0 / 3.0 * M_PI) << std::endl;
462}
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
static constexpr int size()
Returns size of Vector or square Matrix.
Definition matrix.h:248
hila::input - Class for parsing runtime parameter files.
Definition input.h:52
void setup(const CoordinateVector &siz)
General lattice setup.
Definition lattice.cpp:33
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,...
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
constexpr Parity ALL
bit pattern: 011
Checkpointing functions for saving and restoring lattice configurations.
void write_weight_function(string W_function_filename)
Reads the precomputed weight function from run_parameters struct and saves it into a file.
void initialise(const string wfile_name)
Loads parameters and weights for the multicanonical computation.
double weight(double OP)
process 0 interface to "weight function" for the user accessing the weights.
Header for model agnostic implementation of various multicanonical (muca) methods.
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
double gettime()
Definition timing.cpp:167
bool is_rng_seeded()
Check if RNG is seeded already.
Definition random.cpp:241
void finishrun()
Normal, controlled exit - all nodes must call this. Prints timing information and information about c...
void terminate(int status)
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.
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