coherent WaveBurst
Library Reference Guide
home
cWB_docker
git
cWB
library
wat
tutorials
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
22
wavearray<double>
tmp
;
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
}
rnd
TRandom3 rnd(0)
wavearray::rate
virtual void rate(double r)
Definition:
wavearray.hh:141
Duration
int Duration
Definition:
testWDM_2.C:4
wavearray< double >
Rate
int Rate
Definition:
DrawClusterWithSkyplot.C:10
pl
watplot pl
Definition:
DrawSkymapWithSkyplot.C:11
wdm
WDM< double > wdm(nLAYERS, nLAYERS, 6, 10)
i
i drho i
Definition:
cwb_epparameters.C:88
watplot::plot
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:150
ts
wavearray< double > ts(N)
watplot
Definition:
watplot.hh:45
N
int N
Definition:
testWDM_2.C:5
tmp
double * tmp
Definition:
testWDM_5.C:31
WSeries::getLayer
int getLayer(wavearray< DataType_t > &w, double n)
param: n - layer number
Definition:
wseries.cc:193
WSeries< double >
tfmap
WSeries< double > tfmap
Definition:
DrawClusterWithSkyplot.C:26
WDM< double >
WSeries::Forward
void Forward(int n=-1)
param: wavelet - n is number of steps (-1 means full decomposition)
Definition:
wseries.cc:246