00001 #ifndef _fftlib_h_
00002 #define _fftlib_h_
00003
00004 #include <charm++.h>
00005 #include "ckcomplex.h"
00006 #include "rfftw.h"
00007 #include "EachToManyMulticastStrategy.h"
00008 #include "StreamingStrategy.h"
00009 #include "comlib.h"
00010
00011 #define COMPLEX_TO_REAL -11
00012 #define REAL_TO_COMPLEX -12
00013 #define COMPLEX_TO_COMPLEX -13
00014 #define NULL_TO_NULL -14
00015
00016 #define MAX_FFTS 5
00017
00018 class NormalFFTinfo {
00019 public:
00020 NormalFFTinfo(int sDim[2], int dDim[2], int isSrc,
00021 void *dptr, int transType,
00022 int sPlanesPerSlab=1, int dPlanesPerSlab=1) {
00023 init(sDim, dDim, isSrc, dptr, transType, sPlanesPerSlab, dPlanesPerSlab);
00024 }
00025 NormalFFTinfo(NormalFFTinfo &info) {
00026 init(info.srcSize, info.destSize, info.isSrcSlab,
00027 info.dataPtr, info.transformType,
00028 info.srcPlanesPerSlab, info.destPlanesPerSlab);
00029 }
00030
00031 NormalFFTinfo(void) {dataPtr=NULL;}
00032
00033 int srcSize[2], destSize[2];
00034 bool isSrcSlab;
00035 int srcPlanesPerSlab, destPlanesPerSlab;
00036 void *dataPtr;
00037 int transformType;
00038
00039 void pup(PUP::er &p) {
00040 p(srcSize, 2);
00041 p(destSize, 2);
00042 p(isSrcSlab);
00043 p(srcPlanesPerSlab);
00044 p(destPlanesPerSlab);
00045 p|transformType;
00046 if (p.isUnpacking())
00047 dataPtr = NULL;
00048 }
00049
00050 private:
00051 void init( int sDim[2], int dDim[2], int isSrc,
00052 void *dptr, int transT,
00053 int sPlanesPerSlab, int dPlanesPerSlab) {
00054 if (sDim[1] != dDim[1])
00055 ckerr << "WARNING"
00056 << "This configuration of the source and destination "
00057 << "is not consistent, check the dimensions. The program is "
00058 << "likely to misbehave"
00059 << endl;
00060 isSrcSlab = isSrc;
00061 srcPlanesPerSlab = sPlanesPerSlab;
00062 destPlanesPerSlab = dPlanesPerSlab;
00063 dataPtr = dptr;
00064 transformType=transT;
00065 memcpy(srcSize, sDim, 2 * sizeof(int));
00066 memcpy(destSize, dDim, 2 * sizeof(int));
00067 }
00068 };
00069
00070 typedef struct _PencilType{
00071 const static int XLINE = 0;
00072 const static int YLINE = 1;
00073 const static int ZLINE = 2;
00074 }PencilType;
00075
00076 typedef struct _PencilBlock{
00077 const static int FLATBLOCK = 0;
00078 const static int SQUAREBLOCK = 1;
00079 }PencilBlock;
00080
00081
00082 class LineFFTinfo {
00083 public:
00084
00085 LineFFTinfo(int size[3], int _ptype, int _pblock, complex *dptr, int _xPencilsPerSlab=1, int _yPencilsPerSlab=1, int _zPencilsPerSlab=1) {
00086 init(size[0], size[1], size[2], _ptype, _pblock, dptr, _xPencilsPerSlab, _yPencilsPerSlab, _zPencilsPerSlab);
00087 }
00088 LineFFTinfo(LineFFTinfo &info) {
00089 init(info.sizeX, info.sizeY, info.sizeZ, info.ptype, info.pblock, (complex *) NULL, info.xPencilsPerSlab, info.yPencilsPerSlab, info.zPencilsPerSlab);
00090 }
00091 LineFFTinfo(void) {}
00092
00093
00094 void pup(PUP::er &p) {
00095 p|sizeX;
00096 p|sizeY;
00097 p|sizeZ;
00098 p(ptype);
00099 p(pblock);
00100 p(xPencilsPerSlab);
00101 p(yPencilsPerSlab);
00102 p(zPencilsPerSlab);
00103 p(xsquare, 2);
00104 p(ysquare, 2);
00105 p(zsquare, 2);
00106 if (p.isUnpacking())
00107 dataPtr = (complex *) NULL;
00108 }
00109 int sizeX, sizeY, sizeZ;
00110 int ptype;
00111 int pblock;
00112 int xPencilsPerSlab, yPencilsPerSlab, zPencilsPerSlab;
00113 int xsquare[2], ysquare[2], zsquare[2];
00114 complex *dataPtr;
00115 private:
00116 void init(int sizex, int sizey, int sizez, int _ptype, int _pblock, complex *dptr, int _xPencilsPerSlab, int _yPencilsPerSlab, int _zPencilsPerSlab) {
00117 if (sizex != sizey || sizey != sizez)
00118 ckerr << "WARNING"
00119 << "This configuration of the source and destination "
00120 << "is not consistent, check the dimensions. The program is "
00121 << "likely to misbehave"
00122 << endl;
00123 ptype = _ptype;
00124 pblock = _pblock;
00125 xPencilsPerSlab = _xPencilsPerSlab;
00126 yPencilsPerSlab = _yPencilsPerSlab;
00127 zPencilsPerSlab = _zPencilsPerSlab;
00128 dataPtr = dptr;
00129 sizeX = sizex;
00130 sizeY = sizey;
00131 sizeZ = sizez;
00132
00133 CkAssert((sizeX==sizeY) && (sizeX==sizeZ));
00134
00135 if(pblock == PencilBlock::SQUAREBLOCK){
00136 getSquaresize(xPencilsPerSlab, sizeX, xsquare);
00137 getSquaresize(yPencilsPerSlab, sizeX, ysquare);
00138 getSquaresize(zPencilsPerSlab, sizeX, zsquare);
00139
00140 CkAssert((xsquare[1]%ysquare[1]==0) || (ysquare[1]%xsquare[1]==0));
00141 CkAssert((zsquare[0]%ysquare[0]==0) || (ysquare[0]%zsquare[0]==0));
00142 }
00143 else {
00144 xsquare[0]=xPencilsPerSlab; xsquare[1]=1;
00145 ysquare[0]=yPencilsPerSlab; ysquare[1]=1;
00146 zsquare[0]=zPencilsPerSlab; zsquare[1]=1;
00147 CkAssert((yPencilsPerSlab%zPencilsPerSlab==0)
00148 || (zPencilsPerSlab%yPencilsPerSlab==0));
00149 }
00150 }
00151 void getSquaresize(int size, int planesize, int *square) {
00152 int squaresize = (int)sqrt((float)size);
00153 if(size==squaresize*squaresize){
00154 square[0]=squaresize;
00155 square[1]=squaresize;
00156 }
00157 else{
00158 while(squaresize>1 && ((size%squaresize!=0) ||
00159 ((size%squaresize==0)&&(planesize%squaresize!=0||planesize%(size/squaresize)!=0)))){
00160 squaresize--;
00161 }
00162 square[1]=squaresize;
00163 square[0]=size/squaresize;
00164 }
00165 CkAssert(size==square[0]*square[1]);
00166 }
00167 };
00168
00169 #include "fftlib.decl.h"
00170
00171 PUPmarshall(NormalFFTinfo)
00172 PUPmarshall(LineFFTinfo)
00173
00174
00175 typedef struct _SlabArrayInfo{
00176 int count;
00177 NormalFFTinfo info;
00178 ComlibInstanceHandle fftcommInstance;
00179
00180
00181 }SlabArrayInfo;
00182
00183
00184
00185
00186
00187 class SlabArray: public CBase_SlabArray {
00188 public:
00189 SlabArray(CkMigrateMessage *m) {}
00190 SlabArray() {}
00191 virtual ~SlabArray() {}
00192
00193
00194 virtual void doneFFT(int id) {
00195 ckout << "NormalSlabArray finished FFT" << endl;
00196 CkExit();
00197 }
00198 virtual void doneIFFT(int id) {
00199 ckout << "NormalSlabArray finished IFFT" << endl;
00200 CkExit();
00201 }
00202
00203
00204 virtual void doFFT(int, int) = 0;
00205 virtual void doIFFT(int, int) = 0;
00206 protected:
00207 CProxy_SlabArray srcProxy, destProxy;
00208 bool fftuseCommlib;
00209
00210 CkVec<SlabArrayInfo*> infoVec;
00211 };
00212
00213
00214
00215
00216
00217 class NormalSlabArray: public CBase_NormalSlabArray {
00218 public:
00219 NormalSlabArray(CkMigrateMessage *m): SlabArray(m) {CkPrintf("migrate constructor called\n");}
00220 NormalSlabArray() {
00221 #if VERBOSE
00222 CkPrintf("Empty constructor called\n");
00223 #endif
00224 fwd2DPlan = bwd2DPlan = (fftwnd_plan) NULL;
00225 fwd1DPlan = bwd1DPlan = (fftw_plan) NULL;
00226 fftuseCommlib = false;
00227
00228 }
00229
00230 NormalSlabArray(NormalFFTinfo &info,
00231 CProxy_NormalSlabArray src, CProxy_NormalSlabArray dest,
00232 bool useCommlib, ComlibInstanceHandle inst);
00233 ~NormalSlabArray();
00234
00235
00236 void acceptDataForFFT(int, complex *, int, int);
00237 void acceptDataForIFFT(int, complex *, int, int);
00238
00239 void doFFT(int src_id = 0, int dst_id = 0);
00240 void doIFFT(int src_id = 0, int dst_id = 0);
00241
00242 void pup(PUP::er &p);
00243
00244 void setup(NormalFFTinfo &info,
00245 CProxy_NormalSlabArray src, CProxy_NormalSlabArray dest,
00246 bool useCommlib=false,
00247 ComlibInstanceHandle inst=ComlibInstanceHandle());
00248 protected:
00249 fftwnd_plan fwd2DPlan, bwd2DPlan;
00250 fftw_plan fwd1DPlan, bwd1DPlan;
00251
00252 void createPlans(NormalFFTinfo &info);
00253 };
00254
00255 class NormalRealSlabArray: public CBase_NormalRealSlabArray {
00256 public:
00257 NormalRealSlabArray(CkMigrateMessage *m): SlabArray(m) {}
00258 NormalRealSlabArray() {
00259 #if VERBOSE
00260 CkPrintf("Empty constructor called\n");
00261 #endif
00262 tempdataPtr = NULL;
00263 rfwd1DXPlan = rbwd1DXPlan = (rfftw_plan) NULL;
00264 fwd1DYPlan = bwd1DYPlan = (fftw_plan) NULL;
00265 fwd1DZPlan = bwd1DZPlan = (fftw_plan) NULL;
00266 rfwd2DXYPlan = rfwd2DXYPlan = (rfftwnd_plan)NULL;
00267 fftuseCommlib = false;
00268
00269 }
00270
00271 NormalRealSlabArray(NormalFFTinfo &info,
00272 CProxy_NormalRealSlabArray, CProxy_NormalRealSlabArray,
00273 bool useCommlib, ComlibInstanceHandle inst);
00274 ~NormalRealSlabArray();
00275
00276
00277 void acceptDataForFFT(int, complex *, int, int);
00278 void acceptDataForIFFT(int, complex *, int, int);
00279
00280 void doFFT(int src_id = 0, int dst_id = 0);
00281 void doIFFT(int src_id = 0, int dst_id = 0);
00282
00283 void pup(PUP::er &p);
00284
00285 void createPlans(NormalFFTinfo &info);
00286
00287 protected:
00288 rfftwnd_plan rfwd2DXYPlan, rbwd2DXYPlan;
00289 rfftw_plan rfwd1DXPlan, rbwd1DXPlan;
00290 fftw_plan fwd1DYPlan, bwd1DYPlan;
00291 fftw_plan fwd1DZPlan, bwd1DZPlan;
00292
00293 NormalFFTinfo *fftinfos[MAX_FFTS];
00294 bool fftuseCommlib;
00295 ComlibInstanceHandle fftcommInstance;
00296
00297 void setup(NormalFFTinfo &info,
00298 CProxy_NormalRealSlabArray, CProxy_NormalRealSlabArray,
00299 bool useCommlib, ComlibInstanceHandle inst);
00300
00301 private:
00302 complex *tempdataPtr;
00303 };
00304
00305 #if 0
00306
00307
00308
00309
00310
00311
00312
00313
00314
00315 class RunDescriptor {
00316 public:
00317 static const short X = 1;
00318 static const short Y = 2;
00319 static const short Z = 2;
00320
00321
00322 RunDescriptor(int x_, int y_, int z_, int l_, int i_ = 0) {x = x_; y = y_; z = z_; length = l_; inc = i_;}
00323 RunDescriptor(const RunDescriptor &rd) {x = rd.x; y = rd.y; z = rd.z; length = rd.length; inc = rd.inc;}
00324
00325
00326 void pup(PUP::er &p) {p|x; p|y; p|z; p|length; p(dir); p|inc;}
00327
00328
00329
00330 int offset(int i) const { return (inc > 0) ? i : -i; }
00331
00332
00333 int x, y, z;
00334 int length;
00335 short dir;
00336 int inc;
00337 };
00338
00339
00340
00341
00342 class SparseSlabArray: public CBase_SparseSlabArray {
00343 public:
00344 SparseSlabArray(CkMigrateMessage *m): SlabArray(m) {}
00345 SparseSlabArray(): SlabArray() {}
00346 SparseSlabArray(FFTinfo &info);
00347 ~SparseSlabArray();
00348
00349 virtual void getRuns(complex **runs, int *numRuns, int *numPoints) const {*runs = NULL; *numRuns = 0; *numPoints = 0;}
00350
00351 private:
00352 int count;
00353 fftwnd_plan fwd2DPlan, bwd2DPlan;
00354 fftw_plan fwd1DPlan, bwd1DPlan;
00355 FFTinfo info;
00356 int *nonZero;
00357 };
00358 #endif
00359
00360 typedef struct _PencilArrayInfo{
00361 int count;
00362 LineFFTinfo info;
00363 ComlibInstanceHandle fftcommInstance;
00364 }PencilArrayInfo;
00365
00366 class SendFFTMsg: public CMessage_SendFFTMsg {
00367 public:
00368 int size;
00369 int id;
00370 int direction;
00371 int ypos;
00372 int zpos;
00373 complex *data;
00374 };
00375
00376
00377 class NormalLineArray : public CBase_NormalLineArray {
00378 public:
00379 NormalLineArray (CkMigrateMessage *m) {}
00380 NormalLineArray () {
00381 fwdplan = bwdplan =NULL;
00382 fftuseCommlib = false;
00383 id = -1;
00384 line = NULL;
00385 }
00386 NormalLineArray (LineFFTinfo &info, CProxy_NormalLineArray _xProxy, CProxy_NormalLineArray _yProxy, CProxy_NormalLineArray _zProxy, bool useCommlib, ComlibInstanceHandle &inst);
00387 ~NormalLineArray () {}
00388 void setup (LineFFTinfo &info, CProxy_NormalLineArray _xProxy, CProxy_NormalLineArray _yProxy, CProxy_NormalLineArray _zProxy, bool useCommlib, ComlibInstanceHandle &inst);
00389 void doFirstFFT(int id, int direction);
00390 void doSecondFFT(int ypos, complex *val, int size, int id, int direction);
00391 void doThirdFFT(int zpos, int ypos, complex *val, int size, int id, int direction);
00392
00393 void doSecondFFT(SendFFTMsg *msg);
00394 void doThirdFFT(SendFFTMsg *msg);
00395
00396 void doFFT(int id, int direction) {doFirstFFT(id, direction);}
00397 virtual void doneFFT(int id, int direction);
00398 void setInstance(int id_) { id = id_;
00399
00400 contribute(sizeof(int), &id_, CkReduction::sum_int);
00401 }
00402 protected:
00403 complex *line;
00404 fftw_plan fwdplan, bwdplan;
00405 int id;
00406
00407 CProxy_NormalLineArray xProxy, yProxy, zProxy;
00408 bool fftuseCommlib;
00409 CkVec<PencilArrayInfo*> infoVec;
00410 };
00411
00412
00413 #define CAREFUL 1
00414
00415 #endif //_fftlib_h_