Logo Coherent WaveBurst  
Reference Guide
Logo
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
wseries.hh
Go to the documentation of this file.
1 // Wavelet Analysis Tool
2 // universal data container for wavelet transforms
3 // used with DMT and ROOT
4 //
5 //$Id: wseries.hh,v 1.1 2005/05/16 06:05:10 igor Exp $
6 
7 #ifndef WSERIES_HH
8 #define WSERIES_HH
9 
10 #ifndef WAVEARRAY_HH
11 #include "wavearray.hh"
12 #endif
13 
14 #include <vector>
15 #include <list>
16 #include "WaveDWT.hh"
17 #include <complex>
18 #include "wavecomplex.hh"
19 #include "TNamed.h"
20 #include "TH1F.h"
21 
23 typedef std::vector<int> vector_int;
24 
25 
26 template<class DataType_t>
27 class WSeries : public wavearray<DataType_t>
28 {
29 
30  public:
31 
32  // constructors
33 
34  //: Default constructor
35  WSeries();
36 
37  //: Construct WSeries for specific wavelet type
38  //+ default constructor
39  explicit WSeries(const Wavelet &w);
40 
41  //: Construct from wavearray
42  //!param: value - data to initialize the WSeries object
43  explicit WSeries(const wavearray<DataType_t>& value,
44  const Wavelet &w);
45 
46  //: Copy constructor
47  //!param: value - object to copy from
48  WSeries(const WSeries<DataType_t>& value);
49 
50  //: destructor
51  virtual ~WSeries();
52 
53  // operators
54 
57  WSeries<DataType_t>& operator= (const DataType_t);
58 
59  // operator[](const slice &) sets the Slice object of the wavearray class
60  virtual WSeries<DataType_t>& operator[](const std::slice &);
61 
62  //: operators for WSeries objects, which can have different length
63  //: it is required they have the same type of transform (standard or binary)
64  //: and the same size of approximation levels.
65  //: warning: there is no check that approximation levels have the same
66  //: sampling rate.
70 
71  // just to trick ANSI standard
75  virtual WSeries<DataType_t>& operator+=(const DataType_t);
76  virtual WSeries<DataType_t>& operator-=(const DataType_t);
77  virtual WSeries<DataType_t>& operator*=(const DataType_t);
78 
79  // multiply layer by layer this and input wseries
80  void mul(WSeries<DataType_t> &);
81 
82  //: Dump data array to an ASCII file
83  virtual void Dump(const char*, int=0);
84 
85  // accessors
86 
87  //: Get maximum possible level of decompostion
88  int getMaxLevel();
89 
90  //: Get level of decompostion
91  inline int getLevel(){ return pWavelet->m_Level; }
92 
93  //: set level of decompostion
94  inline void setLevel(size_t n){ pWavelet->m_Level = n; }
95 
96  //: Set black pixel probability
97  inline void setbpp(double f) { bpp=f; return; }
98  //: Get black pixel probability
99  inline double getbpp() const { return bpp; }
100 
101  //: Set wavelet rate (0-level)
102  inline void wrate(double r) {wRate=r; return;}
103  //: Get wavelet rate
104  inline double wrate() const { return wRate; }
105 
106  //: Set low frequency boundary
107  inline void setlow(double f) {
108  f_low=f>0. ? f : 0.; return;
109  }
110  //: Get low frequency boundary
111  inline double getlow() const { return f_low; }
112 
113  //: Set high frequency boundary
114  inline void sethigh(double f) {
115  f_high = f; return;
116  }
117  //: get high frequency boundary
118  inline double gethigh() const { return f_high; }
119 
120  //: Get max layer of decompostion
121  inline int maxLayer(){
122  return pWavelet->BinaryTree() ? (1<<getLevel())-1 : getLevel();
123  }
124 
125  // number of samples at zero level
126  inline size_t sizeZero(){return pWavelet->getSlice(0).size();}
127  // number of samples in the original time series
128  inline size_t xsize(){return pWavelet->nSTS;}
129 
130  // maximum sample index in 00 phase
131  inline size_t maxIndex(){return sizeZero()*(maxLayer()+1)-1;}
132 
133  //: Get slice structure for specified layer
134  inline std::slice getSlice(double n){ return pWavelet->getSlice(n); }
135 
136  //: Get wavelet layer frequency resolution
137  inline double resolution(int=0){
138  return frequency(1)-frequency(0);
139  }
140 
141  //: Get central frequency of a layer
142  // l - layer number.
143  // TF map binary dyadic WDM
144  // zero layer Fc dF/2 Fo 0
145  // non-zero Fc +n*dF Fn +n*dF
146  double frequency(int l);
147 
148  //: Get layer index for input frequency
149  // f - frequency Hz
150  // TF map binary dyadic WDM
151  // zero layer Fc dF/2 Fo 0
152  // non-zero Fc +n*dF Fn +n*dF
153  int layer(double f);
154 
155  //: Extract wavelet coefficients from specified layer
156  //!param: n - layer number
157  int getLayer(wavearray<DataType_t> &w, double n);
158 
159  //: replace wavelet data for specified layer with data from Sequence
160  //!param: n - layer number
161  void putLayer(wavearray<DataType_t> &, double n);
162 
163  // extract wavelet amplitude for layer m and wavelet time index n
164  // or sample(n,m)
165  // !param: n - wavelet time index
166  // !param: m - layer number
167  inline DataType_t getSample(int n, double m) {
168  std::slice S = this->getSlice(m);
169  return (n<S.size()) ? this->data[n*S.stride()+S.start()] : 0;
170  }
171 
172  // extract wavelet amplitude for layer m and wavelet time index n
173  // or sample(n,m)
174  // !param: a - new amplitude
175  // !param: n - wavelet time index
176  // !param: m - layer number
177  inline void putSample(DataType_t a, int n, double m) {
178  std::slice S = this->getSlice(m);
179  if(n<S.size()) this->data[n*S.stride()+S.start()]=a;
180  }
181 
182  // mutators
183 
184  virtual void resize(unsigned int);
185  virtual void resample(double, int=6);
186 
187  //: initialize wavelet parameters from Wavelet object
188  void setWavelet(const Wavelet &w);
189  //: return true if WDM transform
190  bool isWDM() {return pWavelet->m_WaveType==WDMT ? true : false;}
191 
192  //: Perform n steps of forward wavelet transform
193  //!param: wavelet - n is number of steps (-1 means full decomposition)
194  // WDM - n=-1 - otrhonormal, n=0 - power, n>0 - upsampled power
195  void Forward(int n = -1);
196  void Forward(wavearray<DataType_t> &, int n = -1);
197  void Forward(wavearray<DataType_t> &, Wavelet &, int n = -1);
198 
199  //: Perform n steps of inverse wavelet transform
200  //!param: n - number of steps (-1 means full reconstruction)
201  void Inverse(int n = -1);
202 
203  // bandpass data and store in TF domain
204  // ts - input time series
205  // flow - low frequence boundary
206  // fhigh - high frequency boundary
207  // n - decomposition parameter
208  void bandpass(wavearray<DataType_t> &ts, double flow, double fhigh, int n=-1);
209 
210  // set wseries coefficients to a in layers between
211  // flow - low frequence boundary
212  // fhigh - high frequency boundary
213  // a - value
214  void bandpass(double flow, double fhigh, double a=0.);
215 
216  // maxEnergy: put maximum energy of delayed samples in this
217  // param: wavearray - input time series
218  // param: wavelet - wavelet used for the transformation
219  // param: double - range of time delays
220  // param: int - downsample factor to obtain coarse TD steps
221  // param: int - clustering mode
222  // returns median energy
223  double maxEnergy(wavearray<DataType_t> &ts, Wavelet &w, double=0, int=1, int=0, TH1F* = NULL);
224 
225  // wdmPacket: converts this to WDM packet series described by pattern
226  // patterns: "/" - chirp, "\" - ringdown, "|" - delta, "*" - pixel
227  // opt = 'e' / 'E' - returns pattern / packet energy
228  // opt = 'l' / 'L' - returns pattern / packet likelihood
229  // opt = 'a' / 'A' - returns packet amplitudes
230  // patterns: "/" - chirp, "\" - ringdown, "|" - delta, "*" - pixel
231  // param: pattern = 0 - "*" single pixel standard search
232  // param: pattern = 1 - "3|" packet
233  // param: pattern = 2 - "3-" packet
234  // param: pattern = 3 - "3/" packet - chirp
235  // param: pattern = 4 - "3\" packet - ringdown
236  // param: pattern = 5 - "5/" packet - chirp
237  // param: pattern = 6 - "5\" packet - ringdown
238  // param: pattern = 7 - "3+" packet
239  // param: pattern = 8 - "3x" cross packet
240  // param: pattern = 9 - "9p" 9-pixel square packet
241  // pattern = else - "*" single pixel standard search
242  // mean of the packet noise distribution is mu=2*K+1, where K is
243  // the effective number of pixels in the packet (K may not be integer)
244  double wdmPacket(int pattern, char opt='L', TH1F* = NULL);
245  double Gamma2Gauss(TH1F* = NULL);
246 
247  // create a wavescan object
248  // produce multi-resolution TF series of input time series x
249  // pws - array of pointers to input time-frequency series
250  // N - number of resolutions
251  // hist - diagnostic histogram
252  void wavescan(WSeries<DataType_t>**, int, TH1F* = NULL);
253 
254 //++++++++++++++ wavelet data conditioning +++++++++++++++++++++++++++
255 
256  //: calculate running medians with window of t seconds
257  //: and subtract the median from this (false key) or
258  //: calculate median for abs(this) and normalize this (true)
259  virtual void median(double t, bool norm=false);
260 
261 
262  //: apply linear predictor to each layer.
263  // param: filter length in seconds
264  // param: filter mode: -1/0/1 - backward/symmetric/forward
265  // param: stride for filter training (0 - train on whole TS)
266  // param: boundary offset to account for wavelet artifacts
267  virtual void lprFilter(double,int=0,double=0.,double=0.);
268 
269  //: tracking of noise non-stationarity and whitening.
270  // param 1 - time window dT. if = 0 - dT=T, where T is wavearray duration
271  // param 2 - mode: 0 - no whitening, 1 - single whitening, >1 - double whitening
272  // mode <0 - whitening using guadrature (WDM wavelet only)
273  // param 3 - boundary offset
274  // param 4 - noise sampling interval (window stride)
275  // the number of measurements is k=int((T-2*offset)/stride)
276  // if stride=0, then stride is set to dT
277  // return: noise array if param2>0, median if param2=0
278  //!what it does: each wavelet layer is devided into k intervals.
279  //!The data for each interval is sorted and the following parameters
280  //!are calculated: median and the amplitude
281  //!corresponding to 31% percentile (wp). Wavelet amplitudes (w) are
282  //!normalized as w' = (w-median(t))/wp(t), where median(t) and wp(t)
283  //!is a linear interpolation between (median,wp) measurements for
284  //!each interval.
285  virtual WSeries<double> white(double,int,double=0.,double=0.);
286 
287  //: whiten TF map by using input noise array
288  virtual bool white(WSeries<double> ws, int mode=0);
289 
290  //: local whitening, works only for binary wavelets.
291  //: returns array of noise rms for wavelet layers
292  //!param: n - number of decomposition steps
293  //!algorithm:
294  //! 1) do forward wavelet transform with n decomposition steps
295  //! 2) whiten wavelet layers and calculate noise rms as
296  //! 1/Sum(1/var)
297  //! 3) do inverse wavelet transform with n reconstruction steps
298  virtual wavearray<double> filter(size_t);
299 
300  //: works only for binary wavelets.
301  //: calculates, corrects and returns noise variability
302  //!param: first - time window to calculate normalization constants
303  //! second - low frequency boundary for correction
304  //! third - high frequency boundary for correction
305  //!algorithm:
306  //! 1) sort wavelet amplitudes with the same time stamp
307  //! 2) calculate left(p) and right(p) amplitudes
308  //! put (right(p)-left(p))/2 into output array
309  //! 3) if first parameter >0 - devide WSeries by average
310  //! variability
311  virtual WSeries<float> variability(double=0., double=-1., double=-1.);
312 
313  //: Selection of a fixed fraction of pixels
314  //: reduced wavelet amplitudes are stored in this
315  //: Returns fraction of non-zero coefficients.
316  //!param: t - sub interval duration. If can not divide on integer
317  // number of sub-intervals then add leftover to the last
318  // one.
319  //!param: f - black pixel fraction
320  //!param: m - mode
321  //!options: f = 0, m = 0 - returns black pixel occupancy
322  //! m = 1 - set threshold f
323  //! m = 2 - random policy
324  //! m = 0 - random pixel selection
325  virtual double fraction(double=0.,double=0.,int=0);
326 
327  //: calculate rank logarithic significance of wavelet pixels
328  //: reduced wavelet amplitudes are stored in this
329  //: Returns pixel occupancy for significance>0.
330  //!param: n - sub-interval duration in seconds
331  //!param: f - black pixel fraction
332  //!options: f = 0 - returns black pixel occupancy
333  virtual double significance(double, double=1.);
334 
335  //: calculate running rank logarithic significance of wavelet pixels
336  //: reduced wavelet amplitudes are stored in this
337  //: Returns pixel occupancy for significance>0.
338  //!param: n - sub-interval duration in domain units
339  //!param: f - black pixel fraction
340  //!options: f = 0 - returns black pixel occupancy
341  virtual double rsignificance(size_t=0, double=1.);
342 
343  //: calculate running rank logarithic significance of wavelet pixels
344  //: reduced wavelet amplitudes are stored in this
345  //: Returns pixel occupancy for significance>0.
346  //!param: T - sliding window duration in seconds
347  //!param: f - black pixel fraction
348  //!param: t - sliding step in seconds
349  //!options: f = 0 - returns black pixel occupancy
350  //!options: t = 0 - sliding step = wavelet time resolution.
351  virtual double rSignificance(double, double=1., double=0.);
352 
353  //: calculate running logarithic significance of wavelet pixels
354  //: reduced wavelet amplitudes are stored in this
355  //: Returns pixel occupancy for significance>0.
356  //!param: T - sliding window duration in seconds
357  //!param: f - black pixel fraction
358  //!param: t - sliding step in seconds
359  //!options: f = 0 - returns black pixel occupancy
360  //!options: t = 0 - sliding step = wavelet time resolution.
361  virtual double gSignificance(double, double=1., double=0.);
362 
363  //: Selection of a fixed fraction of pixels for each wavelet layer
364  //: Returns fraction of non-zero coefficients.
365  //!param: f - black pixel fraction
366  //!param: m - mode
367  //!options: f = 0 - returns black pixel occupancy
368  //! m = 1 - set threshold f, returns percentile amplitudes
369  //! m =-1 - set threshold f, returns wavelet amplitudes
370  //! m > 1 - random policy,returns percentile amplitudes
371  //! m <-1 - random policy,returns wavelet amplitudes
372  //! m = 0 - random pixel selection
373  //! if m<0 return wavelet amplitudes instead of the percentile amplitude
374  virtual double percentile(double=0.,int=0, WSeries<DataType_t>* = NULL);
375 
376  //: clean up single pixels
377  //!param: S - threshold on pixel significance
378  //!return pixel occupancy.
379  virtual double pixclean(double=0.);
380 
381  //: select pixels from *this which satisfy a coincidence rule
382  //: within specified window
383  //!param: WSeries object used for coincidence
384  //!param: coincidence window in seconds
385  //!return pixel occupancy
386  virtual double coincidence(WSeries<DataType_t> &, int=0, int=0, double=0.);
387 
388  //: select pixels from *this which satisfy a coincidence rule
389  //: within specified window w, above threshold T,
390  //!param: WSeries object used for coincidence
391  //!param: coincidence window in seconds
392  //!param: threshold on significance
393  //!return pixel occupancy
394  virtual double Coincidence(WSeries<DataType_t> &, double=0., double=0.);
395 
396  //: calculate calibration coefficients and apply energy calibration
397  //: to wavelet data in this
398  //: AS_Q calibration: R(f)=(1 + gamma*(R*C-1))/alpha*C
399  //: DARM calibration: R(f)=(1 + gamma*(R*C-1))/gamma*C
400  //: input
401  //!param: number of samples in calibration arrays R & C
402  //!param: frequency resolution
403  //!param: pointer to response function R in Fourier domain
404  //!param: pointer to sensing function C in Fourier domain
405  //!param: time dependent calibration coefficient alpha
406  //!param: time dependent calibration coefficient gamma
407  //!param: 0/1 - AS_Q/DARM_ERR calibration, by default is 0
408  //!return array with calibration constants for each wavelet layer
409  virtual WSeries<double> calibrate(size_t, double,
410  d_complex*, d_complex*,
413  size_t ch=0);
414 
415 
416  //: mask WSeries data with the pixel mask defined in wavecluster object
417  //!param: int n
418  //! if n<0, zero pixels defined in mask (regression)
419  //! if n>=0, zero all pixels except ones defined in the mask
420  //!param: bool - if true, set WSeries data to be positive
421  //! if pMask.size()=0, mask(0,true) is equivalent to abs(data)
422  //!return core pixel occupancy
423  //virtual double mask(int=1, bool=false);
424 
425  //: calculate pixel occupancy for clusters passed selection cuts
426  //!param: true - core pccupancy, false - total occupancy;
427  //!return wavearray<double> with occupancy.
428  //virtual wavearray<double> occupancy(bool=true);
429 
430 
431  // print wseries parameters
432  void print(); // *MENU*
433  virtual void Browse(TBrowser *b) {print();}
434 
435 // data members
436 
437  //: parameters of wavelet transform
439  //: whitening mode
440  size_t w_mode;
441  //: black pixel probability
442  double bpp;
443  //: wavelet zero layer rate
444  double wRate;
445  //: low frequency boundary
446  double f_low;
447  //: high frequency boundary
448  double f_high;
449 
450  ClassDef(WSeries,1)
451 
452 }; // class WSeries<DataType_t>
453 
454 
455 #endif // WSERIES_HH
456 
wavearray< double > t(hp.size())
void mul(WSeries< DataType_t > &)
Definition: wseries.cc:117
void sethigh(double f)
Definition: wseries.hh:114
virtual void resize(unsigned int)
Definition: wseries.cc:883
double f_low
Definition: wseries.hh:446
double f_high
Definition: wseries.hh:448
size_t xsize()
Definition: wseries.hh:128
int layer(double f)
Definition: wseries.cc:149
par[0] value
int getMaxLevel()
Definition: wseries.cc:85
tuple f
Definition: cwb_online.py:91
wavearray< double > a(hp.size())
virtual double percentile(double=0., int=0, WSeries< DataType_t > *=NULL)
param: f - black pixel fraction param: m - mode options: f = 0 - returns black pixel occupancy m = 1 ...
Definition: wseries.cc:2064
virtual WSeries< float > variability(double=0., double=-1., double=-1.)
param: first - time window to calculate normalization constants second - low frequency boundary for c...
Definition: wseries.cc:1278
int n
Definition: cwb_net.C:10
virtual WSeries< DataType_t > & operator+=(WSeries< DataType_t > &)
Definition: wseries.cc:778
void bandpass(wavearray< DataType_t > &ts, double flow, double fhigh, int n=-1)
Definition: wseries.cc:295
size_t maxIndex()
Definition: wseries.hh:131
void putSample(DataType_t a, int n, double m)
Definition: wseries.hh:177
std::vector< int > vector_int
Definition: wseries.hh:23
virtual double rSignificance(double, double=1., double=0.)
param: T - sliding window duration in seconds param: f - black pixel fraction param: t - sliding step...
Definition: wseries.cc:1665
void setlow(double f)
Definition: wseries.hh:107
virtual std::slice getSlice() const
Definition: wavearray.hh:129
std::slice getSlice(double n)
Definition: wseries.hh:134
int m
Definition: cwb_net.C:10
virtual WSeries< DataType_t > & operator*=(WSeries< DataType_t > &)
Definition: wseries.cc:753
virtual double significance(double, double=1.)
param: n - sub-interval duration in seconds param: f - black pixel fraction options: f = 0 - returns ...
Definition: wseries.cc:1469
void setWavelet(const Wavelet &w)
Definition: wseries.cc:216
double wRate
Definition: wseries.hh:444
size_t mode
wavearray< double > w
Definition: Test1.C:27
void wrate(double r)
Definition: wseries.hh:102
int getLevel()
Definition: wseries.hh:91
double getlow() const
Definition: wseries.hh:111
double getbpp() const
Definition: wseries.hh:99
size_t w_mode
Definition: wseries.hh:440
virtual double gSignificance(double, double=1., double=0.)
param: T - sliding window duration in seconds param: f - black pixel fraction param: t - sliding step...
Definition: wseries.cc:1819
virtual wavearray< double > filter(size_t)
param: n - number of decomposition steps algorithm: 1) do forward wavelet transform with n decomposit...
Definition: wseries.cc:1242
int pattern
size_t size() const
Definition: wslice.hh:71
double maxEnergy(wavearray< DataType_t > &ts, Wavelet &w, double=0, int=1, int=0, TH1F *=NULL)
Definition: wseries.cc:486
virtual WSeries< DataType_t > & operator-=(WSeries< DataType_t > &)
Definition: wseries.cc:802
bool isWDM()
Definition: wseries.hh:190
virtual double fraction(double=0., double=0., int=0)
param: t - sub interval duration. If can not divide on integer
Definition: wseries.cc:1356
double fhigh
wavecomplex d_complex
Definition: wseries.hh:22
int getLayer(wavearray< DataType_t > &w, double n)
param: n - layer number
Definition: wseries.cc:175
Definition: Wavelet.hh:31
virtual double coincidence(WSeries< DataType_t > &, int=0, int=0, double=0.)
param: WSeries object used for coincidence param: coincidence window in seconds return pixel occupanc...
Definition: wseries.cc:911
void setLevel(size_t n)
Definition: wseries.hh:94
double gethigh() const
Definition: wseries.hh:118
regression r
Definition: Regression_H1.C:44
double flow
virtual ~WSeries()
Definition: wseries.cc:76
virtual WSeries< DataType_t > & operator[](const std::slice &)
Definition: wseries.cc:713
void wavescan(WSeries< DataType_t > **, int, TH1F *=NULL)
Definition: wseries.cc:604
wavearray< double > ts(N)
virtual void Browse(TBrowser *b)
Definition: wseries.hh:433
WSeries()
Definition: wseries.cc:23
virtual void resample(double, int=6)
Definition: wseries.cc:897
void print()
param: int n if n<0, zero pixels defined in mask (regression) if n>=0, zero all pixels except ones de...
Definition: wseries.cc:2334
double resolution(int=0)
Definition: wseries.hh:137
double Gamma2Gauss(TH1F *=NULL)
Definition: wseries.cc:558
virtual double Coincidence(WSeries< DataType_t > &, double=0., double=0.)
param: WSeries object used for coincidence param: coincidence window in seconds param: threshold on s...
Definition: wseries.cc:1002
virtual void lprFilter(double, int=0, double=0., double=0.)
Definition: wseries.cc:1108
void Forward(int n=-1)
param: wavelet - n is number of steps (-1 means full decomposition)
Definition: wseries.cc:228
Meyer< double > S(1024, 2)
int l
Definition: cbc_plots.C:434
virtual void median(double t, bool norm=false)
Definition: wseries.cc:1091
WSeries< DataType_t > & operator=(const wavearray< DataType_t > &)
Definition: wseries.cc:684
void setbpp(double f)
Definition: wseries.hh:97
virtual double pixclean(double=0.)
param: S - threshold on pixel significance return pixel occupancy.
Definition: wseries.cc:1965
DataType_t * data
Definition: wavearray.hh:301
DataType_t getSample(int n, double m)
Definition: wseries.hh:167
WaveDWT< DataType_t > * pWavelet
Definition: wseries.hh:438
double wdmPacket(int pattern, char opt='L', TH1F *=NULL)
Definition: wseries.cc:358
double wrate() const
Definition: wseries.hh:104
size_t stride() const
Definition: wslice.hh:75
void Inverse(int n=-1)
param: n - number of steps (-1 means full reconstruction)
Definition: wseries.cc:273
size_t sizeZero()
Definition: wseries.hh:126
void putLayer(wavearray< DataType_t > &, double n)
param: n - layer number
Definition: wseries.cc:201
virtual void Dump(const char *, int=0)
Definition: wseries.cc:849
double frequency(int l)
Definition: wseries.cc:99
size_t start() const
Definition: wslice.hh:67
int maxLayer()
Definition: wseries.hh:121
virtual double rsignificance(size_t=0, double=1.)
param: n - sub-interval duration in domain units param: f - black pixel fraction options: f = 0 - ret...
Definition: wseries.cc:1572
double bpp
Definition: wseries.hh:442
virtual WSeries< double > calibrate(size_t, double, d_complex *, d_complex *, wavearray< double > &, wavearray< double > &, size_t ch=0)
param: number of samples in calibration arrays R & C param: frequency resolution param: pointer to re...
Definition: wseries.cc:2189
virtual WSeries< double > white(double, int, double=0., double=0.)
what it does: each wavelet layer is devided into k intervals.
Definition: wseries.cc:1128