Logo Coherent WaveBurst  
Reference Guide
Logo
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
wavearray.hh
Go to the documentation of this file.
1 /**********************************************************
2  * Package: Wavelet Analysis Tool
3  * File name: wavearray.h
4  **********************************************************/
5 
6 #ifndef WAVEARRAY_HH
7 #define WAVEARRAY_HH
8 
9 #include <iostream>
10 #include <string.h>
11 #include <math.h>
12 #include <stdio.h>
13 //#include "slice.h" // DMT CVS
14 //#include <valarray>
15 #ifndef __CINT__
16 #include <valarray>
17 #else
18 #include "wslice.hh" // include definition of std::slice for rootcint
19 //namespace std {
20 // class slice; // DMT
21 //}
22 #endif
23 #include "Wavelet.hh"
24 
25 #ifdef _USE_DMT
26  #include "Time.hh"
27  #include "Interval.hh"
28  #include "TSeries.hh"
29 #endif
30 
31 #include "TFFTComplexReal.h"
32 #include "TFFTRealComplex.h"
33 #include "TNamed.h"
34 
35 template<class DataType_t>
36 class wavearray : public TNamed
37 {
38 
39  public:
40 
41  wavearray(int); // Constructor
42 
43  wavearray(); // Default constructor
44 
45  wavearray(const wavearray<DataType_t>&); // copy Constructor
46 
47 // explicit construction from array
48  template <class T>
49  wavearray(const T *, unsigned int, double=0.);
50 
51  virtual ~wavearray(); // Destructor
52 
53 // operators
54 
56 
57 // operator[](const slice &) sets the Slice object of the wavearray class
58 // the operators above do not use Slice object.
60  virtual DataType_t & operator[](const unsigned int);
61 
62 // check if two wavearrays overlap
63  inline virtual size_t limit() const;
64  inline virtual size_t limit(const std::slice &) const;
65  inline virtual size_t limit(const wavearray<DataType_t> &) const;
66 
67 // the Slice object is used in the operators below and set to
68 // slice(0,N,1) when an operator has been executed for both arrays
69 
73  virtual wavearray<DataType_t>& operator<<(wavearray<DataType_t> &);
74 
75  wavearray<DataType_t>& operator= (const DataType_t);
76  virtual wavearray<DataType_t>& operator+=(const DataType_t);
77  virtual wavearray<DataType_t>& operator-=(const DataType_t);
78  virtual wavearray<DataType_t>& operator*=(const DataType_t);
79 
80  virtual char* operator>>(char*);
81 
82 #ifdef _USE_DMT
83  virtual wavearray<DataType_t>& operator= (const TSeries &);
84 #endif
85 
86 // member functions
87 
88  //: Dump data array to an ASCII file
89  virtual void Dump(const char*, int=0);
90 
91  //: Dump data array to a binary file
92  virtual void DumpBinary(const char*, int = 0);
93 
94  //: Dump data array from a binary file as 16bit words
95  virtual void DumpShort(const char*, int=0);
96 
97  //: Dump object array into file root
98  virtual void DumpObject(const char*);
99 
100  //: Read data array from a binary file
101  virtual void ReadBinary(const char*, int=0);
102 
103  //: Read data array from a binary file
104  virtual void ReadShort(const char*);
105 
106  virtual void FFT(int = 1); // fast Fourier transform
107  virtual void FFTW(int = 1); // fast Fourier transform West
108  virtual void resetFFTW(); // release FFTW memory
109 
110  void Resample(const wavearray<DataType_t> &, double, int=6);
111  void resample(const wavearray<DataType_t> &, double, int=6);
112  virtual void resample(double, int=6);
113  virtual void Resample(double); // resample with FFTW
114 
115  //: data time shift based on FFTW
116  //: T = time shift (sec)
117  virtual void delay(double T);
118 
119  virtual void start(double s) {Start=s; if(Rate>0.) Stop=Start+Size/Rate;};
120  virtual double start() const { return Start; };
121  virtual void stop(double s) {Stop = s; };
122  virtual double stop() const { return Stop; };
123  virtual void rate(double r) {Rate = fabs(r); Stop = Start+Size/Rate; };
124  virtual double rate() const { return Rate; };
125  virtual void edge(double s) {Edge = s; };
126  virtual double edge() const { return Edge; };
127  virtual size_t size() const { return Size; };
128  virtual void setSlice(const std::slice &s) { Slice=s; };
129  virtual std::slice getSlice() const { return Slice; };
130 
131  //: return median
132  //: for data between index1 and index2 (including)
133  virtual double median(size_t=0, size_t=0) const;
134 
135  //: calculate running median of data (x) with window of t seconds (par1)
136  //: put result in input array in if specified
137  //: subtract median from x if fl=true otherwise replace x with median
138  //: move running window by n samples
139  virtual void median(double t, wavearray<DataType_t>* in,
140  bool fl=false, size_t n=1);
141 
142  //: return mean
143  virtual double mean() const;
144 
145  //: return f percentile mean for data
146  // f - >0 positive side pecentile, <0 - double-side percentile
147  virtual double mean(double f);
148 
149  //: return mean for wavearray slice
150  virtual double mean(const std::slice&);
151 
152  //: calculate running mean of data (x) with window of t seconds
153  //: put result in input array in if specified
154  //: subtract mean from x if fl=true otherwise replace x with mean
155  //: move running window by n samples
156  virtual void mean(double t, wavearray<DataType_t>* in,
157  bool fl=false, size_t n=1);
158 
159  //: return sqrt<x^2>
160  virtual double rms();
161 
162  //: return sqrt<x^2> for wavearray slice
163  virtual double rms(const std::slice&);
164 
165  //: calculate running rms of data (x) with window of t seconds
166  //: put result in input array in if specified
167  //: subtract rms from x if fl=true otherwise replace x with rms
168  //: move running window by n samples
169  virtual void rms(double t, wavearray<DataType_t>* in,
170  bool fl=false, size_t n=1);
171 
172  //: return maximum for wavearray
173  virtual DataType_t max() const;
174 
175  //: store max elements max(this->data[i],in.data[i]) in this
176  //: param: input wavearray
177  virtual void max(wavearray<DataType_t> &);
178 
179  //: return minimum for wavearray
180  virtual DataType_t min() const;
181 
182  //: take square root of wavearray elements
183  virtual void SQRT() {
184  for(size_t i=0; i<this->size(); i++)
185  if(this->data[i]>0.) this->data[i]=(DataType_t)sqrt(double(this->data[i]));
186  return;
187  }
188 
189  // multiply wavearray by Hann window
190  inline void hann(void);
191 
192  // sort wavearray using quick sorting algorithm
193  // return DataType_t** pp pointer to array of wavearray pointers
194  // sorted from min to max
195  virtual void waveSort(DataType_t** pp, size_t l=0, size_t r=0) const;
196  // sort wavearray using quick sorting algorithm between left (l) and right (r) index
197  // sorted from min to max: this->data[l]/this->data[r] is min/max
198  virtual void waveSort(size_t l=0, size_t r=0);
199 
200  // split input array of pointers pp[i] between l <= i <= r so that:
201  // *pp[i] < *pp[m] if i < m
202  // *pp[i] > *pp[m] if i > m
203  virtual void waveSplit(DataType_t** pp, size_t l, size_t r, size_t m) const;
204  // split wavearray between l and r & at index m
205  virtual DataType_t waveSplit(size_t l, size_t r, size_t m);
206  // split wavearray at index m so that:
207  // *this->data[i] < *this->datap[m] if i < m
208  // *this->data[i] > *this->datap[m] if i > m
209  virtual void waveSplit(size_t m);
210 
211  //: return rank of sample n ranked against samples between
212  //: left (l) and right (r) boundaries
213  //: rank ranges from 1 to r-l+1
214  virtual int getSampleRank(size_t n, size_t l, size_t r) const;
215 
216  //: rank sample n against absolute value of samples between
217  //: left (l) and right (r) boundaries
218  //: rank ranges from 1 to r-l+1
219  virtual int getSampleRankE(size_t n, size_t l, size_t r) const;
220 
221  // calculate rank
222  // input - fraction of laudest samples
223  // output - min amplitude of laudest samples
224  DataType_t rank(double=0.5) const;
225 
226  // Linear Predictor Filter coefficients
227  // calculate autocorrelation function excluding 3% tails and
228  // solve symmetric Yule-Walker problem
229  // param 1: sifter size
230  // param 2: boundary offset
231  virtual wavearray<double> getLPRFilter(size_t, size_t=0);
232 
233  // generate Symmetric Prediction Error with Split Lattice Algorith
234  // param: filter length in seconds
235  // param: window for filter training
236  // param: boundary offset to account for wavelet artifacts
237  virtual void spesla(double,double,double=0.);
238 
239  // apply lpr filter defined in input wavearray<double>
240  virtual void lprFilter(wavearray<double>&);
241 
242  // calculate and apply lpr filter to this
243  // param: filter length in seconds
244  // param: filter mode: -1/0/1 - backward/symmetric/forward
245  // param: stride for filter training (0 - train on whole TS)
246  // param: boundary offset to account for wavelet artifacts
247  virtual void lprFilter(double,int=0,double=0.,double=0.);
248 
249  // normalization by 31% percentile amplitude
250  // param 1 - time window dT. if = 0 - dT=T, where T is wavearray duration
251  // param 2 - 0 - no whitening, 1 - single whitening, >1 - double whitening
252  // param 3 - boundary offset
253  // param 4 - noise sampling interval (window stride)
254  // the number of measurements is k=int((T-2*offset)/stride)
255  // if stride=0, then stride is set to dT
256  // return: noise array if param2>0, median if param2=0
257  virtual wavearray<double> white(double, int=0,double=0.,double=0.) const;
258 
259  // turn wavearray distribution into exponential using rank statistics
260  // input - time window around a sample to rank the sample
261  virtual void exponential(double);
262 
263  // get sample from wavearray
264  // param - sample index or time
265  inline DataType_t get(size_t i) { return data[i]; }
266  inline DataType_t get(double t);
267 
268  inline long uniform(){ return random(); }
269  inline long rand48(long k=1024){ return long(k*drand48()); }
270 
271  double getStatistics(double &mean, double &rms) const;
272 
273  virtual void resize(unsigned int);
274 
275  //: copy, add and subtruct wavearray array from *this
276  void cpf(const wavearray<DataType_t> &, int=0, int=0, int=0);
277  void add(const wavearray<DataType_t> &, int=0, int=0, int=0);
278  void sub(const wavearray<DataType_t> &, int=0, int=0, int=0);
279 
280  //: append to this *this
281  size_t append(const wavearray<DataType_t> &);
282  //: append to this *this
283  size_t append(DataType_t);
284 
285  // count number of pixels above x
286  size_t wavecount(double x, int n=0){
287  size_t N=0;
288  for(int i=n;i<size()-n;i++) if(data[i]>x) N++;
289  return N;
290  }
291 
292  //: cut data on shorter intervals and average
293  double Stack(const wavearray<DataType_t> &, int);
294  double Stack(const wavearray<DataType_t> &, int, int);
295  double Stack(const wavearray<DataType_t> &, double);
296 
297  // print wavearray parameters
298  void print(); // *MENU*
299  virtual void Browse(TBrowser *b) {print();}
300 
301  DataType_t *data; //! data array
302 
303 // private:
304 
305  size_t Size; // number of elements in the data array
306  double Rate; // data sampling rate
307  double Start; // start time
308  double Stop; // end time
309  double Edge; // buffer length in seconds in the beginning and the end
310  std::slice Slice; // the data slice structure
311 
312  TFFTRealComplex* fftw; //! pointer to direct fftw object
313  TFFTComplexReal* ifftw; //! pointer to inverse fftw object
314 
315  inline static int compare(const void *x, const void *y){
316  DataType_t a = *(*(DataType_t**)x) - *(*(DataType_t**)y);
317  if(a > 0) return 1;
318  if(a < 0) return -1;
319  return 0;
320  }
321 
322  ClassDef(wavearray,2)
323 };
324 
325 template<class DataType_t>
327 { return Slice.start() + (Slice.size()-1)*Slice.stride() + 1; }
328 
329 template<class DataType_t>
331 { return s.start() + (s.size()-1)*s.stride() + 1; }
332 
333 template<class DataType_t>
335 {
336  size_t N = a.Slice.size();
337  if(N>Slice.size()) N=Slice.size();
338  return Slice.start() + (N-1)*Slice.stride() + 1;
339 }
340 
341 template<class Tout, class Tin>
342 inline void waveAssign(wavearray<Tout> &aout, wavearray<Tin> &ain)
343 {
344  size_t N = ain.size();
345  aout.rate(ain.rate());
346  aout.start(ain.start());
347  aout.stop(ain.stop());
348  aout.Slice = std::slice(0,N,1);
349  if (N != aout.size()) aout.resize(N);
350  for (unsigned int i=0; i < N; i++) aout.data[i] = (Tout)ain.data[i];
351 }
352 
353 template<class DataType_t>
355 {
356  double phi = 2*3.141592653589793/size();
357  double www = sqrt(2./3.);
358  int nn = size();
359  for(int i=0; i<nn; i++) data[i] *= DataType_t(www*(1.-cos(i*phi)));
360 }
361 
362 template<class DataType_t>
363 inline DataType_t wavearray<DataType_t>::get(double t)
364 {
365  t -= this->Start;
366  if(t<0. || t>Size/Rate) {
367  printf("wavearray<DataType_t>::get(double t): time out of range\n");
368  return 0;
369  }
370  return data[size_t(t*Rate)];
371 }
372 
373 #endif // WAVEARRAY_HH
374 
375 
376 
377 
378 
379 
380 
381 
382 
383 
384 
385 
wavearray< double > t(hp.size())
size_t append(const wavearray< DataType_t > &)
Definition: wavearray.cc:775
virtual size_t size() const
Definition: wavearray.hh:127
tuple f
Definition: cwb_online.py:91
printf("total live time: non-zero lags = %10.1f \n", liveTot)
virtual void rate(double r)
Definition: wavearray.hh:123
long uniform()
Definition: wavearray.hh:268
size_t Size
data array
Definition: wavearray.hh:305
virtual void ReadShort(const char *)
Definition: wavearray.cc:422
wavearray< double > a(hp.size())
virtual double start() const
Definition: wavearray.hh:120
virtual void edge(double s)
Definition: wavearray.hh:125
int n
Definition: cwb_net.C:10
void print()
Definition: wavearray.cc:2203
virtual double edge() const
Definition: wavearray.hh:126
double getStatistics(double &mean, double &rms) const
Definition: wavearray.cc:2105
cout<< endl;cout<< "ts size = "<< ts.size()<< " ts rate = "<< ts.rate()<< endl;tf.Forward(ts, wdm);int levels=tf.getLevel();cout<< "tf size = "<< tf.size()<< endl;double dF=tf.resolution();double dT=1./(2 *dF);cout<< "rate(hz) : "<< RATE<< "\t layers : "<< nLAYERS<< "\t dF(hz) : "<< dF<< "\t dT(ms) : "<< dT *1000.<< endl;int itime=TIME_PIXEL_INDEX;int ifreq=FREQ_PIXEL_INDEX;int index=(levels+1)*itime+ifreq;double time=itime *dT;double freq=(ifreq >0)?ifreq *dF:dF/4;cout<< endl;cout<< "PIXEL TIME = "<< time<< " sec "<< endl;cout<< "PIXEL FREQ = "<< freq<< " Hz "<< endl;cout<< endl;wavearray< double > x
wavearray< DataType_t > & operator=(const wavearray< DataType_t > &)
Definition: wavearray.cc:92
virtual int getSampleRankE(size_t n, size_t l, size_t r) const
Definition: wavearray.cc:1381
double Rate
Definition: wavearray.hh:306
void add(const wavearray< DataType_t > &, int=0, int=0, int=0)
Definition: wavearray.cc:728
virtual void ReadBinary(const char *, int=0)
Definition: wavearray.cc:392
virtual void setSlice(const std::slice &s)
Definition: wavearray.hh:128
double Stop
Definition: wavearray.hh:308
TFFTRealComplex * fftw
Definition: wavearray.hh:312
virtual std::slice getSlice() const
Definition: wavearray.hh:129
virtual wavearray< DataType_t > & operator*=(wavearray< DataType_t > &)
Definition: wavearray.cc:258
virtual ~wavearray()
Definition: wavearray.cc:82
virtual size_t limit() const
Definition: wavearray.hh:326
Long_t size
virtual double median(size_t=0, size_t=0) const
Definition: wavearray.cc:1558
int m
Definition: cwb_net.C:10
virtual void start(double s)
Definition: wavearray.hh:119
virtual void Browse(TBrowser *b)
Definition: wavearray.hh:299
i drho i
double Edge
Definition: wavearray.hh:309
virtual DataType_t min() const
Definition: wavearray.cc:1350
void waveAssign(wavearray< Tout > &aout, wavearray< Tin > &ain)
Definition: wavearray.hh:342
#define N
virtual double mean() const
Definition: wavearray.cc:1053
DataType_t rank(double=0.5) const
Definition: wavearray.cc:1719
virtual double rms()
Definition: wavearray.cc:1188
virtual DataType_t max() const
Definition: wavearray.cc:1316
virtual wavearray< double > getLPRFilter(size_t, size_t=0)
Definition: wavearray.cc:1893
virtual void DumpObject(const char *)
Definition: wavearray.cc:382
float phi
void sub(const wavearray< DataType_t > &, int=0, int=0, int=0)
Definition: wavearray.cc:754
virtual char * operator>>(char *)
Definition: wavearray.cc:2177
size_t wavecount(double x, int n=0)
Definition: wavearray.hh:286
static int compare(const void *x, const void *y)
pointer to inverse fftw object
Definition: wavearray.hh:315
size_t size() const
Definition: wslice.hh:71
TFFTComplexReal * ifftw
pointer to direct fftw object
Definition: wavearray.hh:313
void Resample(const wavearray< DataType_t > &, double, int=6)
Definition: wavearray.cc:622
void hann(void)
Definition: wavearray.hh:354
std::slice Slice
Definition: wavearray.hh:310
virtual double stop() const
Definition: wavearray.hh:122
int k
virtual void resetFFTW()
Definition: wavearray.cc:959
virtual wavearray< DataType_t > & operator+=(wavearray< DataType_t > &)
Definition: wavearray.cc:173
virtual int getSampleRank(size_t n, size_t l, size_t r) const
Definition: wavearray.cc:1361
int Rate
long rand48(long k=1024)
Definition: wavearray.hh:269
virtual void DumpBinary(const char *, int=0)
Definition: wavearray.cc:335
virtual void FFTW(int=1)
Definition: wavearray.cc:878
regression r
Definition: Regression_H1.C:44
s s
Definition: cwb_net.C:137
virtual void SQRT()
Definition: wavearray.hh:183
double T
Definition: testWDM_4.C:11
virtual void lprFilter(wavearray< double > &)
Definition: wavearray.cc:1808
virtual void waveSort(DataType_t **pp, size_t l=0, size_t r=0) const
Definition: wavearray.cc:1403
ifstream in
virtual void delay(double T)
Definition: wavearray.cc:578
virtual void DumpShort(const char *, int=0)
Definition: wavearray.cc:356
virtual wavearray< double > white(double, int=0, double=0., double=0.) const
Definition: wavearray.cc:2002
virtual void stop(double s)
Definition: wavearray.hh:121
virtual wavearray< DataType_t > & operator-=(wavearray< DataType_t > &)
Definition: wavearray.cc:208
double fabs(const Complex &x)
Definition: numpy.cc:37
virtual void spesla(double, double, double=0.)
Definition: wavearray.cc:1752
void resample(const wavearray< DataType_t > &, double, int=6)
Definition: wavearray.cc:485
virtual wavearray< DataType_t > & operator[](const std::slice &)
Definition: wavearray.cc:277
virtual void Dump(const char *, int=0)
Definition: wavearray.cc:301
int l
Definition: cbc_plots.C:434
DataType_t get(size_t i)
Definition: wavearray.hh:265
virtual void exponential(double)
Definition: wavearray.cc:1662
virtual double rate() const
Definition: wavearray.hh:124
DataType_t * data
Definition: wavearray.hh:301
double Stack(const wavearray< DataType_t > &, int)
Definition: wavearray.cc:973
double Start
Definition: wavearray.hh:307
size_t stride() const
Definition: wslice.hh:75
virtual void resize(unsigned int)
Definition: wavearray.cc:445
virtual void waveSplit(DataType_t **pp, size_t l, size_t r, size_t m) const
Definition: wavearray.cc:1481
wavearray< double > y
Definition: Test10.C:31
virtual void FFT(int=1)
Definition: wavearray.cc:814
size_t start() const
Definition: wslice.hh:67
void cpf(const wavearray< DataType_t > &, int=0, int=0, int=0)
Definition: wavearray.cc:699