Logo Coherent WaveBurst  
Reference Guide
Logo
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
testWDM_2.C
Go to the documentation of this file.
1 { // access TF data
2 
3  int Rate = 1024; // sampling rate (Hz)
4  int Duration = 100; // duration (seconds)
5  int N = Rate*Duration; // number of samples
6  wavearray<double> ts(N); //time series container
7  ts.rate(Rate);
8 
9  // time series is filled with white noise data:
10  TRandom3 rnd(0);
11  for(int i=0; i<N; i++) ts[i] = rnd.Gaus();
12 
13  // produce the TF map:
14  WSeries<double> tfmap; // TF map container
15  WDM<double> wdm(32, 64, 4, 8); // define a WDM transform (32 bands)
16  tfmap.Forward(ts, wdm); // apply the WDM to the time series
17  watplot pl; // create canvas
18  pl.plot(ts,const_cast<char*>("APL"),1); // plot time series
19 
20  // extract layer 2
21 
23  tfmap.getLayer(tmp, 2); // extract TF amplitudes of the 3rd band
24  pl.plot(tmp,const_cast<char*>("SAME"),2); // superimpose layer amplitudes
25 
26  // accessing the quadrature bands is done using negative numbers:
27  tfmap.getLayer(tmp, -2); // 3rd band of the quadrature
28  tfmap.getLayer(tmp, -0.5); // 1st band of the quadrature
29 
30 }
TRandom3 rnd(0)
virtual void rate(double r)
Definition: wavearray.hh:123
int Duration
Definition: testWDM_2.C:4
cout<< "skymap size : "<< L<< endl;for(int l=0;l< L;l++) sm.set(l, l);sm > const_cast< char * >("skymap.dat")
watplot pl
WDM< double > wdm(nLAYERS, nLAYERS, 6, 10)
i drho i
void plot(wavearray< double > &, char *=NULL, int=1, double=0., double=0., bool=false, float=0., float=0., bool=false, float=0., bool=false)
Definition: watplot.cc:132
wavearray< double > ts(N)
int N
Definition: testWDM_2.C:5
double * tmp
Definition: testWDM_5.C:31
int getLayer(wavearray< DataType_t > &w, double n)
param: n - layer number
Definition: wseries.cc:175
WSeries< double > tfmap
int Rate
void Forward(int n=-1)
param: wavelet - n is number of steps (-1 means full decomposition)
Definition: wseries.cc:228