CHROMA
syssolver_linop_clover_qphix_iter_refine_w.h
Go to the documentation of this file.
1 // -*- C++ -*-
2 /*! \file
3  * \brief Solve a M*psi=chi linear system by BiCGStab
4  */
5 
6 #ifndef __syssolver_linop_clover_qphix_iter_refine_h__
7 #define __syssolver_linop_clover_qphix_iter_refine_h__
8 
9 #include "chroma_config.h"
10 
11 #ifdef BUILD_QPHIX
12 #include "handle.h"
13 #include "state.h"
14 #include "syssolver.h"
15 #include "linearop.h"
22 #include "io/aniso_io.h"
23 #include <string>
24 
25 #include "chroma_config.h"
26 #include "util/gauge/reunit.h"
27 #include "io/aniso_io.h"
28 
29 // Header files from the dslash package
30 #include "qphix/geometry.h"
31 #include "qphix/qdp_packer.h"
32 #include "qphix/clover.h"
33 #include "qphix/invbicgstab.h"
34 #include "qphix/inv_richardson_multiprec.h"
36 #include "qphix_singleton.h"
37 
38 namespace Chroma
39 {
40 
41  //! Richardson system solver namespace
42  namespace LinOpSysSolverQPhiXCloverIterRefineEnv
43  {
44  //! Register the syssolver
45  bool registerAll();
46 
47  }
48 
49 
50 
51  //! Solve a Clover Fermion System using the QPhiX inverter
52  /*! \ingroup invert
53  *** WARNING THIS SOLVER WORKS FOR Clover FERMIONS ONLY ***
54  */
55  using namespace QPhiXVecTraits;
56  template<typename T, typename U>
57  class LinOpSysSolverQPhiXCloverIterRefine : public LinOpSystemSolver<T>
58  {
59  public:
60  typedef multi1d<U> Q;
61  typedef typename WordType<T>::Type_t REALT;
62  typedef CHROMA_QPHIX_INNER_TYPE InnerReal;
63 
64  typedef typename QPhiX::Geometry<REALT,MixedVecTraits<REALT,InnerReal>::Vec,MixedVecTraits<REALT,InnerReal>::Soa,MixedVecTraits<REALT,InnerReal>::compress12>::FourSpinorBlock QPhiX_Spinor;
65  typedef typename QPhiX::Geometry<REALT,MixedVecTraits<REALT,InnerReal>::Vec,MixedVecTraits<REALT,InnerReal>::Soa,MixedVecTraits<REALT,InnerReal>::compress12>::SU3MatrixBlock QPhiX_Gauge;
66  typedef typename QPhiX::Geometry<REALT,MixedVecTraits<REALT,InnerReal>::Vec,MixedVecTraits<REALT,InnerReal>::Soa,MixedVecTraits<REALT,InnerReal>::compress12>::CloverBlock QPhiX_Clover;
67 
68  typedef typename QPhiX::Geometry<InnerReal,MixedVecTraits<REALT,InnerReal>::VecInner,
69  MixedVecTraits<REALT,InnerReal>::SoaInner,
70  MixedVecTraits<REALT,InnerReal>::compress12>::FourSpinorBlock QPhiX_InnerSpinor;
71 
72  typedef typename QPhiX::Geometry<InnerReal,MixedVecTraits<REALT,InnerReal>::VecInner,
73  MixedVecTraits<REALT,InnerReal>::SoaInner,
74  MixedVecTraits<REALT,InnerReal>::compress12>::SU3MatrixBlock QPhiX_InnerGauge;
75 
76  typedef typename QPhiX::Geometry<InnerReal,MixedVecTraits<REALT,InnerReal>::VecInner,
77  MixedVecTraits<REALT,InnerReal>::SoaInner,
78  MixedVecTraits<REALT,InnerReal>::compress12>::CloverBlock QPhiX_InnerClover;
79 
80 
81  //! Constructor
82  /*!
83  * \param M_ Linear operator ( Read )
84  * \param invParam inverter parameters ( Read )
85  */
86  LinOpSysSolverQPhiXCloverIterRefine(Handle< LinearOperator<T> > A_,
87  Handle< FermState<T,Q,Q> > state_,
88  const SysSolverQPhiXCloverParams& invParam_) :
89  A(A_), invParam(invParam_), clov(new CloverTermT<T, U>()), invclov(new CloverTermT<T, U>())
90  {
91 
92 
93  QDPIO::cout << "LinOpSysSolverQPhiXCloverIterRefine:" << std::endl;
94 
95  if ( toBool( invParam.Delta < 0 ) ) {
96  QDPIO::cout << "Error: Delta parameter for solve should be set" << std::endl;
97  QDP_abort(1);
98  }
99 
100  QDPIO::cout << "AntiPeriodicT is: " << invParam.AntiPeriodicT << std::endl;
101 
102 
103  QDPIO::cout << "Veclen is " << MixedVecTraits<REALT,InnerReal>::Vec << std::endl;
104  QDPIO::cout << "Soalen is " << MixedVecTraits<REALT,InnerReal>::Soa << std::endl;
105  QDPIO::cout << "Inner Veclen is " << MixedVecTraits<REALT,InnerReal>::VecInner << std::endl;
106  QDPIO::cout << "Inner Soalen is " << MixedVecTraits<REALT,InnerReal>::SoaInner << std::endl;
107 
108  if ( MixedVecTraits<REALT,InnerReal>::Soa > MixedVecTraits<REALT,InnerReal>::Vec ) {
109  QDPIO::cerr << "PROBLEM: Soalen > Veclen. Please set soalen appropriately (<=VECLEN) at compile time" << std::endl;
110  QDP_abort(1);
111  }
112 
113  if ( MixedVecTraits<REALT,InnerReal>::SoaInner > MixedVecTraits<REALT,InnerReal>::VecInner ) {
114  QDPIO::cerr << "PROBLEM: Inner Soalen > Inner Veclen. Please set soalen appropriately at compile time" << std::endl;
115  QDP_abort(1);
116  }
117 
118  Q u(Nd);
119  for(int mu=0; mu < Nd; mu++) {
120  u[mu] = state_->getLinks()[mu];
121  }
122 
123  // Set up aniso coefficients
124  multi1d<Real> aniso_coeffs(Nd);
125  for(int mu=0; mu < Nd; mu++) aniso_coeffs[mu] = Real(1);
126 
127  bool anisotropy = invParam.CloverParams.anisoParam.anisoP;
128  if( anisotropy ) {
129  aniso_coeffs = makeFermCoeffs( invParam.CloverParams.anisoParam );
130  }
131 
132  double t_boundary=(double)(1);
133  // NB: In this case, the state will have boundaries applied.
134  // So we only need to apply our own boundaries if Compression is enabled
135  //
136  if (invParam.AntiPeriodicT) {
137  t_boundary=(double)(-1);
138  // Flip off the boundaries -- Dslash expects them off..
139  u[3] *= where(Layout::latticeCoordinate(3) == (Layout::lattSize()[3]-1),
140  Real(t_boundary), Real(1));
141  }
142 
143  cbsize_in_blocks = rb[0].numSiteTable()/MixedVecTraits<REALT,InnerReal>::Soa;
144  const QPhiX::QPhiXCLIArgs& QPhiXParams = TheQPhiXParams::Instance();
145 #ifdef QDP_IS_QDPJIT
146  int pad_xy =0;
147  int pad_xyz = 0;
148 #else
149  int pad_xy = QPhiXParams.getPxy();
150  int pad_xyz = QPhiXParams.getPxyz();
151 #endif
152 
153  // Grab a dslash from which we will get geometry.
154  geom_outer = new QPhiX::Geometry<REALT,
155  MixedVecTraits<REALT,InnerReal>::Vec,
156  MixedVecTraits<REALT,InnerReal>::Soa,
157  MixedVecTraits<REALT,InnerReal>::compress12>(Layout::subgridLattSize().slice(),
158  QPhiXParams.getBy(),
159  QPhiXParams.getBz(),
160  QPhiXParams.getNCores(),
161  QPhiXParams.getSy(),
162  QPhiXParams.getSz(),
163  pad_xy,
164  pad_xyz,
165  QPhiXParams.getMinCt());
166 
167 
168  geom_inner = new QPhiX::Geometry<InnerReal,
169  MixedVecTraits<REALT,InnerReal>::VecInner,
170  MixedVecTraits<REALT,InnerReal>::SoaInner,
171  MixedVecTraits<REALT,InnerReal>::compress12>(Layout::subgridLattSize().slice(),
172  QPhiXParams.getBy(),
173  QPhiXParams.getBz(),
174  QPhiXParams.getNCores(),
175  QPhiXParams.getSy(),
176  QPhiXParams.getSz(),
177 
178  pad_xy,
179  pad_xyz,
180  QPhiXParams.getMinCt());
181 
182 
183 
184 #ifndef QDP_IS_QDPJIT
185  p_even=(QPhiX_Spinor *)geom_outer->allocCBFourSpinor();
186  p_odd=(QPhiX_Spinor *)geom_outer->allocCBFourSpinor();
187  c_even=(QPhiX_Spinor *)geom_outer->allocCBFourSpinor();
188  c_odd=(QPhiX_Spinor *)geom_outer->allocCBFourSpinor();
189  psi_s[0]=p_even;
190  psi_s[1]=p_odd;
191  chi_s[0]=c_even;
192  chi_s[1]=c_odd;
193 #else
194  psi_s[0]=nullptr;
195  psi_s[1]=nullptr;
196  chi_s[0]=nullptr;
197  chi_s[1]=nullptr;
198 #endif
199 
200  // Pack the gauge field
201  QDPIO::cout << "Packing gauge field..." ;
202  // Alloc outer gauge field
203  QPhiX_Gauge* packed_gauge_cb0=(QPhiX_Gauge *)geom_outer->allocCBGauge();
204  QPhiX_Gauge* packed_gauge_cb1=(QPhiX_Gauge *)geom_outer->allocCBGauge();
205 
206  // Alloc inner gauge field
207  QPhiX_InnerGauge* packed_gauge_cb0_inner=(QPhiX_InnerGauge *)geom_inner->allocCBGauge();
208  QPhiX_InnerGauge* packed_gauge_cb1_inner=(QPhiX_InnerGauge *)geom_inner->allocCBGauge();
209 
210  // Pack outer gauge field
211  QPhiX::qdp_pack_gauge<>(u, packed_gauge_cb0,packed_gauge_cb1, *geom_outer);
212  u_packed[0] = packed_gauge_cb0;
213  u_packed[1] = packed_gauge_cb1;
214 
215  // Pack inner gauge field
216  QPhiX::qdp_pack_gauge<>(u, packed_gauge_cb0_inner,packed_gauge_cb1_inner, *geom_inner);
217  u_packed_i[0] = packed_gauge_cb0_inner;
218  u_packed_i[1] = packed_gauge_cb1_inner;
219 
220 
221  QDPIO::cout << "Creating Clover Term" << std::endl;
222  CloverTerm clov_qdp;
223  clov->create(state_, invParam.CloverParams);
224  QDPIO::cout << "Inverting Clover Term" << std::endl;
225  invclov->create(state_, invParam.CloverParams, (*clov));
226  for(int cb=0; cb < 2; cb++) {
227  invclov->choles(cb);
228  }
229  QDPIO::cout << "Done" << std::endl;
230  // Create buffer for outer clover
231  QPhiX_Clover* A_cb0=(QPhiX_Clover *)geom_outer->allocCBClov();
232  QPhiX_Clover* A_cb1=(QPhiX_Clover *)geom_outer->allocCBClov();
233  clov_packed[0] = A_cb0;
234  clov_packed[1] = A_cb1;
235 
236  QPhiX_Clover* A_inv_cb0=(QPhiX_Clover *)geom_outer->allocCBClov();
237  QPhiX_Clover* A_inv_cb1=(QPhiX_Clover *)geom_outer->allocCBClov();
238  invclov_packed[0] = A_inv_cb0;
239  invclov_packed[1] = A_inv_cb1;
240 
241  // Create buffer for inner clover
242  QPhiX_InnerClover* A_cb0_inner=(QPhiX_InnerClover *)geom_inner->allocCBClov();
243  QPhiX_InnerClover* A_cb1_inner=(QPhiX_InnerClover *)geom_inner->allocCBClov();
244  clov_packed_i[0] = A_cb0_inner;
245  clov_packed_i[1] = A_cb1_inner;
246 
247  QPhiX_InnerClover* A_inv_cb0_inner=(QPhiX_InnerClover *)geom_inner->allocCBClov();
248  QPhiX_InnerClover* A_inv_cb1_inner=(QPhiX_InnerClover *)geom_inner->allocCBClov();
249  invclov_packed_i[0] = A_inv_cb0_inner;
250  invclov_packed_i[1] = A_inv_cb1_inner;
251 
252  // Do the actual packing
253  QDPIO::cout << "Packing Clover term..." << std::endl;
254 
255  // Pack outer clover inverse
256  for(int cb=0; cb < 2; cb++) {
257  QPhiX::qdp_pack_clover<>((*invclov), invclov_packed[cb], *geom_outer, cb);
258  }
259 
260  // Pack outer clover
261  for(int cb=0; cb < 2; cb++) {
262  QPhiX::qdp_pack_clover<>((*clov), clov_packed[cb], *geom_outer, cb);
263  }
264 
265  // Pack inner clover inverse
266  for(int cb=0; cb < 2; cb++) {
267  QPhiX::qdp_pack_clover<>((*invclov), invclov_packed_i[cb], *geom_inner, cb);
268  }
269 
270  // Pack inner clover
271  for(int cb=0; cb < 2; cb++) {
272  QPhiX::qdp_pack_clover<>((*clov), clov_packed_i[cb], *geom_inner, cb);
273  }
274  QDPIO::cout << "Done" << std::endl;
275 
276 
277  M_outer=new QPhiX::EvenOddCloverOperator<REALT,
278  MixedVecTraits<REALT,InnerReal>::Vec,
279  MixedVecTraits<REALT,InnerReal>::Soa,
280  MixedVecTraits<REALT,InnerReal>::compress12>(u_packed,
281  clov_packed[1],
282  invclov_packed[0],
283  geom_outer,
284  t_boundary,
285  toDouble(aniso_coeffs[0]),
286  toDouble(aniso_coeffs[3]));
287 
288  M_inner=new QPhiX::EvenOddCloverOperator<InnerReal,
289  MixedVecTraits<REALT,InnerReal>::VecInner,
290  MixedVecTraits<REALT,InnerReal>::SoaInner,
291  MixedVecTraits<REALT,InnerReal>::compress12>(u_packed_i,
292  clov_packed_i[1],
293  invclov_packed_i[0],
294  geom_inner,
295  t_boundary,
296  toDouble(aniso_coeffs[0]),
297  toDouble(aniso_coeffs[3]));
298 
299 
300  bicgstab_inner_solver = new QPhiX::InvBiCGStab<InnerReal,MixedVecTraits<REALT,InnerReal>::VecInner, MixedVecTraits<REALT,InnerReal>::SoaInner, MixedVecTraits<REALT,InnerReal>::compress12>((*M_inner), invParam.MaxIter);
301 
302  mixed_solver = new QPhiX::InvRichardsonMultiPrec<REALT,
303  MixedVecTraits<REALT,InnerReal>::Vec,
304  MixedVecTraits<REALT,InnerReal>::Soa,
305  MixedVecTraits<REALT,InnerReal>::compress12,
306  InnerReal,
307  MixedVecTraits<REALT,InnerReal>::VecInner,
308  MixedVecTraits<REALT,InnerReal>::SoaInner,
309  // NB: &(*) is yukcy.
310  MixedVecTraits<REALT,InnerReal>::compress12>((*M_outer), (*bicgstab_inner_solver), toDouble(invParam.Delta),invParam.MaxIter);
311  }
312 
313 
314 
315  //! Destructor
316  ~LinOpSysSolverQPhiXCloverIterRefine()
317  {
318 
319  // Need to unalloc all the memory...
320 
321  QDPIO::cout << "Destructing" << std::endl;
322 #ifndef QDP_IS_QDPJIT
323  geom_outer->free(p_even);
324  geom_outer->free(p_odd);
325  geom_outer->free(c_even);
326  geom_outer->free(c_odd);
327 #endif
328  psi_s[0] = nullptr;
329  psi_s[1] = nullptr;
330  chi_s[0] = nullptr;
331  chi_s[1] = nullptr;
332 
333  geom_outer->free(invclov_packed[0]);
334  geom_outer->free(invclov_packed[1]);
335  invclov_packed[0] = nullptr;
336  invclov_packed[1] = nullptr;
337 
338 
339  geom_outer->free(clov_packed[0]);
340  geom_outer->free(clov_packed[1]);
341  clov_packed[0] = nullptr;
342  clov_packed[1] = nullptr;
343 
344 
345  geom_outer->free(u_packed[0]);
346  geom_outer->free(u_packed[1]);
347  u_packed[0] = nullptr;
348  u_packed[1] = nullptr;
349 
350 
351  geom_inner->free(invclov_packed_i[0]);
352  geom_inner->free(invclov_packed_i[1]);
353  invclov_packed_i[0] = nullptr;
354  invclov_packed_i[1] = nullptr;
355 
356 
357  geom_inner->free(clov_packed_i[0]);
358  geom_inner->free(clov_packed_i[1]);
359  clov_packed_i[0] = nullptr;
360  clov_packed_i[1] = nullptr;
361 
362 
363  geom_inner->free(u_packed_i[0]);
364  geom_inner->free(u_packed_i[1]);
365  u_packed_i[0] = nullptr;
366  u_packed_i[1] = nullptr;
367 
368  delete geom_inner;
369  delete geom_outer;
370 
371 
372  }
373 
374  //! Return the subset on which the operator acts
375  const Subset& subset() const {return A->subset();}
376 
377 
378 
379  //! Solver the linear system
380  /*!
381  * \param psi solution ( Modify )
382  * \param chi source ( Read )
383  * \return syssolver results
384  */
385  SystemSolverResults_t operator()(T& psi, const T& chi) const
386  {
387  /* Factories here later? */
388  SystemSolverResults_t res;
389 #ifndef QDP_IS_QDPJIT
390  QDPIO::cout << "Packing Spinors" << std::endl << std::flush ;
391 
392  QPhiX::qdp_pack_cb_spinor<>(psi,psi_s[1], (*M_outer).getGeometry(),1);
393  QPhiX::qdp_pack_cb_spinor<>(chi,chi_s[1], (*M_outer).getGeometry(),1);
394 
395  QDPIO::cout << "Done" << std::endl << std::flush;
396 #else
397  psi_s[0] = (QPhiX_Spinor *)( psi.getFjit());
398  psi_s[1] = (QPhiX_Spinor *)( psi.getFjit() ) + cbsize_in_blocks;
399  chi_s[0] = (QPhiX_Spinor *)( chi.getFjit());
400  chi_s[1] = (QPhiX_Spinor *)( chi.getFjit() ) + cbsize_in_blocks;
401 #endif
402  double rsd_final;
403  unsigned long site_flops=0;
404  unsigned long mv_apps=0;
405 
406  QDPIO::cout << "Starting solve" << std::endl << std::flush ;
407  double start = omp_get_wtime();
408  int my_isign=1;
409  (*mixed_solver)(psi_s[1],chi_s[1], toDouble(invParam.RsdTarget), res.n_count, rsd_final, site_flops, mv_apps, my_isign,invParam.VerboseP);
410  double end = omp_get_wtime();
411 
412  QDPIO::cout << "QPHIX_CLOVER_BICGSTAB_ITER_REFINE_SOLVER: " << res.n_count << " iters, rsd_sq_final=" << rsd_final << std::endl;
413 #ifndef QDP_IS_QDPJIT
414  QPhiX::qdp_unpack_cb_spinor<>(psi_s[1], psi, (*M_outer).getGeometry(),1);
415 #endif
416 
417 #if 1
418  // Chi Should now hold the result spinor
419  // Check it against chroma.
420  T r = chi;
421  T tmp;
422  (*A)(tmp, psi, PLUS);
423  r[ A->subset() ] -= tmp;
424 
425  Double r2 = norm2(r,A->subset());
426  Double b2 = norm2(chi, A->subset());
427  Double rel_resid = sqrt(r2/b2);
428 
429  QDPIO::cout << "QPHIX_CLOVER_BICGSTAB_ITER_REFINE_SOLVER: || r || / || b || = " << rel_resid << std::endl;
430 #endif
431 #if 0
432  if ( !toBool ( rel_resid < invParam.RsdTarget*invParam.RsdToleranceFactor ) ) {
433  QDPIO::cout << "SOLVE FAILED" << std::endl;
434  QDP_abort(1);
435  }
436 #endif
437 
438  int num_cb_sites = Layout::vol()/2;
439  unsigned long total_flops = (site_flops + (1320+504+1320+504+48)*mv_apps)*num_cb_sites;
440  double gflops = (double)(total_flops)/(1.0e9);
441 
442  double total_time = end - start;
443  QDPIO::cout << "QPHIX_CLOVER_BICGSTAB_ITER_REFINE_SOLVER: Solver Time="<< total_time <<" (sec) Performance=" << gflops / total_time << " GFLOPS" << std::endl;
444 
445  END_CODE();
446  return res;
447  }
448 
449 
450 
451 
452 
453  private:
454  // Hide default constructor
455  LinOpSysSolverQPhiXCloverIterRefine() {}
456 
457 
458 
459  Handle< LinearOperator<T> > A;
460  const SysSolverQPhiXCloverParams invParam;
461  Handle< CloverTermT<T, U> > clov;
462  Handle< CloverTermT<T, U> > invclov;
463 
464  QPhiX::Geometry<REALT,
465  MixedVecTraits<REALT,InnerReal>::Vec,
466  MixedVecTraits<REALT,InnerReal>::Soa,
467  MixedVecTraits<REALT,InnerReal>::compress12>* geom_outer;
468 
469  QPhiX::Geometry<InnerReal,
470  MixedVecTraits<REALT,InnerReal>::VecInner,
471  MixedVecTraits<REALT,InnerReal>::SoaInner,
472  MixedVecTraits<REALT,InnerReal>::compress12>* geom_inner;
473 
474  Handle< QPhiX::EvenOddCloverOperator<REALT,
475  MixedVecTraits<REALT,InnerReal>::Vec,
476  MixedVecTraits<REALT,InnerReal>::Soa,
477  MixedVecTraits<REALT,InnerReal>::compress12> > M_outer;
478 
479 Handle< QPhiX::EvenOddCloverOperator<InnerReal,
480  MixedVecTraits<REALT,InnerReal>::VecInner,
481  MixedVecTraits<REALT,InnerReal>::SoaInner,
482  MixedVecTraits<REALT,InnerReal>::compress12> > M_inner;
483 
484 
485  Handle< QPhiX::InvBiCGStab<InnerReal,
486  MixedVecTraits<REALT,InnerReal>::VecInner,
487  MixedVecTraits<REALT,InnerReal>::SoaInner,
488  MixedVecTraits<REALT,InnerReal>::compress12> > bicgstab_inner_solver;
489 
490 Handle< QPhiX::InvRichardsonMultiPrec<REALT,
491  MixedVecTraits<REALT,InnerReal>::Vec,
492  MixedVecTraits<REALT,InnerReal>::Soa,
493  MixedVecTraits<REALT,InnerReal>::compress12,
494  InnerReal,
495  MixedVecTraits<REALT,InnerReal>::VecInner,
496  MixedVecTraits<REALT,InnerReal>::SoaInner,
497  MixedVecTraits<REALT,InnerReal>::compress12> > mixed_solver;
498 
499  QPhiX_Clover* invclov_packed[2];
500  QPhiX_Clover* clov_packed[2];
501  QPhiX_Gauge* u_packed[2];
502 
503  QPhiX_InnerClover* invclov_packed_i[2];
504  QPhiX_InnerClover* clov_packed_i[2];
505  QPhiX_InnerGauge* u_packed_i[2];
506 
507  QPhiX_Spinor* p_even;
508  QPhiX_Spinor* p_odd;
509  QPhiX_Spinor* c_even;
510  QPhiX_Spinor* c_odd;
511  mutable QPhiX_Spinor* psi_s[2];
512  mutable QPhiX_Spinor* chi_s[2];
513  size_t cbsize_in_blocks;
514  };
515 
516 
517 } // End namespace
518 
519 #endif // BUILD_QPHIX
520 #endif
521 
Anisotropy parameters.
void create(Handle< FermState< T, multi1d< U >, multi1d< U > > > fs, const CloverFermActParams &param_)
Creation routine.
Include possibly optimized Clover terms.
int mu
Definition: cool.cc:24
Even-odd preconditioned Clover fermion linear operator.
Class for counted reference semantics.
Linear Operators.
Nd
Definition: meslate.cc:74
multi1d< Hadron2PtContraction_t > operator()(const multi1d< LatticeColorMatrix > &u)
Asqtad Staggered-Dirac operator.
Definition: klein_gord.cc:10
static multi1d< LatticeColorMatrix > u
LinOpSysSolverMGProtoClover::Q Q
QDPCloverTerm CloverTerm
Definition: clover_term_w.h:92
LatticeFermion tmp
Definition: mespbg5p_w.cc:36
SystemSolverResults_t InvBiCGStab(const LinearOperator< LatticeFermionF > &A, const LatticeFermionF &chi, LatticeFermionF &psi, const Real &RsdBiCGStab, int MaxBiCGStab, enum PlusMinus isign)
Definition: invbicgstab.cc:222
@ PLUS
Definition: chromabase.h:45
multi1d< LatticeFermion > chi(Ncb)
LatticeFermion psi
Definition: mespbg5p_w.cc:35
QDPCloverTermT< T, U > CloverTermT
Definition: clover_term_w.h:97
A(A, psi, r, Ncb, PLUS)
int cb
Definition: invbicg.cc:120
multi1d< Real > makeFermCoeffs(const AnisoParam_t &aniso)
Make fermion coefficients.
Definition: aniso_io.cc:63
FloatingPoint< double > Double
Definition: gtest.h:7351
Periodic ferm state and a creator.
Fermion action factories.
Reunitarize in place a color matrix to SU(N)
Simple fermionic BC.
Support class for fermion actions and linear operators.
Linear system solvers.
multi1d< LatticeColorMatrix > U
LatticeFermion T
Definition: t_clover.cc:11
multi1d< LatticeColorMatrix > Q
Definition: t_clover.cc:12
Axial gauge fixing.