Logo coherent WaveBurst  
Library Reference Guide
Logo
testWDM_7.C
Go to the documentation of this file.
1 { // residuals with increased precision (4th parameter in the WDM constructor)
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, 6, 12); // define a WDM transform (32 bands, increased precision (12))
16  tfmap.Forward(ts, wdm); // apply the WDM to the time series
17 
18  wavearray<double> ts2(ts);
19 
20  tfmap.Inverse(); // inverse transform
21  tfmap.getLayer(ts2, 0); // move the time domain data in ts2
22  ts -= ts2; // substract from orginal data and plot:
23 
24  TH1F* h = new TH1F("h", "Residuals, P=12", 100, -3e-6, 3e-6);
25  for(int i=0; i<N; ++i)h->Fill(ts[i]);
26  h->Draw();
27 
28  // NOTE: if we call tfmap.Inverse(-2) an inverse transform
29  // that uses the quadrature coefficients is performed
30 }
TRandom3 rnd(0)
virtual void rate(double r)
Definition: wavearray.hh:141
int Duration
Definition: testWDM_7.C:4
WDM< double > wdm(nLAYERS, nLAYERS, 6, 10)
i drho i
wavearray< double > ts(N)
wavearray< double > h
Definition: Regression_H1.C:25
int N
Definition: testWDM_7.C:5
int getLayer(wavearray< DataType_t > &w, double n)
param: n - layer number
Definition: wseries.cc:193
WSeries< double > tfmap
double e
void Forward(int n=-1)
param: wavelet - n is number of steps (-1 means full decomposition)
Definition: wseries.cc:246
void Inverse(int n=-1)
param: n - number of steps (-1 means full reconstruction)
Definition: wseries.cc:291