|
CHROMA
|
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... | |
Support for gauge fixing
| 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.
| ug | gauge field and its axial gauge transform (Modify) |
| v | gauge rotation matrix (Write) |
| decay_dir | time 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.
| ug | gauge field and its axial gauge transform (Modify) |
| decay_dir | time 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().
| 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)!
| u | (gauge fixed) gauge field ( Modify ) |
| n_gf | number of gauge fixing iterations ( Write ) |
| j_decay | direction perpendicular to slices to be gauge fixed ( Read ) |
| GFAccu | desired accuracy for gauge fixing ( Read ) |
| GFMax | maximal number of gauge fixing iterations ( Read ) |
| OrDo | use overrelaxation or not ( Read ) |
| OrPara | overrelaxation parameter ( Read ) |
Definition at line 37 of file coulgauge.cc.
References j_decay, and Chroma::u.
Referenced by main(), and Chroma::InlineCoulGaugeEnv::InlineMeas::operator()().
| 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)!
| u | (gauge fixed) gauge field ( Modify ) |
| g | Gauge transformation matrices (Write) |
| n_gf | number of gauge fixing iterations ( Write ) |
| j_decay | direction perpendicular to slices to be gauge fixed ( Read ) |
| GFAccu | desired accuracy for gauge fixing ( Read ) |
| GFMax | maximal number of gauge fixing iterations ( Read ) |
| OrDo | use overrelaxation or not ( Read ) |
| OrPara | overrelaxation 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().
| 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".
| g | Current (global) gauge transformation matrices ( Modify ) |
| u | original gauge field ( Read ) |
| j_decay | direction perpendicular to slices to be gauge fixed ( Read ) |
| su2_index | SU(2) subgroup index ( Read ) |
| cb | checkerboard index ( Read ) |
| ordo | use overrelaxation or not ( Read ) |
| orpara | overrelaxation 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().
| 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
| c | complex Nc x Nc matrix ( Modify ) on exit it contains the hermitian matrix P |
| v | the projected SU(Nc) Matrix ( Write ) |
| alpha | the phase ( Write ) |
| JacAccu | accuracy in the Jacobi iteration ( Read ) |
| JacMax | maximum 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().
| 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.
| g | Gauge transformation (Write) |
| psi | Input color std::vector field (Read) |
| chi | Output color std::vector field (Write) |
| s_index | color 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().
| 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.
| ug | gauge field and its axial gauge transform (Modify) |
| g | gauge rotation matrix (Write) |
| decay_dir | time 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.