Logo Coherent WaveBurst  
Reference Guide
Logo
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
WDM.hh
Go to the documentation of this file.
1 // Wavelet Analysis Tool:
2 // Sergey Klimenko, University of Florida
3 //--------------------------------------------------------------------
4 // Implementation of
5 // Wilson-Daubechies transform
6 // References:
7 //--------------------------------------------------------------------
8 
9 
10 #ifndef WDM_HH
11 #define WDM_HH
12 
13 #include "SymmArray.hh"
14 #include "SymmArraySSE.hh"
15 #include "SymmObjArray.hh"
16 #include "WaveDWT.hh"
17 #include "wavearray.hh"
18 
19 #define MAXBETA 8
20 
21 template<class DataType_t> class SSeries;
22 
23 template<class DataType_t>
24 class WDM : public WaveDWT<DataType_t>
25 {
26 public:
27  // dummy constructor
28  WDM();
29 
30  // constructor
31  WDM(int, int, int, int);
32 
33  // WDMK constructor
34  WDM(int);
35 
36  //: copy constructors
37  WDM(const WDM &);
38 
39  //: destructor
40  virtual ~WDM();
41 
42  //: Duplicate on heap
43  virtual WDM* Clone() const;
44 
45  //: light-weight duplicate
46  virtual WDM* Init() const;
47 
48  // get maximum possible level of decomposition
49  int getMaxLevel();
50 
51  void forward(int, int);
52  void inverse(int, int);
53 
54  // get the time domain representation of the basis function corresponding to pixel (m,n)
55  // param1: frequency index
56  // param2: time index
57  // param3: where to store it
58  // returns translation constant (time shift not included in w)
59  int getBaseWave(int m, int n, SymmArray<double>& w);
60 
61  // same as above but for the quadrature
62  int getBaseWaveQ(int m, int n, SymmArray<double>& w);
63 
64  // get the time domain representation of the basis function corresponding to pixel j
65  // param1: pixel index
66  // param2: where to store it
67  // param3: Quadrature flag (true for quadrature)
68  // returns translation constant (sampling steps; time shift not included in w)
69  int getBaseWave(int j, wavearray<double>& w, bool Quad=false);
70 
71 
72  //: get array index of the first sample for (level,layer)
73  int getOffset(int, int);
74 
75  //: forward transfom
76  //: param: -1 - orthonormal, 0 - power map, >0 - upsampled map
77  void t2w(int);
78 
79  //: inverse transform (flag == -2 means do inverse of Quadrature)
80  void w2t(int flag);
81  void w2tQ(int);
82  // returns pixel amplitude
83  // param1: frequency index
84  // param2: time index
85  // param3: delay index
86  // param4: 00 phase - true, 90 phase - false
87  double getPixelAmplitude(int, int, int, bool=false);
88  double getPixelAmplitudeSSEOld(int, int, int, bool=false);
89  float getPixelAmplitudeSSE(int m, int n, int dT, bool Quad);
90  void getPixelAmplitudeSSE(int m, int n, int t1, int t2, float* r, bool Quad);
91 
92 
93  double TimeShiftTest(int dt);
94  double TimeShiftTestSSE(int dt);
95 
96  // override getTDamp from WaveDWT
97  // return time-delayed amplitudes for sample n and TD index m:
98  float getTDamp(int n, int m, char c='p');
99 
100 
101  // override getTDvec from WaveDWT
102  // return array of time-delayed amplitudes in the format:
103  // [-n*dt , 0, n*dt, -n*dt, 0, n*dt] - total 2(2n+1) amplitudes
104  // where n = k*LWDM (k - second parameter)
105  // param - index in TF map
106  // param - range of h(t) sample delays k
107  // param - 'a','A' - delayed amplitudes, 'p','P' - delayed power.
108  wavearray<float> getTDvec(int j, int K, char c='p');
109 
110  wavearray<float> getTDvecSSE(int j, int K, char c, SSeries<double>* pss);
111  //wavearray<float> getTDvecSSE(int j, int K, char c);
112 
113  // similar to getTDvec
114  void getTFvec(int j, wavearray<float>& r);
115 
116  // void transform2();
117  // double pixel(int f0, int t0);
118 
119  static double *Cos[MAXBETA], *Cos2[MAXBETA], *SinCos[MAXBETA];
121  static int objCounter;
122 
123  void initFourier();
124 
125  // getFilter(): returns 3 different WDM filters
126  // param1: filter length
127 
129  // set symmetric time-delay filter arrays used to calculate TD filter
130  void setTDFilter(int nCoeffs, int L=1); //L determines fractional increment tau/L
131  wavearray<double> getTDFilter2(int n, int L);
132  wavearray<double> getTDFilter1(int n, int L);
133 
134 
135  // override getSlice() from WaveDWT
136  std::slice getSlice(double n);
137 
138  void SetTFMap();
139  //double PrintTDFilters(int m, int dt, int nC);
140  // return size of global TD array
141  inline size_t Last(int n=0) { return T0.Last(); }
142  // return half size of TD filter
143  virtual size_t getTDFsize() { return T0.Last() ? T0[0].Last() : 0; }
144 
145 
146  int BetaOrder; // beta function order for Meyer
147  int precision; // wavelet precision
148  int KWDM; // K - parameter
149  int LWDM; // unit time delay is tau/LWDM where tau is 1/hot_rate
151 
152  SymmObjArray<SymmArraySSE<float> > T0; // time-delay filters
153  SymmObjArray<SymmArraySSE<float> > Tx; // time-delay filters
155 
156  DataType_t** TFMap00; //! pointer to 0-phase data, by default not initialized
157  DataType_t** TFMap90; //! pointer to 90-phase data, by default not initialized
158  void (*SSE_TDF)(); //!
159  float* td_buffer; //!
160  float* td_data; //!
162 
163 protected:
164  void initSSEPointers();
165 
166  ClassDef(WDM,2)
167 
168 }; // class WDM
169 
170 #endif // WDM_HH
171 
DataType_t ** TFMap00
Definition: WDM.hh:156
int getBaseWave(int m, int n, SymmArray< double > &w)
Definition: WDM.cc:305
void w2t(int flag)
Definition: WDM.cc:1521
wavearray< float > cosTD
Definition: WDM.hh:154
static double SinCosSize[MAXBETA]
Definition: WDM.hh:120
#define MAXBETA
Definition: WDM.hh:19
virtual size_t getTDFsize()
Definition: WDM.hh:143
static double CosSize[MAXBETA]
Definition: WDM.hh:120
wavearray< float > getTDvecSSE(int j, int K, char c, SSeries< double > *pss)
Definition: WDM.cc:1183
unsigned int Last()
Definition: SymmObjArray.hh:22
int n
Definition: cwb_net.C:10
void initFourier()
Definition: WDM.cc:62
static double * SinCos[MAXBETA]
Definition: WDM.hh:119
static double * Cos2[MAXBETA]
Definition: WDM.hh:119
int m
Definition: cwb_net.C:10
virtual WDM * Clone() const
return: Wavelet* - duplicate of *this, allocated on heap
Definition: WDM.cc:274
int j
Definition: cwb_net.C:10
int getMaxLevel()
Definition: WDM.cc:1334
int BetaOrder
Definition: WDM.hh:146
SymmObjArray< SymmArraySSE< float > > Tx
Definition: WDM.hh:153
void setTDFilter(int nCoeffs, int L=1)
Definition: WDM.cc:621
void w2tQ(int)
Definition: WDM.cc:1613
double getPixelAmplitude(int, int, int, bool=false)
Definition: WDM.cc:730
void inverse(int, int)
Definition: WDM.cc:1358
int getOffset(int, int)
Definition: WDM.cc:1342
static int objCounter
Definition: WDM.hh:121
WDM()
Definition: WDM.cc:68
void getTFvec(int j, wavearray< float > &r)
Definition: WDM.cc:1259
wavearray< float > sinTDx
Definition: WDM.hh:154
wavearray< double > w
Definition: Test1.C:27
wavearray< float > getTDvec(int j, int K, char c='p')
Definition: WDM.cc:1141
double getPixelAmplitudeSSEOld(int, int, int, bool=false)
Definition: WDM.cc:815
int getBaseWaveQ(int m, int n, SymmArray< double > &w)
Definition: WDM.cc:356
static double Cos2Size[MAXBETA]
Definition: WDM.hh:120
wavearray< float > sinTD
Definition: WDM.hh:154
static double * Cos[MAXBETA]
Definition: WDM.hh:119
void SetTFMap()
Definition: WDM.cc:429
wavearray< double > getFilter(int n)
Definition: WDM.cc:480
wavearray< double > getTDFilter2(int n, int L)
Definition: WDM.cc:579
int LWDM
Definition: WDM.hh:149
float * td_data
Definition: WDM.hh:160
int KWDM
Definition: WDM.hh:148
void forward(int, int)
Definition: WDM.cc:1350
DataType_t ** TFMap90
pointer to 0-phase data, by default not initialized
Definition: WDM.hh:157
virtual ~WDM()
Definition: WDM.cc:255
void t2w(int)
Definition: WDM.cc:1367
float getPixelAmplitudeSSE(int m, int n, int dT, bool Quad)
Definition: WDM.cc:911
double dt
regression r
Definition: Regression_H1.C:44
Definition: WDM.hh:24
size_t Last(int n=0)
Definition: WDM.hh:141
double TimeShiftTest(int dt)
Definition: WDM.cc:1303
float getTDamp(int n, int m, char c='p')
Definition: WDM.cc:1085
wavearray< double > wdmFilter
Definition: WDM.hh:150
SymmArraySSE< float > td_halo[6]
Definition: WDM.hh:161
virtual WDM * Init() const
Definition: WDM.cc:283
float * td_buffer
Definition: WDM.hh:159
std::slice getSlice(double n)
Definition: WDM.cc:445
void(* SSE_TDF)()
pointer to 90-phase data, by default not initialized
Definition: WDM.hh:158
wavearray< double > getTDFilter1(int n, int L)
Definition: WDM.cc:533
double TimeShiftTestSSE(int dt)
Definition: WDM.cc:1318
TH1 * t1
double dT
Definition: testWDM_5.C:12
int precision
Definition: WDM.hh:147
void initSSEPointers()
Definition: WDM.cc:683
SymmObjArray< SymmArraySSE< float > > T0
Definition: WDM.hh:152