Logo Coherent WaveBurst  
Reference Guide
Logo
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
wavearray.cc
Go to the documentation of this file.
1 /*-------------------------------------------------------
2  * Package: Wavelet Analysis Tool
3  * generic data container to use with DMT and ROOT
4  * File name: wavearray.cc
5  *-------------------------------------------------------
6 */
7 // wavearray class is the base class for wavelet data amd methods .
8 
9 #include <time.h>
10 #include <iostream>
11 //#include <sstream>
12 //#include <strstream>
13 
14 //#include "PConfig.h"
15 #ifdef __GNU_STDC_OLD
16 #include <gnusstream.h>
17 #else
18 #include <sstream>
19 #endif
20 
21 #include "wavearray.hh"
22 #include "wavefft.hh"
23 #include "TComplex.h"
24 #include "TObjArray.h"
25 #include "TObjString.h"
26 #include "TFile.h"
27 
28 extern "C" {
29  typedef int (*qsort_func) (const void*, const void*);
30 }
31 
32 ClassImp(wavearray<double>)
33 
34 using namespace std;
35 
36 //: Default constructor
37 template<class DataType_t>
38 wavearray<DataType_t>::wavearray() :
39  data(NULL),Size(0),Rate(1.),Start(0.),Stop(0.),Edge(0.),fftw(NULL), ifftw(NULL)
40 {
41  Slice = std::slice(0,0,0);
42 }
43 
44 //: allocates a data array with N=n elements.
45 template<class DataType_t>
47 Rate(1.),Start(0.),Edge(0.),fftw(NULL),ifftw(NULL)
48 {
49  if (n <= 0 ) n = 1;
50  data = (DataType_t *)malloc(n*sizeof(DataType_t));
51  Size = n;
52  Stop = double(n);
53  Slice = std::slice(0,n,1);
54  *this = 0;
55 }
56 
57 //: copy constructor
58 template<class DataType_t>
60  data(NULL),Size(0),Rate(1.),Start(0.),Stop(0.),Edge(0.),fftw(NULL), ifftw(NULL)
61 { *this = a; }
62 
63 // explicit construction from array
64 template<class DataType_t> template<class T>
66 wavearray(const T *p, unsigned int n, double r) :
67  data(NULL),Size(0),Rate(1.),Start(0.),Edge(0.),fftw(NULL), ifftw(NULL)
68 {
69  unsigned int i;
70  if(n != 0 && p != NULL){
71  data = (DataType_t *)malloc(n*sizeof(DataType_t));
72  for (i=0; i < n; i++) data[i] = p[i];
73  Size = n;
74  Rate = r;
75  Stop = n/r;
76  }
77  Slice = std::slice(0,n,1);
78 }
79 
80 //: destructor
81 template<class DataType_t>
83 {
84  if(data) free(data);
85  resetFFTW();
86 }
87 
88 //: operators =
89 
90 template<class DataType_t>
93 {
94  if (this==&a) return *this;
95 
96  unsigned int i;
97  unsigned int N = a.Slice.size();
98  unsigned int m = a.Slice.stride();
99 
100  this->Edge = a.Edge;
101  if (N>0 && a.data) {
102  register DataType_t *pm = a.data + a.Slice.start();
104 
105  for (i=0; i < N; i++) { data[i] = *pm; pm+=m; }
106 
107  if(a.rate()>0.) {
108  start(a.start() + a.Slice.start()/a.rate());
109  stop(start()+N*m/a.rate());
110  }
111  else {
112  start(a.start());
113  stop(a.stop());
114 
115  }
116 
117  rate(a.rate()/m);
118  Slice = std::slice(0, size(),1);
119  const_cast<wavearray<DataType_t>&>(a).Slice = std::slice(0,a.size(),1);
120  }
121  else if(data) {
122  free(data);
123  data = NULL;
124  Size = 0;
125  Start = a.start();
126  Stop = a.stop();
127  Rate = a.rate();
128  Slice = std::slice(0,0,0);
129  }
130  this->fftw = NULL;
131  this->ifftw = NULL;
132 
133  return *this;
134 }
135 
136 template<class DataType_t>
139 {
140  unsigned int i;
141  unsigned int N = limit(a);
142  unsigned int n = Slice.stride();
143  unsigned int m = a.Slice.stride();
144  register DataType_t *p = a.data + a.Slice.start();
145 
146  if(size())
147  for (i=Slice.start(); i<N; i+=n){ data[i] = *p; p += m; }
148 
149  Slice = std::slice(0, size(),1);
150  a.Slice = std::slice(0,a.size(),1);
151  return *this;
152 }
153 
154 template<class DataType_t>
156 operator=(const DataType_t c)
157 {
158  unsigned int i;
159  unsigned int n = Slice.stride();
160  unsigned int N = limit();
161 
162  if(size())
163  for (i=Slice.start(); i < N; i+=n) data[i] = c;
164 
165  Slice = std::slice(0, size(),1);
166  return *this;
167 }
168 
169 //: operators +=
170 
171 template<class DataType_t>
174 {
175  unsigned int i;
176  unsigned int N = limit(a);
177  unsigned int n = Slice.stride();
178  unsigned int m = a.Slice.stride();
179  register DataType_t *p = a.data + a.Slice.start();
180 
181  if(size())
182  for (i=Slice.start(); i<N; i+=n){ data[i] += *p; p += m; }
183 
184  Slice = std::slice(0, size(),1);
185  a.Slice = std::slice(0,a.size(),1);
186  return *this;
187 }
188 
189 template<class DataType_t>
191 operator+=(const DataType_t c)
192 {
193  unsigned int i;
194  unsigned int n = Slice.stride();
195  unsigned int N = limit();
196 
197  if(size())
198  for (i=Slice.start(); i < N; i+=n) data[i] += c;
199 
200  Slice = std::slice(0, size(),1);
201  return *this;
202 }
203 
204 //: operators -=
205 
206 template<class DataType_t>
209 {
210  unsigned int i;
211  unsigned int N = limit(a);
212  unsigned int n = Slice.stride();
213  unsigned int m = a.Slice.stride();
214  register DataType_t *p = a.data + a.Slice.start();
215 
216  if(size())
217  for (i=Slice.start(); i<N; i+=n){ data[i] -= *p; p += m; }
218 
219  Slice = std::slice(0, size(),1);
220  a.Slice = std::slice(0,a.size(),1);
221  return *this;
222 }
223 
224 template<class DataType_t>
226 operator-=(const DataType_t c)
227 {
228  unsigned int i;
229  unsigned int n = Slice.stride();
230  unsigned int N = limit();
231 
232  if(size())
233  for (i=Slice.start(); i < N; i+=n) data[i] -= c;
234 
235  Slice = std::slice(0, size(),1);
236  return *this;
237 }
238 
239 //: multiply all elements of data array by constant
240 template<class DataType_t>
242 operator*=(const DataType_t c)
243 {
244  unsigned int i;
245  unsigned int n = Slice.stride();
246  unsigned int N = limit();
247 
248  if(size())
249  for (i=Slice.start(); i < N; i+=n) data[i] *= c;
250 
251  Slice = std::slice(0, size(),1);
252  return *this;
253 }
254 
255 // scalar production
256 template<class DataType_t>
259 {
260  unsigned int i;
261  unsigned int N = limit(a);
262  unsigned int n = Slice.stride();
263  unsigned int m = a.Slice.stride();
264  register DataType_t *p = a.data + a.Slice.start();
265 
266  if(size())
267  for (i=Slice.start(); i<N; i+=n){ data[i] *= *p; p += m; }
268 
269  Slice = std::slice(0, size(),1);
270  a.Slice = std::slice(0,a.size(),1);
271  return *this;
272 }
273 
274 // operator[](const std::slice &)
275 template<class DataType_t>
278 {
279  Slice = s;
280  if(limit() > size()){
281  cout << "wavearray::operator[slice]: Illegal argument "<<limit()<<" "<<size()<<"\n";
282  Slice = std::slice(0,size(),1);
283  }
284  return *this;
285 }
286 
287 // operator[](const std::slice &)
288 template<class DataType_t>
289 DataType_t & wavearray<DataType_t>::
290 operator[](const unsigned int n)
291 {
292  if(n >= size()){
293  cout << "wavearray::operator[int]: Illegal argument\n";
294  return data[0];
295  }
296  return data[n];
297 }
298 
299 //: Dumps data array to file *fname in ASCII format.
300 template<class DataType_t>
301 void wavearray<DataType_t>::Dump(const char *fname, int app)
302 {
303  int n=size();
304  char mode[3] = "w";
305  if (app == 1) strcpy (mode, "a");
306 
307  FILE *fp;
308 
309  if ( (fp = fopen(fname, mode)) == NULL ) {
310  cout << " Dump() error: cannot open file " << fname <<". \n";
311  return;
312  };
313 
314  if(app == 0) {
315  fprintf( fp,"# start time: -start %lf \n", this->Start );
316  fprintf( fp,"# sampling rate: -rate %lf \n", this->Rate );
317  fprintf( fp,"# number of samples: -size %d \n", (int)this->Size );
318  }
319 
320  if(app == 2) {
321  double dt = 1./this->rate();
322  for (int i = 0; i < n; i++) {
323  double time = this->start()+i*dt;
324  fprintf( fp,"%14f %e\n", time, data[i]);
325  }
326  } else {
327  for (int i = 0; i < n; i++) fprintf( fp,"%e ", (float)data[i]);
328  }
329 
330  fclose(fp);
331 }
332 
333 //: Dumps data array to file *fname in binary format.
334 template<class DataType_t>
335 void wavearray<DataType_t>::DumpBinary(const char *fname, int app)
336 {
337  int n = size() * sizeof(DataType_t);
338  char mode[5];
339  strcpy (mode, "w");
340 
341  if (app == 1) strcpy (mode, "a");
342 
343  FILE *fp;
344 
345  if ( (fp=fopen(fname, mode)) == NULL ) {
346  cout << " DumpBinary() error : cannot open file " << fname <<" \n";
347  return ;
348  }
349 
350  fwrite(data, n, 1, fp);
351  fclose(fp);
352 }
353 
354 //: Dumps data array to file *fname in binary format as "short".
355 template<class DataType_t>
356 void wavearray<DataType_t>::DumpShort(const char *fname, int app)
357 {
358  int n = size();
359  char mode[5] = "w";
360  if (app == 1) strcpy (mode, "a");
361 
362  FILE *fp;
363  if ( (fp = fopen(fname, mode)) == NULL ) {
364  cout << " DumpShort() error : cannot open file " << fname <<". \n";
365  return;
366  }
367 
368  short *dtemp;
369  dtemp=new short[n];
370 
371  for ( int i=0; i<n; i++ ) dtemp[i]=short(data[i]) ;
372 
373  n = n * sizeof(short);
374 
375  fwrite(dtemp, n, 1, fp);
376  fclose(fp);
377  delete [] dtemp;
378 }
379 
380 //: Dumps object wavearray to file *fname in root format.
381 template<class DataType_t>
383 {
384  TFile* rfile = new TFile(fname, "RECREATE");
385  this->Write("wavearray"); // write this object
386  rfile->Close();
387  delete rfile;
388 }
389 
390 //: Read data from file in binary format.
391 template<class DataType_t>
393 {
394  int step = sizeof(DataType_t);
395  FILE *fp;
396  double d;
397 
398  if ( (fp=fopen(fname,"r")) == NULL ) {
399  cout << " ReadBinary() error : cannot open file " << fname <<". \n";
400  exit(1);
401  }
402 
403 
404  if(N == 0){ // find the data length
405  while(!feof(fp)){
406  fread(&d,step,1,fp);
407  N++;
408  }
409  N--;
410  rewind(fp);
411  }
412 
413  if(int(size()) != N) resize(N);
414  int n = size() * sizeof(DataType_t);
415 
416  fread(data, n, 1, fp); // Reading binary record
417  fclose(fp);
418 }
419 
420 //: Read data from file as short.
421 template<class DataType_t>
423 {
424  short *dtmp;
425  dtmp = new short[size()];
426  int n = size() * sizeof(short);
427  FILE *fp;
428 
429  if ( (fp=fopen(fname,"r")) == NULL ) {
430  cout << " ReadShort() error : cannot open file " << fname <<". \n";
431  return;
432  }
433 
434  cout << " Reading binary record, size="<< n <<"\n";
435 
436  fread(dtmp, n, 1, fp);
437  for (unsigned int i = 0; i < size(); i++)
438  data[i] = DataType_t(dtmp[i]);
439  fclose(fp);
440  delete [] dtmp;
441 }
442 
443 //: resizes data array to a new length n.
444 template<class DataType_t>
446 {
447  DataType_t *p = NULL;
448  if(n==0){
449  if(data) free(data);
450  data = NULL;
451  Size = 0;
452  Stop = Start;
453  Slice = std::slice(0,0,0);
454  return;
455  }
456 
457  p = (DataType_t *)realloc(data,n*sizeof(DataType_t));
458 
459  if(p){
460  data = p;
461  Size = n;
462  Stop = Start + n/Rate;
463  Slice = std::slice(0,n,1);
464  }
465  else {
466  cout<<"wavearray::resize(): memory allocation failed.\n";
467  exit(0);
468  }
469 
470  resetFFTW();
471 }
472 
473 /************************************************************************
474  * Creates new data set by resampling the original data from "a" *
475  * with new sample frequency "f". Uses polynomial interpolation scheme *
476  * (Lagrange formula) with np-points, "np" must be even number, by *
477  * default np=6 (corresponds to 5-th order polynomial interpolation). *
478  * This function calls wavearray::resize() function to adjust *
479  * current data array if it's size does not exactly match new number *
480  * of points. *
481  ************************************************************************/
482 
483 template<class DataType_t>
485 resample(wavearray<DataType_t> const &a, double f, int nF)
486 {
487  int nP = nF;
488  if(nP<=1) nP = 6;
489  if(int(a.size())<nP) nP = a.size();
490  nP = (nP>>1)<<1;
491 
492  int N;
493  int nP2 = nP/2;
494 
495  const DataType_t *p = a.data;
496 
497  register int i;
498  register int iL;
499  register double x;
500  double *temp=new double[nF];
501 
502  rate(f);
503  double ratio = a.rate() / rate();
504  N = int(a.size()/ratio + 0.5);
505 
506  if ( (int)size() != N ) resize(N);
507 
508 // left border
509 
510  int nL = int(nP2/ratio);
511 
512  for (i = 0; i < nL; i++)
513  data[i] = (DataType_t)Nevill(i*ratio, nP, p, temp);
514 
515 // in the middle of array
516 
517  int nM = int((a.size()-nP2)/ratio);
518  if(nM < nL) nM = nL;
519 
520  if(nM&1 && nM>nL) {
521  x = nL*ratio;
522  iL = int(x) - nP2 + 1 ;
523  data[i] = (DataType_t)Nevill(x-iL, nP, p+iL, temp);
524  nL++;
525  }
526  for (i = nL; i < nM; i+=2) {
527  x = i*ratio;
528  iL = int(x) - nP2 + 1 ;
529  data[i] = (DataType_t)Nevill(x-iL, nP, p+iL, temp);
530  x += ratio;
531  iL = int(x) - nP2 + 1 ;
532  data[i+1] = (DataType_t)Nevill(x-iL, nP, p+iL, temp);
533  }
534 
535 // right border
536 
537  int nR = a.size() - nP;
538  p += nR;
539  for (i = nM; i < N; i++)
540  data[i] = (DataType_t)Nevill(i*ratio-nR, nP, p, temp);
541 
542  delete [] temp;
543 }
544 
545 template<class DataType_t>
546 void wavearray<DataType_t>::resample(double f, int nF)
547 {
549  a = *this;
550  resample(a,f,nF);
551 }
552 
553 template<class DataType_t>
555 {
556  if(f<=0) {
557  cout<<"wavearray::Resample(): error - resample frequency must be >0\n";
558  exit(1);
559  }
560 
561  double rsize = f/this->rate()*this->size(); // new size
562 
563  if(fmod(rsize,2)!=0) {
564  cout<<"wavearray::Resample(): error - resample frequency (" << f << ") not allowed\n";
565  exit(1);
566  }
567 
568  int size = this->size(); // old size
569 
570  this->FFTW(1);
571  this->resize((int)rsize);
572  for(int i=size;i<rsize;i++) this->data[i]=0;
573  this->FFTW(-1);
574  this->rate(f);
575 }
576 
577 template<class DataType_t>
579 {
580  if(T==0) return;
581 
582  // search begin,end of non zero data
583  int ibeg=0; int iend=0;
584  for(int i=0;i<(int)this->size();i++) {
585  if(this->data[i]!=0 && ibeg==0) ibeg=i;
586  if(this->data[i]!=0) iend=i;
587  }
588  int ilen=iend-ibeg+1;
589  // create temporary array for FFTW & add scratch buffer + T
590  int ishift = fabs(T)*this->rate();
591  int isize = 2*ilen+2*ishift;
592  isize = isize + (isize%4 ? 4 - isize%4 : 0); // force to be multiple of 4
593  wavearray<DataType_t> w(isize);
594  w.rate(this->rate()); w=0;
595  // copy this->data !=0 in the middle of w array & set this->data=0
596  for(int i=0;i<ilen;i++) {w[i+ishift+ilen/2]=this->data[ibeg+i];this->data[ibeg+i]=0;}
597 
598  double pi = TMath::Pi();
599  // apply time shift to waveform vector
600  w.FFTW(1);
601  TComplex C;
602  double df = w.rate()/w.size();
603  //cout << "T : " << T << endl;
604  for (int ii=0;ii<(int)w.size()/2;ii++) {
605  TComplex X(w[2*ii],w[2*ii+1]);
606  X=X*C.Exp(TComplex(0.,-2*pi*ii*df*T)); // Time Shift
607  w[2*ii]=X.Re();
608  w[2*ii+1]=X.Im();
609  }
610  w.FFTW(-1);
611 
612  // copy shifted data to input this->data array
613  for(int i=0;i<(int)w.size();i++) {
614  int j=ibeg-(ishift+ilen/2)+i;
615  if((j>=0)&&(j<(int)this->size())) this->data[j]=w[i];
616  }
617 }
618 
619 
620 template<class DataType_t>
622 Resample(wavearray<DataType_t> const &a, double f, int np)
623 {
624  int i1, N, n1, n2, n, np2 = np/2;
625  double s, x, *c, *v;
626  c=new double[np];
627  v=new double[np];
628 
629  rate(f);
630  double ratio = a.rate() / rate();
631  n = a.size();
632  N = int(n/ratio + 0.5);
633 
634  if ( (int)size() != N ) resize(N);
635 
636 // calculate constant filter part c(k) = -(-1)^k/k!/(np-k-1)!
637  for (int i = 0; i < np; i++) {
638  int m = 1;
639  for (int j = 0; j < np; j++) if (j != i) m *= (i - j);
640  c[i] = 1./double(m);
641  }
642 
643  for (int i = 0; i < N; i++)
644  {
645  x = i*ratio;
646 
647  i1 = int(x);
648  x = x - i1 + np2 - 1;
649 
650 // to treat data boundaries we need to calculate critical numbers n1 and n2
651  n1 = i1 - np2 + 1;
652  n2 = i1 + np2 + 1 - n;
653 
654 // here we calculate the sum of products like h(k)*a(k), k=0..np-1,
655 // where h(k) are filter coefficients, a(k) are signal samples
656 // h(k) = c(k)*prod (x-i), i != k, c(k) = -(-1)^k/k!/(np-k-1)!
657 
658 // get signal part multiplied by constant filter coefficients
659  if ( n1 >= 0 )
660  if ( n2 <= 0 )
661 // regular case - far from boundaries
662  for (int j = 0; j < np; j++) v[j] = c[j]*a.data[i1 + j - np2 + 1];
663 
664  else {
665 // right border case
666  x = x + n2;
667  for (int j = 0; j < np; j++) v[j] = c[j]*a.data[n + j - np];
668  }
669  else {
670 // left border case
671  x = x + n1;
672  for (int j = 0; j < np; j++) v[j] = c[j]*a.data[j];
673  }
674 
675 // multiply resulted v[j] by x-dependent factors (x-k), if (k != j)
676  for (int k = 0; k < np; k++) {
677  for (int j = 0; j < np; j++) if (j != k) v[j] *= x;
678  x -= 1.; }
679 
680  s = 0.;
681  for (int j = 0; j < np; j++) s += v[j];
682 
683  data[i] = (DataType_t)s;
684  }
685 
686  delete [] c;
687  delete [] v;
688 }
689 
690 /******************************************************************************
691  * copies data from data array of the object wavearray a to *this
692  * object data array. Procedure starts from the element "a_pos" in
693  * the source array which is copied to the element "pos" in the
694  * destination array. "length" is the number of data elements to copy.
695  * By default "length"=0, which means copy all source data or if destination
696  * array is shorter, copy data until the destination is full.
697  *****************************************************************************/
698 template<class DataType_t> void wavearray<DataType_t>::
699 cpf(const wavearray<DataType_t> &a, int length, int a_pos, int pos)
700 {
701 // if (rate() != a.rate()){
702 // cout << "wavearray::cpf() warning: sample rate mismatch.\n";
703 // cout<<"rate out: "<<rate()<<" rate in: "<<a.rate()<<endl;
704 // }
705 
706  if (length == 0 )
707  length = ((size() - pos) < (a.size() - a_pos))?
708  (size() - pos) : (a.size() - a_pos);
709 
710  if( length > (int)(size() - pos) ) length = size() - pos;
711  if( length > (int)(a.size() - a_pos) ) length = a.size() - a_pos;
712 
713  for (int i = 0; i < length; i++)
714  data[i + pos] = a.data[i + a_pos];
715 
716 // rate(a.rate());
717 }
718 
719 /******************************************************************************
720  * Adds data from data array of the object wavearray a to *this
721  * object data array. Procedure starts from the element "a_pos" in
722  * the source array which is added starting from the element "pos" in the
723  * destination array. "length" is the number of data elements to add.
724  * By default "length"=0, which means add all source data or if destination
725  * array is shorter, add data up to the end of destination.
726  *****************************************************************************/
727 template<class DataType_t> void wavearray<DataType_t>::
728 add(const wavearray<DataType_t> &a, int length, int a_pos, int pos)
729 {
730 // if (rate() != a.rate())
731 // cout << "wavearray::add() warning: sample rate mismatch.\n";
732 
733  if (length == 0 )
734  length = ((size() - pos) < (a.size() - a_pos))?
735  (size() - pos) : (a.size() - a_pos);
736 
737  if( length > (int)( size()- pos) ) length = size() - pos;
738  if( length > (int)(a.size() - a_pos) ) length = a.size() - a_pos;
739 
740  for (int i = 0; i < length; i++)
741  data[i + pos] += a.data[i + a_pos];
742 }
743 
744 
745 /******************************************************************************
746  * Subtracts data array of the object wavearray a from *this
747  * object data array. Procedure starts from the element "a_pos" in
748  * the source array which is subtracted starting from the element "pos" in
749  * the destination array. "length" is number of data elements to subtract.
750  * By default "length"=0, which means subtract all source data or if the
751  * destination array is shorter, subtract data up to the end of destination.
752  ******************************************************************************/
753 template<class DataType_t> void wavearray<DataType_t>::
754 sub(const wavearray<DataType_t> &a, int length, int a_pos, int pos)
755 {
756 // if (rate() != a.rate())
757 // cout << "wavearray::sub() warning: sample rate mismatch.\n";
758 
759  if ( length == 0 )
760  length = ((size() - pos) < (a.size() - a_pos))?
761  (size() - pos) : (a.size() - a_pos);
762 
763  if( length > (int)(size() - pos) ) length = size() - pos;
764  if( length > (int)(a.size() - a_pos) ) length = a.size() - a_pos;
765 
766  for (int i = 0; i < length; i++)
767  data[i + pos] -= a.data[i + a_pos];
768 }
769 
770 /*****************************************************************
771  * append two wavearrays with the same rate ignoring start time
772  * of input array.
773  *****************************************************************/
774 template<class DataType_t> size_t wavearray<DataType_t>::
776 {
777  size_t n = this->size();
778  size_t m = a.size();
779 
780  if (this->rate() != a.rate())
781  cout << "wavearray::append() warning: sample rate mismatch.\n";
782 
783  if(m == 0 ) return this->size();
784  this->resize(n+m);
785  this->cpf(a,m,0,n);
786  this->stop(start()+(n+m/rate()));
787 
788  return n+m;
789 }
790 
791 /*****************************************************************
792  * append a data point
793  *****************************************************************/
794 template<class DataType_t> size_t wavearray<DataType_t>::
795 append(DataType_t a)
796 {
797  size_t n = this->size();
798  this->resize(n+1);
799  this->data[n] = a;
800  this->Stop += 1./rate();
801  return n+1;
802 }
803 
804 
805 /*****************************************************************
806  * Calculates Fourier Transform for real signal using
807  * Fast Fourier Transform (FFT) algorithm. Packs resulting
808  * complex data into original array of real numbers.
809  * Calls wavefft() function which is capable to deal with any
810  * number of data points. FFT(1) means forward transformation,
811  * which is default, FFT(-1) means inverse transformation.
812  *****************************************************************/
813 template<class DataType_t>
814 void wavearray<DataType_t>::FFT(int direction)
815 {
816  double *a, *b;
817  int N = size();
818  int n2 = N/2;
819 
820  a=new double[N];
821  b=new double[N];
822 
823  switch (direction)
824  { case 1:
825 
826 // direct transform
827 
828  for (int i=0; i<N; i++) { a[i] = data[i]; b[i]=0.;}
829 
830  wavefft(a, b, N, N, N, -1);
831 
832 // pack complex numbers to real array
833  for (int i=0; i<n2; i++)
834  { data[2*i] = (DataType_t)a[i]/N;
835  data[2*i+1] = (DataType_t)b[i]/N;
836  }
837 
838 // data[1] is not occupied because imag(F[0]) is always zero and we
839 // store in data[1] the value of F[N/2] which is pure real for even "N"
840  data[1] = (DataType_t)a[n2]/N;
841 
842 // in the case of odd number of data points we store imag(F[N/2])
843 // in the last element of array data[N]
844  if ((N&1) == 1) data[N-1] = (DataType_t)b[n2]/N;
845 
846  break;
847 
848  case -1:
849 // inverse transform
850 
851 // unpack complex numbers from real array
852  for (int i=1;i<n2;i++)
853  { a[i]=data[2*i];
854  b[i]=data[2*i+1];
855  a[N-i]=data[2*i];
856  b[N-i]=-data[2*i+1];
857  }
858 
859  a[0]=data[0];
860  b[0]=0.;
861 
862  if ((N&1) == 1)
863  { a[n2]=data[1]; b[n2]=data[N-1]; } // for odd n
864  else
865  { a[n2]=data[1]; b[n2]=0.; }
866 
867  wavefft(a, b, N, N, N, 1); // call FFT for inverse tranform
868 
869  for (int i=0; i<N; i++)
870  data[i] = (DataType_t)a[i]; // copy the result from array "a"
871  }
872 
873  delete [] b;
874  delete [] a;
875 }
876 
877 template<class DataType_t>
878 void wavearray<DataType_t>::FFTW(int direction)
879 {
880  double *a, *b;
881  int N = size();
882  int n2 = N/2;
883 
884  a=new double[N];
885  b=new double[N];
886 
887  int sign=0,kind=0; // fftw dummy parameters
888 
889  switch (direction)
890  { case 1:
891 
892  if(fftw==NULL) {
893  fftw = new TFFTRealComplex(1,&N,false);
894  fftw->Init("LS",sign,&kind); // EX - optimized, ES - least optimized
895  }
896 
897 // direct transform
898 
899  for (int i=0; i<N; i++) { a[i] = data[i]; b[i]=0.;}
900 
901  fftw->SetPoints(a);
902  fftw->Transform();
903  fftw->GetPointsComplex(a, b);
904 
905 // pack complex numbers to real array
906  for (int i=0; i<n2; i++)
907  { data[2*i] = (DataType_t)a[i]/N;
908  data[2*i+1] = (DataType_t)b[i]/N;
909  }
910 
911 // data[1] is not occupied because imag(F[0]) is always zero and we
912 // store in data[1] the value of F[N/2] which is pure real for even "N"
913  data[1] = (DataType_t)a[n2]/N;
914 
915 // in the case of odd number of data points we store imag(F[N/2])
916 // in the last element of array data[N]
917  if ((N&1) == 1) data[N-1] = (DataType_t)b[n2]/N;
918 
919  break;
920 
921  case -1:
922 // inverse transform
923 
924  if(ifftw==NULL) {
925  ifftw = new TFFTComplexReal(1,&N,false);
926  ifftw->Init("LS",sign,&kind); // EX - optimized, ES - least optimized
927  }
928 
929 // unpack complex numbers from real array
930  for (int i=1;i<n2;i++)
931  { a[i]=data[2*i];
932  b[i]=data[2*i+1];
933  a[N-i]=data[2*i];
934  b[N-i]=-data[2*i+1];
935  }
936 
937  a[0]=data[0];
938  b[0]=0.;
939 
940  if ((N&1) == 1)
941  { a[n2]=data[1]; b[n2]=data[N-1]; } // for odd n
942  else
943  { a[n2]=data[1]; b[n2]=0.; }
944 
945  ifftw->SetPointsComplex(a,b);
946  ifftw->Transform();
947  ifftw->GetPoints(a);
948 
949  for (int i=0; i<N; i++)
950  data[i] = (DataType_t)a[i]; // copy the result from array "a"
951  }
952 
953  delete [] b;
954  delete [] a;
955 }
956 
957 // Release FFTW memory
958 template<class DataType_t>
960 {
961  if(fftw) {delete fftw;fftw=NULL;}
962  if(ifftw) {delete ifftw;ifftw=NULL;}
963 }
964 
965 /*************************************************************************
966  * Stack generates wavearray *this by stacking data from wavearray td.
967  * Input data are devided on subsets with with samples "length"
968  * then they are added together. Average over the all subsets is saved in
969  * *this.
970  *************************************************************************/
971 template<class DataType_t>
974 {
975  double ss, s0, s2;
976  rate(td.rate());
977  int k = td.size()/length;
978  int n = k*length;
979 
980  if (k == 0) {
981  cout <<" Stack() error: data length too short to contain \n"
982  << length << " samples\n";
983  return 0.;
984  }
985 
986  if (size() != (unsigned int)length) resize(length);
987 
988 // sum (stack) all k periods of frequency f to produce 1 cycle
989  s0 = 0.;
990  for (int i = 0; i < length; i++) {
991  ss = 0.;
992  for (int j = i; j < n; j += length) ss += td.data[j];
993  data[i] = (DataType_t)ss/k;
994  s0 += ss;
995  }
996  s0 /= (k*length);
997 
998 // remove constant displacement (zero frequency component)
999  s2 = 0.;
1000  for (int i = 0; i < length; i++) {
1001  data[i] -= (DataType_t)s0;
1002  s2 += data[i]*data[i];
1003  }
1004  s2 /= length;
1005 
1006  return s2; // return stacked signal power (energy/N)
1007 }
1008 
1009 /*************************************************************************
1010  * Another version of Stack:
1011  * Input data (starting from sample "start") are devided on "k" subsets
1012  * with sections "length"
1013  * Average over the all sections is saved in *this.
1014  *************************************************************************/
1015 template<class DataType_t>
1018 {
1019  double avr, rms;
1020  rate(td.rate());
1021 
1022  if(start+length > (int)td.size()) length = td.size()-start;
1023 
1024  int k = (size()<1) ? 0 : length/size();
1025  if (k == 0) {
1026  cout <<" Stack() error: data length too short to contain \n"
1027  << length << " samples\n";
1028  return 0.;
1029  }
1030 
1031  *this = 0;
1032  for (int i = 0; i<k; i++) add(td, size(), start+i*size());
1033  *this *= DataType_t(1./k);
1034  getStatistics(avr,rms);
1035  *this -= DataType_t(avr);
1036 
1037  return rms*rms; // return stacked signal power
1038 }
1039 
1040 // Stack generates wavearray *this by stacking data from wavearray td.
1041 // Input data are devided on subsets with with samples "length"
1042 // then they are added together. Average over the all subsets is saved in
1043 // *this.
1044 template<class DataType_t>
1046 Stack(const wavearray<DataType_t> &td, double window)
1047 {
1048  return this->Stack(td, int(td.rate() * window));
1049 }
1050 
1051 // wavearray mean
1052 template<class DataType_t>
1054 {
1055  register size_t i;
1056  register double x = 0;
1057 
1058  if(!size()) return 0.;
1059 
1060  for(i=0; i<size(); i++) x += data[i];
1061  return x/size();
1062 }
1063 
1064 template<class DataType_t>
1066 {
1067 // return f percentile mean for data
1068 // param - f>0 positive side pecentile, f<0 - double-side percentile
1069 
1070  if(!size()) return 0;
1071  double ff = fabs(f);
1072  size_t nn = size_t(this->Edge*rate());
1073  if(ff > 1) ff = 1.;
1074  if(!nn || 2*nn>=size()-2) {return mean();}
1075  if(ff==0.) return median(nn,size()-nn-1);
1076 
1077  size_t i;
1078  size_t N = size()-2*nn;
1079  size_t mL = f<0. ? (N-int(N*ff))/2 : 0; // left boundary
1080  size_t mR = f<0. ? (N+int(N*ff))/2-1 : int(N*ff)-1; // right boundary
1081 
1082  double x = 0.;
1083 
1084  wavearray<DataType_t> wtmp = *this;
1085  if(f>0) wtmp *= wtmp;
1086  DataType_t *p = wtmp.data;
1087  DataType_t **pp = (DataType_t **)malloc(N*sizeof(DataType_t*));
1088  for(i=nn; i<N+nn; i++) pp[i-nn] = p + i;
1089 
1090  if(mL>0) waveSplit(pp,0,N-1,mL); // split on left side
1091  if(mR<N-2) waveSplit(pp,0,N-1,mR); // split on right side
1092  for(i=mL; i<=mR; i++) x += this->data[pp[i]-p];
1093 
1094  free(pp);
1095  return x/(mR-mL+1.);
1096 }
1097 
1098 template<class DataType_t>
1100 {
1101  register size_t i;
1102  register double x = 0.;
1103  register DataType_t *p = data + s.start();
1104  size_t N = s.size();
1105  size_t m = (s.stride()<=0) ? 1 : s.stride();
1106  if(size()<limit(s)) N = (limit(s) - s.start() - 1)/m;
1107  for(i=0; i<N; i++) { x += *p; p += m; }
1108  return (N==0) ? 0. : x/N;
1109 }
1110 
1111 // running mean
1112 template<class DataType_t>
1115  bool clean,
1116  size_t skip)
1117 {
1118 // calculate running mean of data (x) with window of t seconds
1119 // put result in input array *pm if specified
1120 // subtract median from x if clean=true otherwise replace x with median
1121 // move running window by skip samples
1122 
1123  DataType_t* p=NULL;
1124  DataType_t* q=NULL;
1125  DataType_t* xx;
1126  double sum = 0.;
1127 
1128  size_t i,last;
1129  size_t step = Slice.stride();
1130  size_t N = Slice.size(); // number of samples in wavearray
1131  size_t n = size_t(t*rate()/step); // # of samples in the window
1132 
1133  if(n<4) {
1134  cout<<"wavearray<DataType_t>::mean() short time window"<<endl;
1135  return;
1136  }
1137 
1138  if(n&1) n--; // # of samples in the window - 1
1139 
1140  size_t nM = n/2; // index of median sample
1141  size_t nL = N-nM-1;
1142 
1143  if(pm){
1144  pm->resize(N/skip);
1145  pm->start(start());
1146  pm->rate(rate());
1147  }
1148 
1149  xx = (DataType_t *)malloc((n+1)*sizeof(DataType_t));
1150 
1151  p = data+Slice.start();
1152  q = data+Slice.start();
1153  for(i=0; i<=n; i++) {
1154  xx[i] = *p;
1155  sum += xx[i];
1156  p += step;
1157  }
1158  last = 0;
1159 
1160  for(i=0; i<N; i++){
1161 
1162  if(pm) {
1163  pm->data[i/skip] = DataType_t(sum/(n+1.));
1164  if(clean) q[i*step] -= DataType_t(sum/(n+1.));
1165  }
1166  else {
1167  if(clean) q[i*step] -= DataType_t(sum/(n+1.));
1168  else q[i*step] = DataType_t(sum/(n+1.));
1169  }
1170 
1171  if(i>=nM && i<nL) { // copy next sample into last
1172  sum -= xx[last];
1173  sum += *p;
1174  xx[last++] = *p;
1175  p += step;
1176  }
1177 
1178  if(last>n) last = 0;
1179 
1180  }
1181 
1182  free(xx);
1183  return;
1184 }
1185 
1186 
1187 template<class DataType_t>
1189 {
1190  register size_t i;
1191  register double x = 0.;
1192  register double y = 0.;
1193  size_t N = (size()>>2)<<2;
1194  register DataType_t *p = data + size() - N;
1195 
1196  if(!size()) return 0.;
1197 
1198  for(i=0; i<size()-N; i++) { x += data[i]; y += data[i]*data[i]; }
1199  for(i=0; i<N; i+=4){
1200  x += p[i] + p[i+1] + p[i+2] + p[i+3];
1201  y += p[i]*p[i] + p[i+1]*p[i+1] + p[i+2]*p[i+2] + p[i+3]*p[i+3];
1202  }
1203  x /= size();
1204  return sqrt(y/size()-x*x);
1205 }
1206 
1207 // running 50% percentile rms
1208 template<class DataType_t>
1211  bool clean,
1212  size_t skip)
1213 {
1214 
1215  DataType_t* p=NULL;
1216  DataType_t* q=NULL;
1217  DataType_t** pp;
1218  DataType_t* xx;
1219  DataType_t rm = 1;
1220 
1221  size_t i,last;
1222  size_t step = Slice.stride();
1223  size_t N = Slice.size(); // number of samples in wavearray
1224  size_t n = size_t(t*rate()/step);
1225 
1226  if(n<4) {
1227  cout<<"wavearray<DataType_t>::median() short time window"<<endl;
1228  return;
1229  }
1230 
1231  if(n&1) n--; // # of samples - 1
1232 
1233  size_t nM = n/2; // index of median sample
1234  size_t nL = N-nM-1;
1235 
1236  if(pm){
1237  pm->resize(N/skip);
1238  pm->start(start());
1239  pm->rate(rate());
1240  }
1241 
1242  pp = (DataType_t **)malloc((n+1)*sizeof(DataType_t*));
1243  xx = (DataType_t *)malloc((n+1)*sizeof(DataType_t));
1244 
1245  p = data+Slice.start();
1246  q = data+Slice.start();
1247  for(i=0; i<=n; i++) {
1248  xx[i] = *p>0 ? *p : -(*p);
1249  pp[i] = xx+i;
1250  p += step;
1251  }
1252  last = 0;
1253 
1254  for(i=0; i<N; i++){
1255 
1256  if(i==(i/skip)*skip) {
1257  waveSplit(pp,0,n,nM); // median split
1258  rm=*pp[nM];
1259  }
1260 
1261  if(pm) {
1262  pm->data[i/skip] = DataType_t(rm/0.6745);
1263  if(clean) q[i*step] *= DataType_t(0.6745/rm);
1264  }
1265  else {
1266  if(clean) q[i*step] *= DataType_t(0.6745/rm);
1267  else q[i*step] = DataType_t(rm/0.6745);
1268  }
1269 
1270  if(i>=nM && i<nL) { // copy next sample into last
1271  xx[last++] = *p>0 ? *p : -(*p);
1272  p += step;
1273  }
1274 
1275  if(last>n) last = 0;
1276 
1277  }
1278 
1279 
1280  free(pp);
1281  free(xx);
1282  return;
1283 }
1284 
1285 
1286 template<class DataType_t>
1288 {
1289  register size_t i;
1290  register double a = 0.;
1291  register double x = 0.;
1292  register double y = 0.;
1293  register DataType_t *p = data + s.start();
1294  size_t n = s.size();
1295  size_t m = (s.stride()<=0) ? 1 : s.stride();
1296 
1297  if(size()<limit(s)) n = (limit(s) - s.start() - 1)/m;
1298  if(!n) return 0.;
1299  size_t N = (n>>2)<<2;
1300 
1301  for(i=0; i<n-N; i++) {
1302  a = *p; x += a; y += a*a; p += m;
1303  }
1304 
1305  for(i=0; i<N; i+=4) {
1306  a = *p; x += a; y += a*a; p += m;
1307  a = *p; x += a; y += a*a; p += m;
1308  a = *p; x += a; y += a*a; p += m;
1309  a = *p; x += a; y += a*a; p += m;
1310  }
1311  x /= N;
1312  return sqrt(y/N-x*x);
1313 }
1314 
1315 template<class DataType_t>
1316 DataType_t wavearray<DataType_t>::max() const
1317 {
1318  if(!size()) return 0;
1319  register unsigned int i;
1320  register DataType_t x = data[0];
1321  for(i=1; i<size(); i++) { if(x<data[i]) x=data[i]; }
1322  return x;
1323 }
1324 
1325 template<class DataType_t>
1327 {
1328  if(!size() ||
1329  Slice.size()!=x.Slice.size() ||
1330  Slice.start()!=x.Slice.start() ||
1331  Slice.stride()!=x.Slice.stride()) {
1332  cout<<"wavearray::max(): illegal imput array\n";
1333  return;
1334  }
1335 
1336  size_t n;
1337  size_t K = Slice.stride();
1338  size_t I = Slice.start();
1339  size_t N = Slice.size();
1340 
1341  DataType_t* pin = x.data+I;
1342  DataType_t* pou = this->data+I;
1343  for(n=0; n<N; n+=K) {
1344  if(pou[n]<pin[n]) pou[n]=pin[n];
1345  }
1346  return;
1347 }
1348 
1349 template<class DataType_t>
1350 DataType_t wavearray<DataType_t>::min() const
1351 {
1352  if(!size()) return 0;
1353  register size_t i;
1354  register DataType_t x = data[0];
1355  for(i=1; i<size(); i++) { if(x>data[i]) x=data[i]; }
1356  return x;
1357 }
1358 
1359 
1360 template<class DataType_t>
1361 int wavearray<DataType_t>::getSampleRank(size_t n, size_t l, size_t r) const
1362 {
1363  DataType_t v;
1364  register int i = l-1; // save left boundary
1365  register int j = r; // save right boundary
1366 
1367  v = data[n]; // pivot
1368  data[n]=data[r]; data[r]=v; // store pivot
1369 
1370  while(i<j)
1371  {
1372  while(data[++i]<v && i<j);
1373  while(data[--j]>v && i<j);
1374  }
1375  data[r]=data[n]; data[n]=v; // put pivot back
1376 
1377  return i-int(l)+1; // rank ranges from 1 to r-l+1
1378 }
1379 
1380 template<class DataType_t>
1381 int wavearray<DataType_t>::getSampleRankE(size_t n, size_t l, size_t r) const
1382 {
1383  DataType_t v,vv;
1384  register int i = l-1; // save left boundary
1385  register int j = r; // save right boundary
1386 
1387  v = data[n]; // pivot
1388  data[n]=data[r]; data[r]=v; // store pivot
1389  vv = v>0 ? v : -v; // sort absolute value
1390 
1391  while(i<j)
1392  {
1393  while((data[++i]>0 ? data[i] : -data[i])<vv && i<j);
1394  while((data[--j]>0 ? data[j] : -data[j])>vv && i<j);
1395  }
1396  data[r]=data[n]; data[n]=v; // put pivot back
1397 
1398  return i-int(l)+1; // rank ranges from 1 to r-l+1
1399 }
1400 
1401 
1402 template<class DataType_t>
1403 void wavearray<DataType_t>::waveSort(DataType_t** pin, size_t l, size_t r) const
1404 {
1405 // sort data array using quick sorting algorithm between left (l) and right (r) index
1406 // DataType_t** pin is pointer to array of data pointers
1407 // sorted from min to max: *pin[l]/*pin[r] points to min/max
1408 
1409  size_t k;
1410 
1411  register DataType_t v;
1412  register DataType_t* p;
1413  register DataType_t** pp = pin;
1414 
1415  if(pp==NULL || !this->size()) return;
1416  if(!r) r = this->size()-1;
1417  if(l>=r) return;
1418 
1419  register size_t i = (r+l)/2; // median
1420  register size_t j = r-1; // pivot storage index
1421 
1422 // sort l,i,r
1423  if(*pp[l] > *pp[i]) {p=pp[l]; pp[l]=pp[i]; pp[i]=p;}
1424  if(*pp[l] > *pp[r]) {p=pp[l]; pp[l]=pp[r]; pp[r]=p;}
1425  if(*pp[i] > *pp[r]) {p=pp[i]; pp[i]=pp[r]; pp[r]=p;}
1426  if(r-l < 3) return; // all sorted
1427 
1428  v = *pp[i]; // pivot
1429  p=pp[i]; pp[i]=pp[j]; pp[j]=p; // store pivot
1430  i = l;
1431 
1432  for(;;)
1433  {
1434  while(*pp[++i] < v);
1435  while(*pp[--j] > v);
1436  if(j<i) break;
1437  p=pp[i]; pp[i]=pp[j]; pp[j]=p; // swap i,j
1438  }
1439 
1440  p=pp[i]; pp[i++]=pp[r-1]; pp[r-1]=p; // return pivot back
1441 
1442  if(j-l > 2) waveSort(pp,l,j);
1443  else if(j>l){ // sort l,k,j
1444  k = l+1;
1445  if(*pp[l] > *pp[k]) {p=pp[l]; pp[l]=pp[k]; pp[k]=p;}
1446  if(*pp[l] > *pp[j]) {p=pp[l]; pp[l]=pp[j]; pp[j]=p;}
1447  if(*pp[k] > *pp[j]) {p=pp[k]; pp[k]=pp[j]; pp[j]=p;}
1448  }
1449 
1450  if(r-i > 2) waveSort(pp,i,r);
1451  else if(r>i){ // sort i,k,r
1452  k = i+1;
1453  if(*pp[i] > *pp[k]) {p=pp[i]; pp[i]=pp[k]; pp[k]=p;}
1454  if(*pp[i] > *pp[r]) {p=pp[i]; pp[i]=pp[r]; pp[r]=p;}
1455  if(*pp[k] > *pp[r]) {p=pp[k]; pp[k]=pp[r]; pp[r]=p;}
1456  }
1457 
1458  return;
1459 }
1460 
1461 template<class DataType_t>
1463 {
1464 // sort wavearray using quick sorting algorithm between left (l) and right (r) index
1465 // sorted from min to max: this->data[l]/this->data[r] is min/max
1466  size_t N = this->size();
1467  DataType_t *pd = (DataType_t *)malloc(N*sizeof(DataType_t));
1468  DataType_t **pp = (DataType_t **)malloc(N*sizeof(DataType_t*));
1469  for(size_t i=0; i<N; i++) {
1470  pd[i] = this->data[i];
1471  pp[i] = pd + i;
1472  }
1473  waveSort(pp,l,r);
1474  for(size_t i=0; i<N; i++) this->data[i] = *pp[i];
1475 
1476  free(pd);
1477  free(pp);
1478 }
1479 
1480 template<class DataType_t>
1481 void wavearray<DataType_t>::waveSplit(DataType_t** pp, size_t l, size_t r, size_t m) const
1482 {
1483 // split input array of pointers pp[i] between l <= i <= r so that:
1484 // *pp[i] < *pp[m] if i < m
1485 // *pp[i] > *pp[m] if i > m
1486  DataType_t v;
1487  DataType_t* p;
1488  size_t i = (r+l)/2;
1489  size_t j = r-1;
1490 
1491 // sort l,i,r
1492  if(*pp[l] > *pp[i]) {p=pp[l]; pp[l]=pp[i]; pp[i]=p;}
1493  if(*pp[l] > *pp[r]) {p=pp[l]; pp[l]=pp[r]; pp[r]=p;}
1494  if(*pp[i] > *pp[r]) {p=pp[i]; pp[i]=pp[r]; pp[r]=p;}
1495  if(r-l < 3) return; // all sorted
1496 
1497  v = *pp[i]; // pivot
1498  p=pp[i]; pp[i]=pp[j]; pp[j]=p; // store pivot
1499  i = l;
1500 
1501  for(;;)
1502  {
1503  while(*pp[++i] < v);
1504  while(*pp[--j] > v);
1505  if(j<i) break;
1506  p=pp[i]; pp[i]=pp[j]; pp[j]=p; // swap i,j
1507  }
1508  p=pp[i]; pp[i]=pp[r-1]; pp[r-1]=p; // put pivot
1509 
1510  if(i > m) waveSplit(pp,l,i,m);
1511  else if(i < m) waveSplit(pp,i,r,m);
1512 
1513  return;
1514 }
1515 
1516 template<class DataType_t>
1517 DataType_t wavearray<DataType_t>::waveSplit(size_t l, size_t r, size_t m)
1518 {
1519 // split wavearray between l and r & at index m so that:
1520 // *this->data[i] < *this->datap[m] if i < m
1521 // *this->data[i] > *this->datap[m] if i > m
1522  size_t N = this->size();
1523  DataType_t *pd = (DataType_t *)malloc(N*sizeof(DataType_t));
1524  DataType_t **pp = (DataType_t **)malloc(N*sizeof(DataType_t*));
1525  for(size_t i=0; i<N; i++) {
1526  pd[i] = this->data[i];
1527  pp[i] = pd + i;
1528  }
1529  waveSplit(pp,l,r,m);
1530  for(size_t i=0; i<N; i++) this->data[i] = *pp[i];
1531 
1532  free(pd);
1533  free(pp);
1534  return this->data[m];
1535 }
1536 
1537 template<class DataType_t>
1539 {
1540 // split wavearray at index m so that:
1541 // *this->data[i] < *this->datap[m] if i < m
1542 // *this->data[i] > *this->datap[m] if i > m
1543  size_t N = this->size();
1544  DataType_t *pd = (DataType_t *)malloc(N*sizeof(DataType_t));
1545  DataType_t **pp = (DataType_t **)malloc(N*sizeof(DataType_t*));
1546  for(size_t i=0; i<N; i++) {
1547  pd[i] = this->data[i];
1548  pp[i] = pd + i;
1549  }
1550  waveSplit(pp,0,N-1,m);
1551  for(size_t i=0; i<N; i++) this->data[i] = *pp[i];
1552 
1553  free(pd);
1554  free(pp);
1555 }
1556 
1557 template<class DataType_t>
1558 double wavearray<DataType_t>::median(size_t l, size_t r) const
1559 {
1560  if(!r) r = size()-1;
1561  if(r<=l) return 0.;
1562 
1563  size_t i;
1564  size_t N = r-l+1;
1565  size_t m = N/2+(N&1); // median
1566 
1567  double x = 0.;
1568 
1569  DataType_t **pp = (DataType_t **)malloc(N*sizeof(DataType_t*));
1570  for(i=l; i<=r; i++) pp[i-l] = data + i;
1571 
1572  waveSplit(pp,0,N-1,m);
1573  x = *pp[m];
1574 
1575  free(pp);
1576  return x;
1577 }
1578 
1579 
1580 template<class DataType_t>
1583  bool clean,
1584  size_t skip)
1585 {
1586 // calculate running mean of data (x) with window of t seconds
1587 // put result in input array pm if specified
1588 // subtract median from x if clean=true otherwise replace x with median
1589 // move running window by skip samples
1590 
1591  DataType_t* p=NULL;
1592  DataType_t* q=NULL;
1593  DataType_t** pp;
1594  DataType_t* xx;
1595  DataType_t am=0;
1596 
1597  size_t i,last;
1598  size_t step = Slice.stride();
1599  size_t N = Slice.size(); // number of samples in wavearray
1600  size_t n = size_t(t*rate()/step); // # of samples in running window
1601 
1602  if(n<4) {
1603  cout<<"wavearray<DataType_t>::median() short time window"<<endl;
1604  return;
1605  }
1606 
1607  if(n&1) n--; // # of samples - 1
1608 
1609  size_t nM = n/2; // index of median sample
1610  size_t nL = N-nM-1;
1611 
1612  if(pm){
1613  pm->resize(N/skip);
1614  pm->start(start());
1615  pm->rate(rate()/skip);
1616  }
1617 
1618  pp = (DataType_t **)malloc((n+1)*sizeof(DataType_t*));
1619  xx = (DataType_t *)malloc((n+1)*sizeof(DataType_t));
1620 
1621  p = data+Slice.start();
1622  q = data+Slice.start();
1623  for(i=0; i<=n; i++) {
1624  xx[i] = *p;
1625  pp[i] = xx+i;
1626  p += step;
1627  }
1628  last = 0;
1629 
1630  for(i=0; i<N; i++){
1631 
1632  if(i==(i/skip)*skip) {
1633  waveSplit(pp,0,n,nM); // median split
1634  am=*pp[nM];
1635  }
1636 
1637  if(pm) {
1638  pm->data[i/skip] = am;
1639  if(clean) q[i*step] -= am;
1640  }
1641  else {
1642  if(clean) q[i*step] -= am;
1643  else q[i*step] = am;
1644  }
1645 
1646  if(i>=nM && i<nL) { // copy next sample into last
1647  xx[last++] = *p;
1648  p += step;
1649  }
1650 
1651  if(last>n) last = 0;
1652 
1653  }
1654 
1655  free(pp);
1656  free(xx);
1657  return;
1658 }
1659 
1660 
1661 template<class DataType_t>
1663 {
1664 
1665  DataType_t* p=NULL;
1666  DataType_t* q=NULL;
1667 
1668  size_t i;
1669  double r;
1670  size_t last,next;
1671  size_t N = Slice.size(); // number of samples in wavearray
1672  size_t step = Slice.stride();
1673  size_t n = size_t(t*rate()/step);
1674 
1675  if(n<4) {
1676  cout<<"wavearray<DataType_t>::median() short time window"<<endl;
1677  return;
1678  }
1679 
1680  if(n&1) n--; // # of samples in running window
1681 
1682  size_t nM = n/2; // index of median sample
1683  size_t nL = N-nM-1;
1684 
1685  DataType_t** pp = (DataType_t **)malloc((n+1)*sizeof(DataType_t*));
1687 
1688  p = data+Slice.start();
1689  q = data+Slice.start();
1690  for(i=0; i<=n; i++) {
1691  xx.data[i] = *p;
1692  pp[i] = xx.data+i;
1693  p += step;
1694  }
1695  last = 0;
1696  next = 0;
1697 
1698  for(i=0; i<N; i++){
1699 
1700  r = (xx.getSampleRank(next,0,n)-double(nM)-1.)/(nM+1.);
1701  q[i*step] = DataType_t(r>0. ? -log(1.-r) : log(1.+r));
1702 
1703  if(i>=nM && i<nL) { // copy next sample into last
1704  xx.data[last++] = *p; p+=step;
1705  }
1706 
1707  next++;
1708  if(next>n) next = 0;
1709  if(last>n) last = 0;
1710 
1711  }
1712 
1713  free(pp);
1714  return;
1715 }
1716 
1717 
1718 template<class DataType_t>
1719 DataType_t wavearray<DataType_t>::rank(double f) const
1720 {
1721  int i;
1722  int n = size();
1723  DataType_t out = 0;
1724  DataType_t **pp = NULL;
1725  if(f<0.) f = 0.;
1726  if(f>1.) f = 1.;
1727 
1728  if(n)
1729  pp = (DataType_t **)malloc(n*sizeof(DataType_t*));
1730  else
1731  return out;
1732 
1733  for(i=0; i<n; i++) pp[i] = data + i;
1734  qsort(pp, n, sizeof(DataType_t *),
1736 
1737  i =int((1.-f)*n);
1738  if(i==0) out = *(pp[0]);
1739  else if(i>=n-1) out = *(pp[n-1]);
1740  else out = (*(pp[i])+*(pp[i+1]))/2;
1741 
1742  for(i=0; i<n; i++) *(pp[i]) = n-i;
1743 
1744  free(pp);
1745  return out;
1746 }
1747 
1748 
1749 // symmetric prediction error signal produced with
1750 // adaptive split lattice algoritm
1751 template<class DataType_t>
1752 void wavearray<DataType_t>::spesla(double T, double w, double oFFset)
1753 {
1754  int k,j;
1755  double p,q,xp,xm;
1756 
1757  int K = int(rate()*T+0.5); // filter duration
1758  int M = int(rate()*w+0.5); // integration window
1759  if(M&1) M++; // make M even
1760  int m = M/2;
1761 
1762  int offset = int(oFFset*rate()); // data offset
1763  int N = size(); // data size
1764  int nf = offset+m;
1765  int nl = N-nf;
1766 
1769 
1770  x0 = *this; // previous X iteration
1771  x1 = *this;
1772  x1.add(x0,x0.size()-1,0,1); // current X iteration
1773  x1 *= DataType_t(0.5);
1774 
1775 // cout<<nf<<" "<<nl<<" "<<offset<<" "<<K<<" "<<M<<" "<<x0.rms()<<" "<<x1.rms();
1776 
1777  for(k=1; k<K; k++) {
1778 
1779  p = q = 0.;
1780  for(j=offset; j<offset+M; j++) {
1781  p += x1.data[j]*x1.data[j];
1782  q += x0.data[j]*x0.data[j];
1783  }
1784 
1785  for(j=1; j<N; j++) {
1786  if(j>=nf && j<nl) { // update p & q
1787  xp = x1.data[j+m];
1788  xm = x1.data[j-m];
1789  p += (xp-xm)*(xp+xm);
1790  xp = x0.data[j+m];
1791  xm = x0.data[j-m];
1792  q += (xp-xm)*(xp+xm);
1793  }
1794  xp = k==1 ? 2*p/q : p/q;
1795  data[j] = x1.data[j]+x1.data[j-1]-DataType_t(x0.data[j-1]*xp);
1796  }
1797 
1798  x0 = x1;
1799  x1 = *this;
1800 
1801  }
1802  return;
1803 }
1804 
1805 
1806 // apply filter
1807 template<class DataType_t>
1809 {
1810  int i,j;
1811  int N = size();
1812  int m = w.size();
1814  x = *this;
1815 
1816  for(i=0; i<N; i++) {
1817  for(j=1; j<m && (i-j)>=0; j++) {
1818  data[i] += DataType_t(w.data[j]*x.data[i-j]);
1819  }
1820  }
1821  return;
1822 }
1823 
1824 
1825 
1826 // calculate and apply lpr filter
1827 template<class DataType_t>
1828 void wavearray<DataType_t>::lprFilter(double T, int mode, double stride, double oFFset)
1829 {
1830  int i,j,l,n,nf,nl;
1831  double duration = stop()-start()-2*oFFset;
1832  double a=0.;
1833  double b=1.;
1834 
1835  int N = size();
1836  int M = int(rate()*T+0.5); // filter duration
1837  int k = stride>0 ? int(duration/stride) : 1; // number of intervals
1838 
1839  if(!k) k++;
1840  int m = int((N-2.*oFFset*rate())/k);
1841  if(m&1) m--; // make m even
1842 
1843  int offset = int(oFFset*rate()+0.5); // data offset
1844  if(offset&1) offset++;
1845 
1847  wavearray<DataType_t> x = *this;
1849  w.rate(rate());
1850 
1851  for(l=0; l<k; l++) {
1852 
1853  n = l*m+(N-k*m)/2;
1854  w.cpf(x,m,n);
1855  f = w.getLPRFilter(M,0);
1856 
1857  nf = l==0 ? 0 : n;
1858  nl = l==k-1 ? N : n+m;
1859  n = offset+M;
1860 
1861  if(mode == 0 || mode == 1) { // forward LP
1862  for(i=nf; i<nl; i++) {
1863  if(i < n) continue;
1864  b = (!mode && i<N-n) ? 0.5 : 1.; // symmetric LP
1865  for(j=1; j<M; j++) {
1866  a = double(f.data[j]*x.data[i-j]);
1867  data[i] += DataType_t(a*b);
1868  }
1869  }
1870  }
1871 
1872  if(mode == 0 || mode == -1) { // backward LP
1873  for(i=nf; i<nl; i++) {
1874  if(i >= N-n) continue;
1875  b = (!mode && i>=n) ? 0.5 : 1.; // symmetric LP
1876  for(j=1; j<M; j++) {
1877  a = double(f.data[j]*x.data[i+j]);
1878  data[i] += DataType_t(a*b);
1879  }
1880  }
1881  }
1882 
1883  }
1884 
1885  return;
1886 }
1887 
1888 //**************************************************************
1889 // calculate autocorrelation function and
1890 // solve symmetric Yule-Walker problem
1891 //**************************************************************
1892 template<class DataType_t>
1894 {
1895  size_t i,m;
1896  double f = 0.03; // exclude tail
1897  size_t N = size()-2*offset;
1898  size_t nL = size_t(f*N+0.5); // left percentile
1899  size_t nR = N-nL-1; // right percentile
1900  size_t nn = N/2; // median
1901  size_t n = size()-offset;
1902 
1903  if(size()<=offset) {
1904  cout<<"wavearray<DataType_t>::getLPRFilter() invalid input parameters\n";
1905  wavearray<double> a(1);
1906  return a;
1907  }
1908 
1909  wavearray<double> r(M);
1910  wavearray<double> a(M);
1911  wavearray<double> b(M);
1912 
1913  wavearray<DataType_t> x = *this;
1914 
1915  DataType_t ** pp = (DataType_t **)malloc(N*sizeof(DataType_t*));
1916 
1917  for(i=offset; i<n; i++) pp[i-offset] = x.data + i;
1918 
1919  waveSplit(pp,0,N-1,nn);
1920  waveSplit(pp,0,nn,nL);
1921  waveSplit(pp,nn,N-1,nR);
1922 
1923  x -= *pp[nn]; // subtract median
1924  for(i=0; i<nL; i++) *pp[i] = 0; // exclude tails
1925  for(i=nR; i<N; i++) *pp[i] = 0; // exclude tails
1926 
1927 // autocorrelation
1928 
1929  offset += M;
1930  n = size()-offset;
1931 
1932  for (m=0; m<M; m++) {
1933  r.data[m] = 0.;
1934  for (i=offset; i<n; i++) {
1935  r.data[m] += x.data[i]*(x.data[i-m] + x.data[i+m])/2.;
1936  }
1937  }
1938 
1939 
1940 // Levinson algorithm: P.Delsarte & Y.Genin, IEEE, ASSP-35, #5 May 1987
1941 
1942  double p,q;
1943  double s = r.data[0];
1944 
1945  a = 0.; a.data[0] = 1.;
1946 
1947  for (m=1; m<M; m++) {
1948 
1949  q = 0.;
1950  for (i=0; i<m; i++) q -= a.data[i] * r.data[m-i];
1951 
1952  p = q/s; // reflection coefficient
1953  s -= q*p;
1954 
1955  for (i=1; i<=m; i++) b.data[i] = a.data[i] + p*a.data[m-i];
1956  for (i=1; i<=m; i++) a.data[i] = b.data[i];
1957 
1958  }
1959 
1960 /*
1961  double tmp, num, den;
1962  M--;
1963 
1964  a.data[1] = - r.data[1] / r.data[0];
1965 
1966  for (m=1; m<M; m++) {
1967 
1968  num = r.data[m + 1];
1969  den = r.data[0];
1970 
1971  for (i=1; i<=m; i++) {
1972  num += a.data[i] * r.data[m+1-i];
1973  den += a.data[i] * r.data[i];
1974  }
1975 
1976  a.data[m+1] = - num / den;
1977 
1978  for (i=1; i <= ((m+1)>>1); i++) {
1979  tmp = a.data[i] + a.data[m+1] * a.data[m+1-i];
1980  a.data[m+1-i] = a.data[m+1-i] + a.data[m+1] * a.data[i];
1981  a.data[i] = tmp;
1982  }
1983 
1984  }
1985  a.data[0] = 1;
1986 */
1987 
1988  free(pp);
1989  return a;
1990 }
1991 
1992 
1993 //**************************************************************
1994 // normalize data by noise variance
1995 // offset | stride |
1996 // |*****|*X********X********X********X********X********X*|*****|
1997 // ^ ^ ^
1998 // measurement |<- interval ->|
1999 //**************************************************************
2000 template<class DataType_t>
2002 wavearray<DataType_t>::white(double t, int mode, double oFFset, double stride) const
2003 {
2004  int i,j;
2005  int N = size();
2006 
2007  double segT = size()/rate(); // segment duration
2008  if(t<=0.) t = segT-2.*oFFset;
2009  int offset = int(oFFset*rate()+0.5); // offset in samples
2010  if(offset&1) offset--; // make offset even
2011 
2012  if(stride > t || stride<=0.) stride = t;
2013  int K = int((segT-2.*oFFset)/stride); // number of noise measurement minus 1
2014  if(!K) K++; // special case
2015 
2016  int n = N-2*offset; // total number of samples used for noise estimate
2017  int k = n/K; // shift in samples
2018  if(k&1) k--; // make k even
2019 
2020  int m = int(t*rate()+0.5); // number of samples in one interval
2021  int mm = m/2; // median m
2022  int mL = int(0.15865*m+0.5); // -sigma index (CL=0.31732)
2023  int mR = m-mL-1; // +sigma index
2024  int jL = (N-k*K)/2; // array index for first measurement
2025  int jR = N-offset-m; // array index for last measurement
2026  int jj = jL-mm; // start of the first interval
2027 
2028  wavearray<double> meDIan(K+1);
2029  wavearray<double> norm50(K+1);
2030 
2031  meDIan.rate(rate()/k);
2032  meDIan.start(start()+jL/rate());
2033  meDIan.stop(start()+(N-offset)/rate());
2034 
2035  norm50.rate(rate()/k);
2036  norm50.start(start()+jL/rate());
2037  norm50.stop(start()+(N-offset)/rate());
2038 
2039  mode = abs(mode);
2040 
2041  if(m<3 || mL<2 || mR>m-2) {
2042  cout<<"wavearray::white(): too short input array."<<endl;
2043  return mode!=0 ? norm50 : meDIan;
2044  }
2045 
2046  register DataType_t *p = data;
2048  double x;
2049  double r;
2050 
2051  DataType_t ** pp = (DataType_t **)malloc(m*sizeof(DataType_t*));
2052 
2053  for(j=0; j<=K; j++) {
2054 
2055  if(jj < offset) p = data + offset; // left boundary
2056  else if(jj >= jR) p = data + jR; // right boundary
2057  else p = data + jj;
2058  jj += k; // update jj
2059 
2060  if(p+m>data+N) {cout<<"wavearray::white(): error1\n"; exit(1);}
2061 
2062  for(i=0; i<m; i++) pp[i] = p + i;
2063  waveSplit(pp,0,m-1,mm);
2064  waveSplit(pp,0,mm,mL);
2065  waveSplit(pp,mm,m-1,mR);
2066  meDIan[j] = mode ? *pp[mm] : sqrt((*pp[mm])*0.7191);
2067  norm50[j] = (*pp[mR] - *pp[mL])/2.;
2068 
2069  }
2070 
2071  p = data;
2072 
2073  if(mode) {
2074 
2075  mm = jL;
2076  for(i=0; i<mm; i++){
2077  x = double(*p)-meDIan.data[0];
2078  r = norm50.data[0];
2079  *(p++) = mode==1 ? DataType_t(x/r) : DataType_t(x/r/r);
2080  }
2081 
2082  for(j=0; j<K; j++) {
2083  for(i=0; i<k; i++) {
2084  x = double(*p)-(meDIan.data[j+1]*i + meDIan.data[j]*(k-i))/k;
2085  r = (norm50.data[j+1]*i + norm50.data[j]*(k-i))/k;
2086  *(p++) = mode==1 ? DataType_t(x/r) : DataType_t(x/r/r);
2087  }
2088  }
2089 
2090  mm = (data+N)-p;
2091  for(i=0; i<mm; i++){
2092  x = double(*p)-meDIan.data[K];
2093  r = norm50.data[K];
2094  *(p++) = mode==1 ? DataType_t(x/r) : DataType_t(x/r/r);
2095  }
2096 
2097  }
2098 
2099  free(pp);
2100  return mode!=0 ? norm50 : meDIan;
2101 }
2102 
2103 //: returns mean and root mean square of the signal.
2104 template<class DataType_t>
2105 double wavearray<DataType_t>::getStatistics(double &m, double &r) const
2106 {
2107  register size_t i;
2108  register double a;
2109  register double b;
2110  DataType_t *p = const_cast<DataType_t *>(data);
2111  double y = 0.;
2112  size_t N = size() - 1 + size_t(size()&1);
2113 
2114  if(!size()) return 0.;
2115 
2116  m = p[0];
2117  r = p[0]*p[0];
2118  if(N < size()){
2119  m += p[N];
2120  r += p[N]*p[N];
2121  y += p[N]*p[N-1];
2122  }
2123 
2124  for(i=1; i<N; i+=2) {
2125  a = p[i]; b = p[i+1];
2126  m += a + b;
2127  r += a*a + b*b;
2128  y += a*(p[i-1]+b);
2129  }
2130 
2131  N = size();
2132  m = m/double(N);
2133  r = r/double(N) - m*m;
2134 
2135  y = 4.*(y/N - m*m + m*(p[0]+p[i]-m)/N);
2136  y /= 4.*r - 2.*((p[0]-m)*(p[0]-m)+(p[i]-m)*(p[i]-m))/N;
2137  r = sqrt(r);
2138 
2139  a = (fabs(y) < 1) ? sqrt(0.5*(1.-fabs(y))) : 0.;
2140 
2141  return y>0 ? -a : a;
2142 }
2143 
2144 template <class DataType_t>
2145 void wavearray<DataType_t>::Streamer(TBuffer &R__b)
2146 {
2147  // Stream an object of class wavearray<DataType_t>.
2148  UInt_t R__s, R__c;
2149  if (R__b.IsReading()) {
2150  Version_t R__v = R__b.ReadVersion(&R__s, &R__c); if (R__v) { }
2151  TNamed::Streamer(R__b);
2152  R__b >> Size;
2153  R__b >> Rate;
2154  R__b >> Start;
2155  R__b >> Stop;
2156  if(R__v>1) R__b >> Edge;
2157  R__b.StreamObject(&(Slice),typeid(slice));
2158  this->resize(Size);
2159  R__b.ReadFastArray((char*)data,Size*sizeof(DataType_t));
2160  R__b.CheckByteCount(R__s, R__c, wavearray<DataType_t>::IsA());
2161  } else {
2162  R__c = R__b.WriteVersion(wavearray<DataType_t>::IsA(), kTRUE);
2163  TNamed::Streamer(R__b);
2164  R__b << Size;
2165  R__b << Rate;
2166  R__b << Start;
2167  R__b << Stop;
2168  R__b << Edge;
2169  R__b.StreamObject(&(Slice),typeid(slice));
2170  R__b.WriteFastArray((char*)data,Size*sizeof(DataType_t));
2171  R__b.SetByteCount(R__c, kTRUE);
2172  }
2173 }
2174 
2175 //: save to file operator >>
2176 template <class DataType_t>
2178 {
2179  TString name = fname;
2180  name.ReplaceAll(" ","");
2181  TObjArray* token = TString(fname).Tokenize(TString("."));
2182  TObjString* ext_tok = (TObjString*)token->At(token->GetEntries()-1);
2183  TString ext = ext_tok->GetString();
2184  if(ext=="dat") {
2185  //: Dump data array to a binary file
2186  DumpBinary(fname,0);
2187  } else
2188  if(ext=="txt") {
2189  //: Dump data array to an ASCII file
2190  Dump(fname,0);
2191  } else
2192  if(ext=="root") {
2193  //: Dump object into root file
2194  DumpObject(fname);
2195  } else {
2196  cout << "wavearray operator (>>) : file type " << ext.Data() << " not supported" << endl;
2197  }
2198  return fname;
2199 }
2200 
2201 //: save to file operator >>
2202 template <class DataType_t>
2204 {
2205  cout << endl << endl;
2206  cout.precision(14);
2207  cout << "Size\t= " << this->Size << " samples" << endl;
2208  cout << "Rate\t= " << this->Rate << " Hz" << endl;
2209  cout << "Start\t= " << this->Start << " sec" << endl;
2210  cout << "Stop\t= " << this->Stop << " sec" << endl;
2211  cout << "Length\t= "<< this->Stop-this->Start << " sec" << endl;
2212  cout << "Edge\t= " << this->Edge << " sec" << endl;
2213  cout << "Mean\t= " << this->mean() << endl;
2214  cout << "RMS\t= " << this->rms() << endl;
2215  cout << endl;
2216 
2217  return;
2218 }
2219 
2220 #ifdef _USE_DMT
2221 
2222 template<class DataType_t>
2224 operator=(const TSeries &ts)
2225 {
2226  Interval ti = ts.getTStep();
2227  double Tsample = ti.GetSecs();
2228 
2229  unsigned int n=ts.getNSample();
2230  if(n != size()) resize(n);
2231 
2232  if ( Tsample > 0. )
2233  rate(double(int(1./Tsample + 0.5)));
2234  else {
2235  cout <<" Invalid sampling interval = 0 sec.\n";
2236  }
2237 
2238  start(double(ts.getStartTime().totalS()));
2239 
2240  TSeries r(ts.getStartTime(), ti, ts.getNSample());
2241  r = ts;
2242  float *vec_ref;
2243  vec_ref= (float*)(r.refData());
2244 
2245  for ( unsigned int i=0; i<n; i++ )
2246  data[i] = (DataType_t) (vec_ref[i]);
2247  return *this;
2248 }
2249 
2250 #endif
2251 
2252 // instantiations
2253 
2254 #define CLASS_INSTANTIATION(class_) template class wavearray< class_ >;
2255 
2257 CLASS_INSTANTIATION(short)
2258 CLASS_INSTANTIATION(long)
2259 CLASS_INSTANTIATION(long long)
2260 CLASS_INSTANTIATION(unsigned int)
2261 CLASS_INSTANTIATION(float)
2262 CLASS_INSTANTIATION(double)
2263 //CLASS_INSTANTIATION(std::complex<float>)
2264 //CLASS_INSTANTIATION(std::complex<double>)
2265 
2266 #undef CLASS_INSTANTIATION
2267 
2268 
2269 #if !defined (__SUNPRO_CC)
2270 template wavearray<double>::
2271 wavearray(const double *, unsigned int, double);
2272 
2273 template wavearray<double>::
2274 wavearray(const float *, unsigned int, double );
2275 
2276 template wavearray<double>::
2277 wavearray(const short *, unsigned int, double );
2278 
2279 template wavearray<float>::
2280 wavearray(const double *, unsigned int, double );
2281 
2282 template wavearray<float>::
2283 wavearray(const float *, unsigned int, double );
2284 
2285 template wavearray<float>::
2286 wavearray(const short *, unsigned int, double );
2287 
2288 #else
2289 // FAKE CALLS FOR SUN CC since above instatinations are
2290 // not recognized
2291 static void fake_instatination_SUN_CC ()
2292 {
2293  double x;
2294  float y;
2295  short s;
2296  wavearray<double> wvdd (&x, 1, 0);
2297  wavearray<double> wvdf (&y, 1, 0);
2298  wavearray<double> wvds (&s, 1, 0);
2299  wavearray<float> wvfd (&x, 1, 0);
2300  wavearray<float> wvff (&y, 1, 0);
2301  wavearray<float> wvfs (&s, 1, 0);
2302 }
2303 
2304 #endif
2305 
2306 
2307 
2308 
2309 
2310 
2311 
2312 
2313 
wavearray< double > t(hp.size())
size_t append(const wavearray< DataType_t > &)
Definition: wavearray.cc:775
virtual size_t size() const
Definition: wavearray.hh:127
static const double C
Definition: GNGen.cc:10
double duration
tuple f
Definition: cwb_online.py:91
fprintf(stdout,"start=%f duration=%f rate=%f\n", x.start(), x.size()/x.rate(), x.rate())
int offset
Definition: TestSTFT_2.C:19
virtual void rate(double r)
Definition: wavearray.hh:123
virtual void ReadShort(const char *)
Definition: wavearray.cc:422
size_t Size
data array
Definition: wavearray.hh:305
#define np
wavearray< double > a(hp.size())
WSeries< float > v[nIFO]
Definition: cwb_net.C:62
par[0] name
int n
Definition: cwb_net.C:10
void print()
Definition: wavearray.cc:2203
TString("c")
double getStatistics(double &mean, double &rms) const
Definition: wavearray.cc:2105
cout<< endl;cout<< "ts size = "<< ts.size()<< " ts rate = "<< ts.rate()<< endl;tf.Forward(ts, wdm);int levels=tf.getLevel();cout<< "tf size = "<< tf.size()<< endl;double dF=tf.resolution();double dT=1./(2 *dF);cout<< "rate(hz) : "<< RATE<< "\t layers : "<< nLAYERS<< "\t dF(hz) : "<< dF<< "\t dT(ms) : "<< dT *1000.<< endl;int itime=TIME_PIXEL_INDEX;int ifreq=FREQ_PIXEL_INDEX;int index=(levels+1)*itime+ifreq;double time=itime *dT;double freq=(ifreq >0)?ifreq *dF:dF/4;cout<< endl;cout<< "PIXEL TIME = "<< time<< " sec "<< endl;cout<< "PIXEL FREQ = "<< freq<< " Hz "<< endl;cout<< endl;wavearray< double > x
wavearray< DataType_t > & operator=(const wavearray< DataType_t > &)
Definition: wavearray.cc:92
virtual int getSampleRankE(size_t n, size_t l, size_t r) const
Definition: wavearray.cc:1381
double Rate
Definition: wavearray.hh:306
double median(std::vector< double > &vec)
Definition: wavegraph.cc:467
void add(const wavearray< DataType_t > &, int=0, int=0, int=0)
Definition: wavearray.cc:728
virtual void ReadBinary(const char *, int=0)
Definition: wavearray.cc:392
double Stop
Definition: wavearray.hh:308
STL namespace.
virtual wavearray< DataType_t > & operator*=(wavearray< DataType_t > &)
Definition: wavearray.cc:258
virtual ~wavearray()
Definition: wavearray.cc:82
double Nevill(const double x0, int n, DataType_t *p, double *q)
Definition: watfun.hh:38
Long_t size
virtual double median(size_t=0, size_t=0) const
Definition: wavearray.cc:1558
#define M
Definition: UniqSLagsList.C:3
int m
Definition: cwb_net.C:10
virtual void start(double s)
Definition: wavearray.hh:119
int j
Definition: cwb_net.C:10
i drho i
int isize
double Edge
Definition: wavearray.hh:309
virtual DataType_t min() const
Definition: wavearray.cc:1350
#define N
virtual double mean() const
Definition: wavearray.cc:1053
tuple ff
Definition: cwb_online.py:394
double pi
Definition: TestChirp.C:18
cout<< "Selected Pixels : "<< nPix<< endl;wc.cluster(1, 1);SSeries< double > ss
size_t mode
DataType_t rank(double=0.5) const
Definition: wavearray.cc:1719
virtual double rms()
Definition: wavearray.cc:1188
wavearray< double > w
Definition: Test1.C:27
virtual DataType_t max() const
Definition: wavearray.cc:1316
ofstream out
Definition: cwb_merge.C:196
pTF[i] DumpBinary(file)
virtual wavearray< double > getLPRFilter(size_t, size_t=0)
Definition: wavearray.cc:1893
tlive_fix Write()
virtual void DumpObject(const char *)
Definition: wavearray.cc:382
void wavefft(double a[], double b[], int ntot, int n, int nspan, int isn)
Definition: wavefft.cc:23
void sub(const wavearray< DataType_t > &, int=0, int=0, int=0)
Definition: wavearray.cc:754
virtual char * operator>>(char *)
Definition: wavearray.cc:2177
double time[6]
Definition: cbc_plots.C:435
wavearray< double > xx
Definition: TestFrame1.C:11
x FFTW(-1)
size_t size() const
Definition: wslice.hh:71
i() int(T_cor *100))
double Pi
x resize(N)
bool log
Definition: WaveMDC.C:41
void Resample(const wavearray< DataType_t > &, double, int=6)
Definition: wavearray.cc:622
TIter next(twave->GetListOfBranches())
char fname[1024]
std::slice Slice
Definition: wavearray.hh:310
int k
virtual void resetFFTW()
Definition: wavearray.cc:959
void Dump(WSeries< double > &w, double t1, double t2, const char *fname)
virtual wavearray< DataType_t > & operator+=(wavearray< DataType_t > &)
Definition: wavearray.cc:173
r add(tfmap,"hchannel")
virtual int getSampleRank(size_t n, size_t l, size_t r) const
Definition: wavearray.cc:1361
TObjArray * token
int Rate
virtual wavearray< DataType_t > & operator<<(wavearray< DataType_t > &)
Definition: wavearray.cc:138
double dt
virtual void DumpBinary(const char *, int=0)
Definition: wavearray.cc:335
virtual void FFTW(int=1)
Definition: wavearray.cc:878
regression r
Definition: Regression_H1.C:44
Double_t x0
s s
Definition: cwb_net.C:137
double T
Definition: testWDM_4.C:11
virtual void lprFilter(wavearray< double > &)
Definition: wavearray.cc:1808
virtual void waveSort(DataType_t **pp, size_t l=0, size_t r=0) const
Definition: wavearray.cc:1403
wavearray< double > ts(N)
virtual void delay(double T)
Definition: wavearray.cc:578
virtual void DumpShort(const char *, int=0)
Definition: wavearray.cc:356
virtual wavearray< double > white(double, int=0, double=0., double=0.) const
Definition: wavearray.cc:2002
virtual void stop(double s)
Definition: wavearray.hh:121
virtual wavearray< DataType_t > & operator-=(wavearray< DataType_t > &)
Definition: wavearray.cc:208
double fabs(const Complex &x)
Definition: numpy.cc:37
virtual void spesla(double, double, double=0.)
Definition: wavearray.cc:1752
void resample(const wavearray< DataType_t > &, double, int=6)
Definition: wavearray.cc:485
virtual wavearray< DataType_t > & operator[](const std::slice &)
Definition: wavearray.cc:277
strcpy(RunLabel, RUN_LABEL)
virtual void Dump(const char *, int=0)
Definition: wavearray.cc:301
double df
int l
Definition: cbc_plots.C:434
virtual void exponential(double)
Definition: wavearray.cc:1662
int(* qsort_func)(const void *, const void *)
Definition: wavearray.cc:29
#define CLASS_INSTANTIATION(class_)
Definition: wavearray.cc:2254
DataType_t * data
Definition: wavearray.hh:301
double Stack(const wavearray< DataType_t > &, int)
Definition: wavearray.cc:973
fclose(ftrig)
size_t stride() const
Definition: wslice.hh:75
for(int i=0;i< 101;++i) Cos2[2][i]=0
virtual void resize(unsigned int)
Definition: wavearray.cc:445
double length
Definition: TestBandPass.C:18
virtual void waveSplit(DataType_t **pp, size_t l, size_t r, size_t m) const
Definition: wavearray.cc:1481
wavearray< double > y
Definition: Test10.C:31
char pm[8]
virtual void FFT(int=1)
Definition: wavearray.cc:814
size_t start() const
Definition: wslice.hh:67
exit(0)
void cpf(const wavearray< DataType_t > &, int=0, int=0, int=0)
Definition: wavearray.cc:699
double avr