CHROMA
axgauge.cc
Go to the documentation of this file.
1 /*! \file
2  * \brief Axial gauge fixing
3  */
4 
5 #include "chromabase.h"
6 #include "meas/gfix/axgauge.h"
7 
8 namespace Chroma
9 {
10 
11  //! Axial gauge fixing
12  /*!
13  * \ingroup gfix
14  *
15  * Transfroms a gauge configuration, in place, into axial
16  * gauge with special direction decay_dir.
17  *
18  * Note: The non-unity time-like gauge fields from the last time slice
19  * will be copied to all other time slices.
20  *
21  * \param ug gauge field and its axial gauge transform (Modify)
22  * \param v gauge rotation matrix (Write)
23  * \param decay_dir time direction (Read)
24  */
25 
26 
27  void axGauge(multi1d<LatticeColorMatrix>& ug, int decay_dir)
28  {
29  START_CODE();
30 
31  START_CODE();
32 
33  int lsizet = Layout::lattSize()[decay_dir];
34 
35  LatticeColorMatrix v = 1;
36  LatticeInteger t_coord = Layout::latticeCoordinate(decay_dir);
37 
38  /* Transform the "time-like" links to unity, slice by slice except for the */
39  /* last slice and thereby construct the gauge transformation V. */
40  for(int t = 1; t < lsizet; ++t)
41  {
42  LatticeBoolean btmp = t_coord == t;
43 
44  LatticeColorMatrix tmp_1 = shift(ug[decay_dir], BACKWARD, decay_dir);
45  copymask(v, btmp, tmp_1);
46  LatticeColorMatrix tmp_2 = tmp_1 * ug[decay_dir];
47  copymask(ug[decay_dir], btmp, tmp_2);
48  }
49 
50  /* Now do the gauge transformation on the space-like links */
51  for(int mu = 0; mu < Nd; ++mu)
52  {
53  if ( mu != decay_dir )
54  {
55  LatticeColorMatrix tmp_2 = ug[mu] * shift(adj(v), FORWARD, mu);
56  ug[mu] = v * tmp_2;
57  }
58  }
59 
60  /* Finally "broadcast" the t-link from the last time slice to all others */
61  for(int t = lsizet-2; t >= 0; --t)
62  {
63  LatticeBoolean btmp = t_coord == t;
64 
65  copymask(ug[decay_dir], btmp, LatticeColorMatrix(shift(ug[decay_dir], FORWARD, decay_dir)));
66  }
67 
68 
69  END_CODE();
70  }
71 
72 } // end namespace Chroma
Axial gauge fixing.
Primary include file for CHROMA library code.
int mu
Definition: cool.cc:24
void axGauge(multi1d< LatticeColorMatrix > &ug, int decay_dir)
Axial gauge fixing.
Definition: axgauge.cc:27
LatticeColorMatrix tmp_2
Definition: meslate.cc:51
int t
Definition: meslate.cc:37
LatticeColorMatrix tmp_1
Definition: meslate.cc:50
Nd
Definition: meslate.cc:74
Asqtad Staggered-Dirac operator.
Definition: klein_gord.cc:10
START_CODE()
#define FORWARD
Definition: primitives.h:82
#define BACKWARD
Definition: primitives.h:83
copymask(lcoord, lbit, ltmp_1, REPLACE)