QUDA  v1.1.0
A library for QCD on GPUs
transfer.h
Go to the documentation of this file.
1 #ifndef _TRANSFER_H
2 #define _TRANSFER_H
3 
13 #include <color_spinor_field.h>
14 #include <vector>
15 
16 namespace quda {
17 
29  class Transfer {
30 
31  private:
32 
34  const std::vector<ColorSpinorField*> &B;
35 
37  const int Nvec;
38 
40  const int NblockOrtho;
41 
43  const QudaPrecision null_precision;
44 
46  mutable ColorSpinorField *V_h;
47 
49  mutable ColorSpinorField *V_d;
50 
52  mutable ColorSpinorField *fine_tmp_h;
53 
55  mutable ColorSpinorField *fine_tmp_d;
56 
58  mutable ColorSpinorField *coarse_tmp_h;
59 
61  mutable ColorSpinorField *coarse_tmp_d;
62 
64  int *geo_bs;
65 
70  mutable int *fine_to_coarse_h;
71 
75  mutable int *coarse_to_fine_h;
76 
81  mutable int *fine_to_coarse_d;
82 
86  mutable int *coarse_to_fine_d;
87 
89  int spin_bs;
90 
92  int **spin_map;
93 
95  const int nspin_fine;
96 
98  QudaSiteSubset site_subset;
99 
101  QudaParity parity;
102 
104  mutable bool enable_gpu;
105 
107  mutable bool enable_cpu;
108 
111  mutable bool use_gpu;
112 
115  mutable QudaTransferType transfer_type;
116 
121  void createV(QudaFieldLocation location) const;
122 
127  void createTmp(QudaFieldLocation location) const;
128 
133  void createGeoMap(int *geo_bs);
134 
139  void createSpinMap(int spin_bs);
140 
145  void initializeLazy(QudaFieldLocation location) const;
146 
150  mutable double flops_;
151 
156  TimeProfile &profile;
157 
158  public:
171  Transfer(const std::vector<ColorSpinorField *> &B, int Nvec, int NblockOrtho, int *geo_bs, int spin_bs,
172  QudaPrecision null_precision, const QudaTransferType transfer_type, TimeProfile &profile);
173 
175  virtual ~Transfer();
176 
180  void reset();
181 
187  void P(ColorSpinorField &out, const ColorSpinorField &in) const;
188 
194  void R(ColorSpinorField &out, const ColorSpinorField &in) const;
195 
200  {
201  return location == QUDA_CUDA_FIELD_LOCATION ? null_precision : std::max(B[0]->Precision(), QUDA_SINGLE_PRECISION);
202  }
203 
210  if (location == QUDA_INVALID_FIELD_LOCATION) {
211  // if not set then we return the memory space where the input vectors are stored
212  return B[0]->Location() == QUDA_CUDA_FIELD_LOCATION ? *V_d : *V_h;
213  } else {
214  return location == QUDA_CUDA_FIELD_LOCATION ? *V_d : *V_h;
215  }
216  }
217 
222  int nvec() const {return Nvec;}
223 
228  int Spin_bs() const {return spin_bs;}
229 
234  const int *Geo_bs() const {return geo_bs;}
235 
240  QudaTransferType getTransferType() const { return transfer_type; }
241 
246  { return location == QUDA_CPU_FIELD_LOCATION ? fine_to_coarse_h : fine_to_coarse_d; }
247 
252  { return location == QUDA_CPU_FIELD_LOCATION ? coarse_to_fine_h : coarse_to_fine_d; }
253 
258  void setTransferGPU(bool use_gpu) const { this->use_gpu = use_gpu; }
259 
268  void setSiteSubset(QudaSiteSubset site_subset, QudaParity parity);
269 
274  double flops() const;
275  };
276 
290  void BlockOrthogonalize(ColorSpinorField &V, const std::vector<ColorSpinorField *> &B, const int *fine_to_coarse,
291  const int *coarse_to_fine, const int *geo_bs, const int spin_bs, const int n_block_ortho);
292 
304  int Nvec, const int *fine_to_coarse, const int * const *spin_map,
306 
318  int Nvec, const int *fine_to_coarse, const int *coarse_to_fine, const int * const *spin_map,
320 
329  void StaggeredProlongate(ColorSpinorField &out, const ColorSpinorField &in, const int *fine_to_coarse,
330  const int *const *spin_map, int parity = QUDA_INVALID_PARITY);
331 
340  void StaggeredRestrict(ColorSpinorField &out, const ColorSpinorField &in, const int *fine_to_coarse,
341  const int *const *spin_map, int parity = QUDA_INVALID_PARITY);
342 
343 } // namespace quda
344 #endif // _TRANSFER_H
QudaFieldLocation Location() const
int nvec() const
Definition: transfer.h:222
void reset()
for resetting the Transfer when the null vectors have changed
Definition: transfer.cpp:226
virtual ~Transfer()
Definition: transfer.cpp:251
void R(ColorSpinorField &out, const ColorSpinorField &in) const
Definition: transfer.cpp:412
Transfer(const std::vector< ColorSpinorField * > &B, int Nvec, int NblockOrtho, int *geo_bs, int spin_bs, QudaPrecision null_precision, const QudaTransferType transfer_type, TimeProfile &profile)
Definition: transfer.cpp:21
const int * Geo_bs() const
Definition: transfer.h:234
void setSiteSubset(QudaSiteSubset site_subset, QudaParity parity)
Sets whether the transfer operator is to act on full fields or single parity fields,...
Definition: transfer.cpp:273
double flops() const
Definition: transfer.cpp:478
void setTransferGPU(bool use_gpu) const
Definition: transfer.h:258
const int * coarseToFine(QudaFieldLocation location=QUDA_CPU_FIELD_LOCATION) const
Definition: transfer.h:251
QudaPrecision NullPrecision(QudaFieldLocation location) const
The precision of the packed null-space vectors.
Definition: transfer.h:199
const int * fineToCoarse(QudaFieldLocation location=QUDA_CPU_FIELD_LOCATION) const
Definition: transfer.h:245
QudaTransferType getTransferType() const
Definition: transfer.h:240
int Spin_bs() const
Definition: transfer.h:228
void P(ColorSpinorField &out, const ColorSpinorField &in) const
Definition: transfer.cpp:350
const ColorSpinorField & Vectors(QudaFieldLocation location=QUDA_INVALID_FIELD_LOCATION) const
Definition: transfer.h:209
quda::mgarray< int > n_block_ortho
int V
Definition: host_utils.cpp:37
QudaParity parity
Definition: covdev_test.cpp:40
enum QudaPrecision_s QudaPrecision
@ QUDA_CUDA_FIELD_LOCATION
Definition: enum_quda.h:326
@ QUDA_CPU_FIELD_LOCATION
Definition: enum_quda.h:325
@ QUDA_INVALID_FIELD_LOCATION
Definition: enum_quda.h:327
enum QudaTransferType_s QudaTransferType
@ QUDA_INVALID_PARITY
Definition: enum_quda.h:284
enum QudaSiteSubset_s QudaSiteSubset
enum QudaFieldLocation_s QudaFieldLocation
@ QUDA_SINGLE_PRECISION
Definition: enum_quda.h:64
enum QudaParity_s QudaParity
void StaggeredRestrict(ColorSpinorField &out, const ColorSpinorField &in, const int *fine_to_coarse, const int *const *spin_map, int parity=QUDA_INVALID_PARITY)
Apply the unitary "restriction" operator for Kahler-Dirac preconditioning.
void Restrict(ColorSpinorField &out, const ColorSpinorField &in, const ColorSpinorField &v, int Nvec, const int *fine_to_coarse, const int *coarse_to_fine, const int *const *spin_map, int parity=QUDA_INVALID_PARITY)
Apply the restriction operator.
void BlockOrthogonalize(ColorSpinorField &V, const std::vector< ColorSpinorField * > &B, const int *fine_to_coarse, const int *coarse_to_fine, const int *geo_bs, const int spin_bs, const int n_block_ortho)
Block orthogonnalize the matrix field, where the blocks are defined by lookup tables that map the fin...
void StaggeredProlongate(ColorSpinorField &out, const ColorSpinorField &in, const int *fine_to_coarse, const int *const *spin_map, int parity=QUDA_INVALID_PARITY)
Apply the unitary "prolongation" operator for Kahler-Dirac preconditioning.
void Prolongate(ColorSpinorField &out, const ColorSpinorField &in, const ColorSpinorField &v, int Nvec, const int *fine_to_coarse, const int *const *spin_map, int parity=QUDA_INVALID_PARITY)
Apply the prolongation operator.