QUDA  1.0.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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 
44 
47 
50 
53 
56 
59 
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 
99 
102 
104  mutable bool enable_gpu;
105 
107  mutable bool enable_cpu;
108 
111  mutable bool use_gpu;
112 
117  void createV(QudaFieldLocation location) const;
118 
123  void createTmp(QudaFieldLocation location) const;
124 
129  void createGeoMap(int *geo_bs);
130 
135  void createSpinMap(int spin_bs);
136 
141  void initializeLazy(QudaFieldLocation location) const;
142 
146  mutable double flops_;
147 
153 
154  public:
167  Transfer(const std::vector<ColorSpinorField *> &B, int Nvec, int NblockOrtho, int *geo_bs, int spin_bs,
168  QudaPrecision null_precision, TimeProfile &profile);
169 
171  virtual ~Transfer();
172 
176  void reset();
177 
183  void P(ColorSpinorField &out, const ColorSpinorField &in) const;
184 
190  void R(ColorSpinorField &out, const ColorSpinorField &in) const;
191 
196  {
197  return location == QUDA_CUDA_FIELD_LOCATION ? null_precision : std::max(B[0]->Precision(), QUDA_SINGLE_PRECISION);
198  }
199 
206  if (location == QUDA_INVALID_FIELD_LOCATION) {
207  // if not set then we return the memory space where the input vectors are stored
208  return B[0]->Location() == QUDA_CUDA_FIELD_LOCATION ? *V_d : *V_h;
209  } else {
210  return location == QUDA_CUDA_FIELD_LOCATION ? *V_d : *V_h;
211  }
212  }
213 
218  int nvec() const {return Nvec;}
219 
224  int Spin_bs() const {return spin_bs;}
225 
230  const int *Geo_bs() const {return geo_bs;}
231 
236  { return location == QUDA_CPU_FIELD_LOCATION ? fine_to_coarse_h : fine_to_coarse_d; }
237 
242  { return location == QUDA_CPU_FIELD_LOCATION ? coarse_to_fine_h : coarse_to_fine_d; }
243 
248  void setTransferGPU(bool use_gpu) const { this->use_gpu = use_gpu; }
249 
258  void setSiteSubset(QudaSiteSubset site_subset, QudaParity parity);
259 
264  double flops() const;
265  };
266 
280  void BlockOrthogonalize(ColorSpinorField &V, const std::vector<ColorSpinorField *> &B, const int *fine_to_coarse,
281  const int *coarse_to_fine, const int *geo_bs, const int spin_bs, const int n_block_ortho);
282 
294  int Nvec, const int *fine_to_coarse, const int * const *spin_map,
296 
307  void Restrict(ColorSpinorField &out, const ColorSpinorField &in, const ColorSpinorField &v,
308  int Nvec, const int *fine_to_coarse, const int *coarse_to_fine, const int * const *spin_map,
310 
311 
312 } // namespace quda
313 #endif // _TRANSFER_H
bool enable_gpu
Definition: transfer.h:104
TimeProfile & profile
Definition: transfer.h:152
void createTmp(QudaFieldLocation location) const
Allocate temporaries used when applying transfer operators.
Definition: transfer.cpp:134
const int Nvec
Definition: transfer.h:37
Transfer(const std::vector< ColorSpinorField *> &B, int Nvec, int NblockOrtho, int *geo_bs, int spin_bs, QudaPrecision null_precision, TimeProfile &profile)
Definition: transfer.cpp:19
enum QudaPrecision_s QudaPrecision
void initializeLazy(QudaFieldLocation location) const
Lazy allocation of the transfer operator in a given location.
Definition: transfer.cpp:154
const QudaPrecision null_precision
Definition: transfer.h:43
void R(ColorSpinorField &out, const ColorSpinorField &in) const
Definition: transfer.cpp:344
bool enable_cpu
Definition: transfer.h:107
int * fine_to_coarse_d
Definition: transfer.h:81
void P(ColorSpinorField &out, const ColorSpinorField &in) const
Definition: transfer.cpp:305
ColorSpinorField * fine_tmp_h
Definition: transfer.h:52
double flops_
Definition: transfer.h:146
int nvec() const
Definition: transfer.h:218
void createV(QudaFieldLocation location) const
Allocate V field.
Definition: transfer.cpp:101
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.
Definition: transfer.cpp:227
const int * fineToCoarse(QudaFieldLocation location=QUDA_CPU_FIELD_LOCATION) const
Definition: transfer.h:235
QudaParity parity
Definition: transfer.h:101
ColorSpinorField * coarse_tmp_d
Definition: transfer.h:61
const int * coarseToFine(QudaFieldLocation location=QUDA_CPU_FIELD_LOCATION) const
Definition: transfer.h:241
int * coarse_to_fine_h
Definition: transfer.h:75
QudaSiteSubset site_subset
Definition: transfer.h:98
void createGeoMap(int *geo_bs)
Creates the map between fine and coarse grids.
Definition: transfer.cpp:248
cpuColorSpinorField * in
ColorSpinorField * coarse_tmp_h
Definition: transfer.h:58
const int NblockOrtho
Definition: transfer.h:40
int Spin_bs() const
Definition: transfer.h:224
ColorSpinorField * V_h
Definition: transfer.h:46
enum QudaParity_s QudaParity
int V
Definition: test_util.cpp:27
void Restrict(Arg arg)
Definition: restrictor.cuh:90
const int * Geo_bs() const
Definition: transfer.h:230
int ** spin_map
Definition: transfer.h:92
enum QudaSiteSubset_s QudaSiteSubset
int * fine_to_coarse_h
Definition: transfer.h:70
QudaPrecision NullPrecision(QudaFieldLocation location) const
The precision of the packed null-space vectors.
Definition: transfer.h:195
enum QudaFieldLocation_s QudaFieldLocation
cpuColorSpinorField * out
const int nspin_fine
Definition: transfer.h:95
ColorSpinorField * V_d
Definition: transfer.h:49
void setTransferGPU(bool use_gpu) const
Definition: transfer.h:248
ColorSpinorField * fine_tmp_d
Definition: transfer.h:55
void createSpinMap(int spin_bs)
Creates the map between fine spin and parity to coarse spin dimensions.
Definition: transfer.cpp:289
const std::vector< ColorSpinorField * > & B
Definition: transfer.h:34
int n_block_ortho[QUDA_MAX_MG_LEVEL]
Definition: test_util.cpp:1673
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.
Definition: prolongator.cu:296
double flops() const
Definition: transfer.cpp:387
const ColorSpinorField & Vectors(QudaFieldLocation location=QUDA_INVALID_FIELD_LOCATION) const
Definition: transfer.h:205
int * coarse_to_fine_d
Definition: transfer.h:86
int * geo_bs
Definition: transfer.h:64
virtual ~Transfer()
Definition: transfer.cpp:205
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 reset()
for resetting the Transfer when the null vectors have changed
Definition: transfer.cpp:182