CHROMA
Functions

Functions

void Chroma::axGauge (multi1d< LatticeColorMatrix > &ug, int decay_dir)
 Axial gauge fixing. More...
 
void Chroma::coulGauge (multi1d< LatticeColorMatrix > &u, int &n_gf, int j_decay, const Real &GFAccu, int GFMax, bool OrDo, const Real &OrPara)
 Coulomb (and Landau) gauge fixing. More...
 
void Chroma::coulGauge (multi1d< LatticeColorMatrix > &u, LatticeColorMatrix &g, int &n_gf, int j_decay, const Real &GFAccu, int GFMax, bool OrDo, const Real &OrPara)
 Coulomb (and Landau) gauge fixing. More...
 
void Chroma::grelax (LatticeColorMatrix &g, const multi1d< LatticeColorMatrix > &u, int j_decay, int su2_index, int cb, bool ordo, const Real &orpara)
 Perform a single gauge fixing iteration. More...
 
void Chroma::polar_dec (LatticeColorMatrix &c, LatticeColorMatrix &v, LatticeReal &alpha, const Real &JacAccu, int JacMax)
 Decompose a complex matrix as C = exp(i\alpha) V P. More...
 
void Chroma::rot_colvec (LatticeColorMatrix &g, const LatticeColorVector &psi, LatticeColorVector &chi, int s_index)
 Rotate a color std::vector. More...
 
void Chroma::temporalGauge (multi1d< LatticeColorMatrix > &ug, LatticeColorMatrix &g, int decay_dir)
 Temporal gauge fixing. More...
 

Detailed Description

Support for gauge fixing

Function Documentation

◆ axGauge()

void Chroma::axGauge ( multi1d< LatticeColorMatrix > &  ug,
int  decay_dir 
)

Axial gauge fixing.

Transfroms a gauge configuration, in place, into axial gauge with special direction decay_dir.

Note: The non-unity time-like gauge fields from the last time slice will be copied to all other time slices.

Parameters
uggauge field and its axial gauge transform (Modify)
vgauge rotation matrix (Write)
decay_dirtime direction (Read)

Transfroms a gauge configuration, in place, into axial gauge with special direction decay_dir.

Note: The non-unity time-like gauge fields from the last time slice will be copied to all other time slices.

Parameters
uggauge field and its axial gauge transform (Modify)
decay_dirtime direction (Read)

Definition at line 27 of file axgauge.cc.

References BACKWARD, copymask(), Chroma::END_CODE(), FORWARD, mu, Nd, Chroma::START_CODE(), t, tmp_1, and tmp_2.

Referenced by Chroma::wilslp().

◆ coulGauge() [1/2]

void Chroma::coulGauge ( multi1d< LatticeColorMatrix > &  u,
int &  n_gf,
int  j_decay,
const Real &  GFAccu,
int  GFMax,
bool  OrDo,
const Real &  OrPara 
)

Coulomb (and Landau) gauge fixing.

Driver for gauge fixing to Coulomb gauge in slices perpendicular to the direction "j_decay". If j_decay >= Nd: fix to Landau gauge. Note: as written this works only for SU(2) and SU(3)!

Parameters
u(gauge fixed) gauge field ( Modify )
n_gfnumber of gauge fixing iterations ( Write )
j_decaydirection perpendicular to slices to be gauge fixed ( Read )
GFAccudesired accuracy for gauge fixing ( Read )
GFMaxmaximal number of gauge fixing iterations ( Read )
OrDouse overrelaxation or not ( Read )
OrParaoverrelaxation parameter ( Read )

Definition at line 37 of file coulgauge.cc.

References j_decay, and Chroma::u.

Referenced by main(), and Chroma::InlineCoulGaugeEnv::InlineMeas::operator()().

◆ coulGauge() [2/2]

void Chroma::coulGauge ( multi1d< LatticeColorMatrix > &  u,
LatticeColorMatrix &  g,
int &  n_gf,
int  j_decay,
const Real &  GFAccu,
int  GFMax,
bool  OrDo,
const Real &  OrPara 
)

Coulomb (and Landau) gauge fixing.

Driver for gauge fixing to Coulomb gauge in slices perpendicular to the direction "j_decay". If j_decay >= Nd: fix to Landau gauge. Note: as written this works only for SU(2) and SU(3)!

Parameters
u(gauge fixed) gauge field ( Modify )
gGauge transformation matrices (Write)
n_gfnumber of gauge fixing iterations ( Write )
j_decaydirection perpendicular to slices to be gauge fixed ( Read )
GFAccudesired accuracy for gauge fixing ( Read )
GFMaxmaximal number of gauge fixing iterations ( Read )
OrDouse overrelaxation or not ( Read )
OrParaoverrelaxation parameter ( Read )

Definition at line 68 of file coulgauge.cc.

References Chroma::cb, Chroma::END_CODE(), FORWARD, Chroma::grelax(), j_decay, mu, Nd, norm, Chroma::pop(), Chroma::push(), Chroma::reunit(), Chroma::START_CODE(), su2_index, sum, Chroma::tDir(), Chroma::u, Chroma::write(), and Chroma::xi_0().

◆ grelax()

void Chroma::grelax ( LatticeColorMatrix &  g,
const multi1d< LatticeColorMatrix > &  u,
int  j_decay,
int  su2_index,
int  cb,
bool  ordo,
const Real &  orpara 
)

Perform a single gauge fixing iteration.

Performs one gauge fixing 'iteration', one checkerboard and SU(2) subgroup only, for gauge fixing to Coulomb gauge in slices perpendicular to the direction "j_decay".

Parameters
gCurrent (global) gauge transformation matrices ( Modify )
uoriginal gauge field ( Read )
j_decaydirection perpendicular to slices to be gauge fixed ( Read )
su2_indexSU(2) subgroup index ( Read )
cbcheckerboard index ( Read )
ordouse overrelaxation or not ( Read )
orparaoverrelaxation parameter ( Read )

Definition at line 37 of file grelax.cc.

References Chroma::a, Chroma::anisoP(), BACKWARD, Chroma::cb, Chroma::END_CODE(), FORWARD, j_decay, mu, Nd, Chroma::r, Chroma::START_CODE(), su2_index, Chroma::su2Extract(), Chroma::sunFill(), Chroma::tDir(), Chroma::twopi, Chroma::u, and Chroma::xi_0().

Referenced by Chroma::coulGauge().

◆ polar_dec()

void Chroma::polar_dec ( LatticeColorMatrix &  c,
LatticeColorMatrix &  v,
LatticeReal &  alpha,
const Real &  JacAccu,
int  JacMax 
)

Decompose a complex matrix as C = exp(i\alpha) V P.

Decompose a complex matrix as C = exp(i\alpha) V P with V SU(Nc) and P = (C^\dagger C)^{1/2} positive hermitian

Parameters
ccomplex Nc x Nc matrix ( Modify ) on exit it contains the hermitian matrix P
vthe projected SU(Nc) Matrix ( Write )
alphathe phase ( Write )
JacAccuaccuracy in the Jacobi iteration ( Read )
JacMaxmaximum number of Jacobi iterations ( Read )

Definition at line 25 of file polar_dec.cc.

References Chroma::StagPhases::alpha(), Chroma::c, copymask(), Chroma::END_CODE(), Chroma::i, j, Chroma::k, numbad, Chroma::pop(), Chroma::push(), Chroma::reunit(), Chroma::REUNITARIZE_LABEL, Chroma::START_CODE(), sum, and Chroma::write().

Referenced by Chroma::HisqFermAct::createState().

◆ rot_colvec()

void Chroma::rot_colvec ( LatticeColorMatrix &  g,
const LatticeColorVector &  psi,
LatticeColorVector &  chi,
int  s_index 
)

Rotate a color std::vector.

Rotate a color std::vector into the form where the component with index s_index is real and the component with larger index are zero. We do this by a series of SU(2) gauge rotations.

Parameters
gGauge transformation (Write)
psiInput color std::vector field (Read)
chiOutput color std::vector field (Write)
s_indexcolor index (Read)

Definition at line 25 of file rot_colvec.cc.

References Chroma::a, Chroma::chi(), Chroma::END_CODE(), Chroma::i, Chroma::psi, Chroma::QDP_error_exit(), Chroma::START_CODE(), and Chroma::sunFill().

◆ temporalGauge()

void Chroma::temporalGauge ( multi1d< LatticeColorMatrix > &  ug,
LatticeColorMatrix &  g,
int  decay_dir 
)

Temporal gauge fixing.

Transfroms a gauge configuration, in place, into axial gauge with special direction decay_dir.

Note: The non-unity time-like gauge fields from the last time slice will be copied to all other time slices.

Parameters
uggauge field and its axial gauge transform (Modify)
ggauge rotation matrix (Write)
decay_dirtime direction (Read)

Definition at line 26 of file temporal_gauge.cc.

References BACKWARD, copymask(), Chroma::END_CODE(), FORWARD, mu, Nd, Chroma::START_CODE(), t, Chroma::tmp, and Chroma::u.