QUDA  0.9.0
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 <dirac_quda.h>
15 #include <vector>
16 
17 namespace quda {
18 
30  class Transfer {
31 
32  private:
33 
35  const std::vector<ColorSpinorField*> &B;
36 
38  const int Nvec;
39 
42 
45 
48 
51 
54 
57 
59  int *geo_bs;
60 
66 
71 
77 
82 
84  int spin_bs;
85 
87  int *spin_map;
88 
91 
94 
96  bool enable_gpu;
97 
100  mutable bool use_gpu;
101 
105  void fillV(ColorSpinorField &V);
106 
111  void createGeoMap(int *geo_bs);
112 
119  //template <class Complex, class FieldOrder>
120  //void blockOrderV(Complex *out, const FieldOrder &in, int *geo_bs, int spin_bs);
121 
128  //template <class FieldOrder, class Complex>
129  //void undoblockOrderV(FieldOrder &out, Complex *in, int *geo_bs, int spin_bs);
130 
138  //template <class Complex>
139  //void blockGramSchmidt(Complex *v, int nBlocks, int Nc, int blockSize);
140 
145  void createSpinMap(int spin_bs);
146 
150  mutable double flops_;
151 
157 
158  public:
159 
170  Transfer(const std::vector<ColorSpinorField*> &B, int Nvec, int *geo_bs, int spin_bs,
172 
174  virtual ~Transfer();
175 
181  void P(ColorSpinorField &out, const ColorSpinorField &in) const;
182 
188  void R(ColorSpinorField &out, const ColorSpinorField &in) const;
189 
196  return (location == QUDA_CUDA_FIELD_LOCATION) ? *V_d : *V_h;
197  }
198 
203  int nvec() const {return Nvec;}
204 
209  int Spin_bs() const {return spin_bs;}
210 
215  const int *Geo_bs() const {return geo_bs;}
216 
221  void setTransferGPU(bool use_gpu) const { this->use_gpu = use_gpu; }
222 
239 
244  double flops() const;
245  };
246 
253  void FillV(ColorSpinorField &V, const std::vector<ColorSpinorField*> &B, int Nvec);
254 
266  void BlockOrthogonalize(ColorSpinorField &V, int Nvec, const int *geo_bs,
267  const int *fine_to_coarse, int spin_bs);
268 
279  void Prolongate(ColorSpinorField &out, const ColorSpinorField &in, const ColorSpinorField &v,
280  int Nvec, const int *fine_to_coarse, const int *spin_map,
282 
293  void Restrict(ColorSpinorField &out, const ColorSpinorField &in, const ColorSpinorField &v,
294  int Nvec, const int *fine_to_coarse, const int *coarse_to_fine, const int *spin_map,
296 
297 
298 } // namespace quda
299 #endif // _TRANSFER_H
bool enable_gpu
Definition: transfer.h:96
TimeProfile & profile
Definition: transfer.h:156
const int Nvec
Definition: transfer.h:38
void R(ColorSpinorField &out, const ColorSpinorField &in) const
Definition: transfer.cpp:300
int * fine_to_coarse_d
Definition: transfer.h:76
void P(ColorSpinorField &out, const ColorSpinorField &in) const
Definition: transfer.cpp:262
ColorSpinorField * fine_tmp_h
Definition: transfer.h:47
double flops_
Definition: transfer.h:150
int nvec() const
Definition: transfer.h:203
void setSiteSubset(QudaSiteSubset site_subset, QudaParity parity)
Sets whether the transfer operator is to act on full fields or single parity fields, and if single-parity which parity. If site_subset is QUDA_FULL_SITE_SUBSET, the transfer operator can still be applied to single-parity fields, however, if site_subset is QUDA_PARITY_SITE_SUBSET, then the transfer operator cannot be applied to full fields, and setSiteSubset will need to be called first to reset to QUDA_FULL_SITE_SUBSET. This method exists to reduce GPU memory overhead - if only transfering single-parity fine fields then we only store a single-parity copy of the null space components on the device.
Definition: transfer.cpp:155
QudaParity parity
Definition: transfer.h:93
ColorSpinorField * coarse_tmp_d
Definition: transfer.h:56
void FillV(ColorSpinorField &V, const std::vector< ColorSpinorField *> &B, int Nvec)
int * coarse_to_fine_h
Definition: transfer.h:70
const ColorSpinorField & Vectors(QudaFieldLocation location=QUDA_CPU_FIELD_LOCATION) const
Definition: transfer.h:195
QudaSiteSubset site_subset
Definition: transfer.h:90
void createGeoMap(int *geo_bs)
Definition: transfer.cpp:212
cpuColorSpinorField * in
ColorSpinorField * coarse_tmp_h
Definition: transfer.h:53
int V
Definition: test_util.cpp:28
void Restrict(ColorSpinorField &out, const ColorSpinorField &in, const ColorSpinorField &v, int Nvec, const int *fine_to_coarse, const int *coarse_to_fine, const int *spin_map, int parity=QUDA_INVALID_PARITY)
Apply the restriction operator.
Definition: restrictor.cu:509
int Spin_bs() const
Definition: transfer.h:209
int * spin_map
Definition: transfer.h:87
ColorSpinorField * V_h
Definition: transfer.h:41
enum QudaParity_s QudaParity
const int * Geo_bs() const
Definition: transfer.h:215
enum QudaSiteSubset_s QudaSiteSubset
int * fine_to_coarse_h
Definition: transfer.h:65
enum QudaFieldLocation_s QudaFieldLocation
cpuColorSpinorField * out
void fillV(ColorSpinorField &V)
Definition: transfer.cpp:197
Transfer(const std::vector< ColorSpinorField *> &B, int Nvec, int *geo_bs, int spin_bs, bool enable_gpu, TimeProfile &profile)
Definition: transfer.cpp:21
ColorSpinorField * V_d
Definition: transfer.h:44
void setTransferGPU(bool use_gpu) const
Definition: transfer.h:221
ColorSpinorField * fine_tmp_d
Definition: transfer.h:50
void createSpinMap(int spin_bs)
Definition: transfer.cpp:253
void BlockOrthogonalize(ColorSpinorField &V, int Nvec, const int *geo_bs, const int *fine_to_coarse, int spin_bs)
Block orthogonnalize the matrix field, where the blocks are defined by lookup tables that map the fin...
const std::vector< ColorSpinorField * > & B
Definition: transfer.h:35
QudaParity parity
Definition: covdev_test.cpp:53
double flops() const
Definition: transfer.cpp:342
int * coarse_to_fine_d
Definition: transfer.h:81
int * geo_bs
Definition: transfer.h:59
virtual ~Transfer()
Definition: transfer.cpp:137
void Prolongate(ColorSpinorField &out, const ColorSpinorField &in, const ColorSpinorField &v, int Nvec, const int *fine_to_coarse, const int *spin_map, int parity=QUDA_INVALID_PARITY)
Apply the prolongation operator.
Definition: prolongator.cu:284