CHROMA
temporal_gauge.cc
Go to the documentation of this file.
1 /*! \file
2  * \brief Axial gauge fixing
3  */
4 
5 #include "chromabase.h"
7 
8 namespace Chroma
9 {
10 
11  //! Temporal 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 g gauge rotation matrix (Write)
23  * \param decay_dir time direction (Read)
24  */
25 
26  void temporalGauge(multi1d<LatticeColorMatrix>& ug, LatticeColorMatrix& g, int decay_dir)
27  {
28  START_CODE();
29 
30  // Sanity check
31  if (decay_dir < 0 || decay_dir >= ug.size())
32  {
33  QDPIO::cerr << __func__ << ": invalid decay_dir= " << decay_dir << std::endl;
34  QDP_abort(1);
35  }
36 
37  // Initialize
38  int N_t = Layout::lattSize()[decay_dir];
39  LatticeInteger t_coord = Layout::latticeCoordinate(decay_dir);
40 
41  {
42  // Upcast to double precision
43  LatticeColorMatrixD u = ug[decay_dir];
44  LatticeColorMatrixD gt = 1; // Unitm matrix
45 
46  // U_prev holds U_{t-1}
47  LatticeColorMatrixD U_prev = shift(u, BACKWARD, decay_dir);
48  LatticeColorMatrixD G_prev = gt; // Not shifting becuse it is unit at this point
49 
50  for(int t = 1; t < N_t; ++t) {
51  LatticeBoolean btmp = (t_coord == t);
52  LatticeColorMatrixD t1 = G_prev*U_prev ; // DO whole lattice... Sucks
53  copymask(gt, btmp, t1);
54 
55  // Reuse t1, to compute G_prev for next t
56  t1 = shift(gt, BACKWARD, decay_dir);
57  G_prev = t1;
58  }
59 
60  g=gt; // Save gauge transform matrix (downcast to base precision)
61  }
62 
63  /* Now do the gauge transformation on all the links */
64  for(int mu = 0; mu < Nd; ++mu) {
65  LatticeColorMatrix tmp = ug[mu] * shift(adj(g), FORWARD, mu);
66  ug[mu] = g * tmp;
67  }
68 
69  // Check temporal links (except last timeslice are unit
70  {
71  // Take a copy
72  LatticeColorMatrix u = ug[decay_dir];
73  LatticeColorMatrix g_unit = 1; // For copying into the last timeslice
74 
75  // Blast last non-unit timeslice with a unit one...
76  LatticeBoolean btmp = (t_coord == (N_t-1) );
77  copymask(u, btmp, g_unit); // Replace last timeslice with units...
78 
79  // U should now be indistinguishable from unit gauge
80  QDPIO::cout << "Norm of Unit-Gauge violation / link = " << sqrt( norm2(u-g_unit) )/ Layout::vol() << std::endl;
81 
82  }
83 
84  END_CODE();
85  }
86 
87 } // end namespace Chroma
Primary include file for CHROMA library code.
int mu
Definition: cool.cc:24
void temporalGauge(multi1d< LatticeColorMatrix > &ug, LatticeColorMatrix &g, int decay_dir)
Temporal gauge fixing.
int t
Definition: meslate.cc:37
Nd
Definition: meslate.cc:74
Asqtad Staggered-Dirac operator.
Definition: klein_gord.cc:10
static multi1d< LatticeColorMatrix > u
LatticeFermion tmp
Definition: mespbg5p_w.cc:36
START_CODE()
#define FORWARD
Definition: primitives.h:82
#define BACKWARD
Definition: primitives.h:83
copymask(lcoord, lbit, ltmp_1, REPLACE)
Axial gauge fixing.