Logo Coherent WaveBurst  
Reference Guide
Logo
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
testWDM_5.C
Go to the documentation of this file.
1 { // using Time Delay Filters to shift a sine-gaussian in the TF domain
2 
3  int Rate = 1024; // sampling rate (Hz)
4  int Duration = 20; // duration (seconds)
5  int N = Rate*Duration; // number of samples
6  wavearray<double> ts(N); //time series container
7  ts.rate(Rate);
9 
10  // time series: sine gaussian @ 100Hz, Q=8, t = 10s
11  double Q = 8;
12  double dT = 1./Rate;
13  double omega = TMath::TwoPi()*100;
14  double aux = omega/Q;
15  aux *= aux/2;
16  for(int i=0; i<N; i++){
17  double tt = (i*dT-10);
18  ts[i] = sin(omega*tt)*exp(-tt*tt*aux);
19  }
20 
21  pl.plot(ts, const_cast<char*>("APL"), 4, 9.9, 10.1);
22 
23  WDM<double> wdm(32, 64, 4, 8); // define a WDM transform (32 bands)
24  WSeries<double> tfmap; // TF map container
25  tfmap.Forward(ts, wdm);
26 
27 
29  pWDM->setTDFilter(12); // computes TD filters
30 
31  double* tmp = new double[pWDM->nWWS/2]; // will store time delayed coefficients
32  for(int j=9*32; j<11*32; ++j) // time indeces between 9-11s
33  for(int i=0;i<=32; ++i) // frequency bands
34  tmp[j*33+i] = (double)pWDM->getPixelAmplitude(i, j, 13, false); //delay by 13 samples
35 
36  // overwrite original tmap:
37  double* TFMap = pWDM->pWWS; // all TF coefficients are stored in TFMap
38  for(int j=9*32; j<11*32; ++j) // time indeces between 9-11s
39  for(int i=0;i<=32; ++i) // frequency bands
40  TFMap[j*33+i] = tmp[j*33+i];
41 
42  // inverse to time domain and plot shifted signal on top of original (ZOOM around 10s to see)
43  tfmap.Inverse();
44  tfmap.getLayer(ts, 0);
45  pl.plot(ts, const_cast<char*>("SAME"), 2);
46 
47  delete [] tmp;
48 }
double aux
Definition: testWDM_5.C:14
virtual void rate(double r)
Definition: wavearray.hh:123
cout<< "skymap size : "<< L<< endl;for(int l=0;l< L;l++) sm.set(l, l);sm > const_cast< char * >("skymap.dat")
double omega
Definition: testWDM_5.C:13
int Duration
Definition: testWDM_5.C:4
DataType_t * pWWS
Definition: WaveDWT.hh:123
int j
Definition: cwb_net.C:10
i drho i
void setTDFilter(int nCoeffs, int L=1)
Definition: WDM.cc:621
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
double getPixelAmplitude(int, int, int, bool=false)
Definition: WDM.cc:730
watplot pl
Definition: testWDM_5.C:8
WSeries< double > tfmap
Definition: testWDM_5.C:24
wavearray< double > ts(N)
double * tmp
Definition: testWDM_5.C:31
int getLayer(wavearray< DataType_t > &w, double n)
param: n - layer number
Definition: wseries.cc:175
int N
Definition: testWDM_5.C:5
int Rate
unsigned long nWWS
pointer to wavelet work space
Definition: WaveDWT.hh:124
WDM< double > * pWDM
Definition: testWDM_5.C:28
double Q
Definition: testWDM_5.C:11
void Forward(int n=-1)
param: wavelet - n is number of steps (-1 means full decomposition)
Definition: wseries.cc:228
WDM< double > wdm(32, 64, 4, 8)
WaveDWT< DataType_t > * pWavelet
Definition: wseries.hh:438
double dT
Definition: testWDM_5.C:12
void Inverse(int n=-1)
param: n - number of steps (-1 means full reconstruction)
Definition: wseries.cc:273