CHROMA
tprec_wilson_linop_w.cc
Go to the documentation of this file.
1 /*! \file
2  * \brief Even-odd preconditioned Wilson linear operator
3  */
4 
5 #include "chromabase.h"
6 #include "actions/ferm/linop/prec_wilson_linop_w.h"
7 
8 using namespace QDP::Hints;
9 namespace Chroma
10 {
11  //! Creation routine
12  /*!
13  * \param u_ gauge field (Read)
14  * \param Mass_ fermion kappa (Read)
15  */
16  void EvenOddPrecWilsonLinOp::create(Handle< FermState<T,P,Q> > fs,
17  const Real& Mass_)
18  {
19  AnisoParam_t anisoParam;
20  create(fs, Mass_, anisoParam);
21  }
22 
23 
24  //! Creation routine with Anisotropy
25  /*!
26  * \param u_ gauge field (Read)
27  * \param Mass_ fermion kappa (Read)
28  * \param aniso anisotropy struct (Read)
29  */
30  void EvenOddPrecWilsonLinOp::create(Handle< FermState<T,P,Q> > fs,
31  const Real& Mass_,
32  const AnisoParam_t& anisoParam)
33  {
34  START_CODE();
35 
36  D.create(fs,anisoParam);
37 
38  Mass = Mass_;
39  Real ff = where(anisoParam.anisoP, anisoParam.nu / anisoParam.xi_0, Real(1));
40  fact = 1 + (Nd-1)*ff + Mass;
41  invfact = 1/fact;
42 
43  END_CODE();
44  }
45 
46 
47  //! Apply even-odd linop component
48  /*!
49  * The operator acts on the entire even sublattice
50  *
51  * \param chi Pseudofermion field (Write)
52  * \param psi Pseudofermion field (Read)
53  * \param isign Flag ( PLUS | MINUS ) (Read)
54  */
55  void
56  EvenOddPrecWilsonLinOp::evenOddLinOp(LatticeFermion& chi,
57  const LatticeFermion& psi,
58  enum PlusMinus isign) const
59  {
60  START_CODE();
61 
62  Real mhalf = -0.5;
63 
64  D.apply(chi, psi, isign, 0);
65  chi[rb[0]] *= mhalf;
66 
67  END_CODE();
68  }
69 
70  //! Apply odd-even linop component
71  /*!
72  * The operator acts on the entire odd sublattice
73  *
74  * \param chi Pseudofermion field (Write)
75  * \param psi Pseudofermion field (Read)
76  * \param isign Flag ( PLUS | MINUS ) (Read)
77  */
78  void
79  EvenOddPrecWilsonLinOp::oddEvenLinOp(LatticeFermion& chi,
80  const LatticeFermion& psi,
81  enum PlusMinus isign) const
82  {
83  START_CODE();
84 
85  Real mhalf = -0.5;
86 
87  D.apply(chi, psi, isign, 1);
88  chi[rb[1]] *= mhalf;
89 
90  END_CODE();
91  }
92 
93 
94 
95  //! Override inherited one with a few more funkies
96  void EvenOddPrecWilsonLinOp::operator()(LatticeFermion & chi,
97  const LatticeFermion& psi,
98  enum PlusMinus isign) const
99  {
100  START_CODE();
101 
102  LatticeFermion tmp1, tmp2, tmp3; // if an array is used here,
103 
104  moveToFastMemoryHint(tmp1);
105  moveToFastMemoryHint(tmp2);
106  moveToFastMemoryHint(tmp3);
107 
108 
109  Real mquarterinvfact = -0.25*invfact;
110 
111  // tmp1[0] = D_eo psi[1]
112  D.apply(tmp1, psi, isign, 0);
113 
114  // tmp2[1] = D_oe tmp1[0]
115  D.apply(tmp2, tmp1, isign, 1);
116 
117  // Now we have tmp2[1] = D_oe D_eo psi[1]
118 
119  // now scale tmp2[1] with (-1/4)/fact = (-1/4)*(1/(Nd + m))
120  // with a vscale -- using tmp2 on both sides should be OK, but
121  // just to be safe use tmp3
122  tmp3[rb[1]] = mquarterinvfact*tmp2;
123 
124  // now tmp3[1] should be = (-1/4)*(1/(Nd + m) D_oe D_eo psi[1]
125 
126  // Now get chi[1] = fact*psi[1] + tmp3[1]
127  // with a vaxpy3
128  // chi[1] = (Nd + m) - (1/4)*(1/(Nd + m)) D_oe D_eo psi[1]
129  //
130  // in this order, this last job could be replaced with a
131  // vaxpy3_norm if we wanted the || M psi ||^2 over the subset.
132  chi[rb[1]] = fact*psi + tmp3;
133 
134  getFermBC().modifyF(chi, rb[1]);
135 
136  END_CODE();
137  }
138 
139 
140  //! Derivative of even-odd linop component
141  void
142  EvenOddPrecWilsonLinOp::derivEvenOddLinOp(multi1d<LatticeColorMatrix>& ds_u,
143  const LatticeFermion& chi,
144  const LatticeFermion& psi,
145  enum PlusMinus isign) const
146  {
147  START_CODE();
148 
149  ds_u.resize(Nd);
150 
151  D.deriv(ds_u, chi, psi, isign, 0);
152  for(int mu=0; mu < Nd; mu++) {
153  ds_u[mu] *= Real(-0.5);
154  }
155 
156  END_CODE();
157  }
158 
159 
160  //! Derivative of odd-even linop component
161  void
162  EvenOddPrecWilsonLinOp::derivOddEvenLinOp(multi1d<LatticeColorMatrix>& ds_u,
163  const LatticeFermion& chi,
164  const LatticeFermion& psi,
165  enum PlusMinus isign) const
166  {
167  START_CODE();
168 
169  ds_u.resize(Nd);
170 
171  D.deriv(ds_u, chi, psi, isign, 1);
172  for(int mu=0; mu < Nd; mu++) {
173  ds_u[mu] *= Real(-0.5);
174  }
175  END_CODE();
176  }
177 
178  //! Return flops performed by the operator()
179  unsigned long EvenOddPrecWilsonLinOp::nFlops() const
180  {
181  unsigned long cbsite_flops = 2*D.nFlops()+6*Nc*Ns;
182  return cbsite_flops*(Layout::sitesOnNode()/2);
183  }
184 
185 } // End Namespace Chroma
Primary include file for CHROMA library code.
#define END_CODE()
Definition: chromabase.h:65
#define START_CODE()
Definition: chromabase.h:64
int mu
Definition: cool.cc:24
Real Mass
Nd
Definition: meslate.cc:74
Double tmp2
Definition: mesq.cc:30
Double tmp3
Definition: mesq.cc:31
multi1d< Hadron2PtContraction_t > operator()(const multi1d< LatticeColorMatrix > &u)
Asqtad Staggered-Dirac operator.
Definition: klein_gord.cc:10
chi
Definition: pade_trln_w.cc:24
psi
Definition: pade_trln_w.cc:191