37 std::string config_file;
54template <
typename group>
56 for (
int i = 0; i < 2 * NDIM; ++i) {
59 int tpar = 1 + (tdp % 2);
78template <
typename group>
88 staples_timer.start();
98#ifdef SUN_OVERRELAX_dFJ
99 suN_overrelax_dFJ(U[d][X], staples[X], p.beta);
125template <
typename group>
128 for (
int n = 0; n < p.n_update; n++) {
129 for (
int i = 0; i <= p.n_overrelax; i++) {
142template <
typename group>
146 static bool first =
true;
152 auto plaq = measure_s_wplaq(U) / (lattice.volume() * NDIM * (NDIM - 1) / 2);
154 hila::out0 << string_format(
"MEAS % 0.6e % 0.6e % 0.6e", plaq, poly.real(), poly.imag())
162template <
typename group>
163void checkpoint(
const GaugeField<group> &U,
int trajectory,
const parameters &p) {
166 std::string config_file =
167 p.config_file +
"_" + std::to_string(
abs((trajectory + 1) / p.n_save) % 2);
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';
178 outf <<
"config name " << config_file <<
'\n';
181 std::stringstream msg;
183 hila::timestamp(msg.str().c_str());
186template <
typename group>
192 if (status.
open(
"run_status",
false,
false)) {
194 trajectory = status.
get(
"trajectory");
195 seed = status.
get(
"seed");
196 p.time_offset = status.
get(
"time");
198 std::string config_file = status.
get(
"config name");
201 U.config_read(config_file);
205 in.open(p.config_file, std::ios::in | std::ios::binary);
209 U.config_read(p.config_file);
222int main(
int argc,
char **argv) {
240 hila::out0 <<
"Using floating point epsilon: " << fp<ftype>::epsilon <<
"\n";
248 lsize = par.get(
"lattice size");
250 p.beta = par.get(
"beta");
252 p.n_traj = par.get(
"number of trajectories");
254 p.n_update = par.get(
"updates in trajectory");
256 p.n_overrelax = par.get(
"overrelax steps");
258 p.n_therm = par.get(
"thermalization trajs");
260 p.gflow_freq = par.get(
"gflow freq");
262 p.gflow_max_l = par.get(
"gflow max lambda");
264 p.gflow_l_step = par.get(
"gflow lambda step");
266 p.gflow_a_accu = par.get(
"gflow abs. accuracy");
268 p.gflow_r_accu = par.get(
"gflow rel. accuracy");
270 p.num_csflows = par.get(
"number of csflows");
272 p.csflow_t_step = par.get(
"csflow t step");
274 p.csflow_a_accu = par.get(
"csflow abs. accuracy");
276 p.csflow_r_accu = par.get(
"csflow rel. accuracy");
278 long seed = par.get(
"random seed");
280 p.n_save = par.get(
"trajs/saved");
282 p.config_file = par.get(
"config name");
287 lattice.
setup(lsize);
297 int start_traj = -p.n_therm;
299 if (!restore_checkpoint(U, start_traj, p)) {
310 for (
int trajectory = start_traj; trajectory < p.n_traj; ++trajectory) {
314 update_timer.start();
319 hila::synchronize_threads();
322 measure_timer.start();
324 hila::out0 <<
"Measure_start " << trajectory <<
'\n';
328 hila::out0 <<
"Measure_end " << trajectory <<
'\n';
330 measure_timer.stop();
332 if (trajectory >= 0) {
334 if (p.gflow_freq > 0 && trajectory % p.gflow_freq == 0) {
337 int gtrajectory = trajectory / p.gflow_freq;
338 if (p.gflow_l_step > 0) {
341 for (
int ics = 0; ics < p.num_csflows; ++ics) {
347 do_cs_flow_adapt(U, E, p.csflow_t_step, p.csflow_a_accu, p.csflow_r_accu);
355 int nflow_steps = (int)(p.gflow_max_l / p.gflow_l_step);
358 hila::out0 <<
"Gflow_start " << gtrajectory <<
'\n';
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);
369 for (
int i = 1; i < nflow_steps; ++i) {
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);
375 measure_gradient_flow_stuff(V, (i + 1) * p.gflow_l_step, t_step);
378 hila::out0 <<
"Gflow_end " << gtrajectory <<
" time " << std::setprecision(3)
388 if (p.n_save > 0 && (trajectory + 1) % p.n_save == 0) {
389 checkpoint(U, trajectory, p);
The field class implements the standard methods for accessing Fields. Hilapp replaces the parity acce...
void reunitarize_gauge()
Reunitarize Gauge Field consisting of matrices.
void config_write(const std::string &filename) const
config_write writes the gauge field to file, with additional "verifying" header
static constexpr int size()
Returns size of Vector or square Matrix.
void setup(const CoordinateVector &siz)
General lattice setup.
T abs(const Complex< T > &a)
Return absolute value of Complex number.
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)
Direction
Enumerator for direction that assigns integer to direction to be interpreted as unit vector.
constexpr Parity ALL
bit pattern: 011
double random()
Real valued uniform random number generator.
int myrank()
rank of this node
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...
T broadcast(T &var, int rank=0)
Broadcast the value of var to all MPI ranks from rank (default=0).
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.
void staplesum(const GaugeField< T > &U, Field< T > &staples, Direction d1, Parity par=ALL)
Sum the staples of link matrices to direction dir.
void update(GaugeField< group > &U, const parameters &p, bool relax)
Wrapper update function.
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.
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