Logo Coherent WaveBurst  
Reference Guide
Logo
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
gwavearray.cc
Go to the documentation of this file.
1 #include "gwavearray.hh"
2 
3 ClassImp(gwavearray<double>)
4 
5 //______________________________________________________________________________
6 /* Begin_Html
7 <center><h2>gwavearray class</h2></center>
8 This class gwavearray is derived from the wavearray class
9 It is used to produce time,freq and time-freq plots of a signal stored in a wavearray structure.
10 
11 <p>
12 <h3><a name="example">Example</a></h3>
13 <p>
14 The macro <a href="./tutorials/gwat/DrawGWaveArray.C.html">DrawGWaveArray.C</a> is an example which shown how to use the gWaveArray class.<br>
15 The pictures below gives the macro output plots.<br>
16 <p>
17 
18 End_Html
19 
20 Begin_Macro
21 DrawGWaveArray.C
22 End_Macro */
23 
24 using namespace std;
25 
26 
27 //______________________________________________________________________________
28 // destructor
29 template<class DataType_t>
30 gwavearray<DataType_t>::~gwavearray()
31 {
32  if(wextern) this->data = NULL;
33 
34  if(stft!=NULL) delete stft;
35  if(pts!=NULL) delete pts;
36 }
37 
38 //______________________________________________________________________________
39 template<class DataType_t>
41 wavearray<DataType_t>(w) {Init();wextern=false;}
42 
43 //______________________________________________________________________________
44 template<class DataType_t>
46 {
47  Init();
48  wextern = true;
49 
50  this->data = w->data;
51  this->Rate = w->rate();
52  this->Size = w->size();
53  this->Start = w->start();
54  this->Stop = w->start() + w->size()/w->rate();
55  this->Slice = std::slice(0,w->size(),1);
56 }
57 
58 //______________________________________________________________________________
59 template<class DataType_t>
61 
62  gRandom->SetSeed(0);
63  tMin=0;
64  tMax=0;
65 
66  stft=NULL;pts=NULL;
67 }
68 
69 //______________________________________________________________________________
70 template<class DataType_t>
72 //
73 // Draw waveform in time/frequency domain
74 //
75 // Input:
76 // type - GWAT_TIME, GWAT_FFT, GWAT_TF
77 // options - graphic options
78 // color - option for GWAT_TIME, GWAT_FFT
79 //
80 
81  return Draw(NULL,type,options,color);
82 }
83 
84 //______________________________________________________________________________
85 template<class DataType_t>
87 //
88 // Draw waveform in time/frequency domain
89 //
90 // Input:
91 // wavearray - data array
92 // type - GWAT_TIME, GWAT_FFT, GWAT_TF
93 // options - graphic options
94 // color - option for GWAT_TIME, GWAT_FFT
95 //
96 
97  if(type==GWAT_TIME) DrawTime(x,options,color);
98  if(type==GWAT_FFT) DrawFFT(x,options,color);
99  if(type==GWAT_TF) DrawTF(x,options);
100  if(type==GWAT_SG) cout << "gwavearray<DataType_t>::Draw - Error : draw type "
101  << "GWAT_SG" << " not allowed" << endl;
102 }
103 
104 //______________________________________________________________________________
105 template<class DataType_t>
107 //
108 // Draw waveform in time domain
109 //
110 // Input:
111 // options - draw options (same as TGraph)
112 // if contains FULL the full interval is showed
113 //
114 // return pointer to watplot object
115 //
116 
117  return DrawTime(NULL,options,color);
118 }
119 
120 //______________________________________________________________________________
121 template<class DataType_t>
123 //
124 // Draw waveform in time domain
125 //
126 // Input:
127 // wavearray - data array
128 // options - draw options (same as TGraph)
129 // if contains FULL the full interval is showed
130 //
131 // return pointer to watplot object
132 //
133 
134  if(this->rate()<=0) {
135  cout << "gwavearray::DrawTime : Error - rate must be >0" << endl;
136  exit(1);
137  }
138 
139  options.ToUpper();
140  if(stft!=NULL) delete stft;
141  if(options.Contains("SAME")&&(pts!=NULL)) {
142  } else {
143  if(pts!=NULL) delete pts;
144  char name[32];sprintf(name,"TIME-gID:%d",int(gRandom->Rndm(13)*1.e9));
145  pts = new watplot(const_cast<char*>(name),200,20,800,500);
146  }
147  double tStart,tStop;
148  if(options.Contains("FULL")) {
149  options.ReplaceAll("FULL","");
150  tStart=this->start();
151  tStop=this->start()+this->size()/this->rate();
152  tStop-=this->start();tStart-=this->start();
153  } else if(options.Contains("CUSTOM")) {
154  options.ReplaceAll("CUSTOM","");
155  tStart=this->start();
156  tStop=this->start()+this->size()/this->rate();
157  tStart = tMin>tStart&&tMin<tStop ? tMin : tStart;
158  tStop = tMax>tStart&&tMax<tStop ? tMax : tStop;
159  tStop-=this->start();tStart-=this->start();
160  } else {
161  GetTimeRange(tStart, tStop);
162  }
164  if(x!=NULL) waveAssign(tf,*x);
165  else waveAssign(tf,*this);
166  char title[256];
167  sprintf(title,"START : %10.3f (gps) - LENGHT : %g (sec) - RATE : %g (hz)",
168  tf.start()+tStart,tStop-tStart,tf.rate());
169  if(!options.Contains("SAME")) {
170  tSave=tf.start(); // save start time of first plot
171  tf.start(0);
172  } else {
173  tf.start(tf.start()-tSave);
174  }
175  if(!options.Contains("SAME")) options=options+TString(" ALP");
176  options.ReplaceAll(" ","");
177  pts->plot(tf, const_cast<char*>(options.Data()), color, tStart, tStop);
178  pts->graph[pts->graph.size()-1]->SetTitle(title);
179 
180  return pts;
181 }
182 
183 //______________________________________________________________________________
184 template<class DataType_t>
186 //
187 // Draw waveform spectrum
188 //
189 // Input:
190 // wavearray - data array
191 // options - draw options (same as TGraph)
192 // if contains FULL the full interval is showed
193 // if contains NOLOGX/NOLOGY X/Y axis are linear
194 //
195 // return pointer to watplot object
196 //
197 
198  return DrawFFT(NULL,options,color);
199 }
200 
201 //______________________________________________________________________________
202 template<class DataType_t>
204 //
205 // Draw waveform spectrum
206 //
207 // Input:
208 // wavearray - data array
209 // options - draw options (same as TGraph)
210 // if contains FULL the full interval is showed
211 // if contains NOLOGX/NOLOGY X/Y axis are linear
212 //
213 // return pointer to watplot object
214 //
215 
216  options.ToUpper();
217  if(stft!=NULL) delete stft;
218  if(options.Contains("SAME")&&(pts!=NULL)) {
219  } else {
220  if(pts!=NULL) delete pts;
221  char name[32];sprintf(name,"FFT-gID:%d",int(gRandom->Rndm(13)*1.e9));
222  pts = new watplot(const_cast<char*>(name),200,20,800,500);
223  }
224 
225  double fLow = 4.;
226  double fHigh = this->rate()/2.;
227  double tStart,tStop;
228  bool psd = options.Contains("PSD") ? true : false;
229  options.ReplaceAll("PSD","");
230  if(options.Contains("FULL")) {
231  options.ReplaceAll("FULL","");
232  tStart=this->start();
233  tStop=this->start()+this->size()/this->rate();
234  tStop-=this->start();tStart-=this->start();
235  } else if(options.Contains("CUSTOM")) {
236  options.ReplaceAll("CUSTOM","");
237  tStart=this->start();
238  tStop=this->start()+this->size()/this->rate();
239  tStart = tMin>tStart&&tMin<tStop ? tMin : tStart;
240  tStop = tMax>tStart&&tMax<tStop ? tMax : tStop;
241  tStop-=this->start();tStart-=this->start();
242  } else {
243  GetTimeRange(tStart, tStop);
244  }
245  bool logx = true;
246  if(options.Contains("NOLOGX")) {logx=false;options.ReplaceAll("NOLOGX","");}
247  bool logy = true;
248  if(options.Contains("NOLOGY")) {logy=false;options.ReplaceAll("NOLOGY","");}
250  if(x!=NULL) waveAssign(tf,*x);
251  else waveAssign(tf,*this);
252  // use only the selected range
253  int nb=tStart*tf.rate();
254  int ne=tStop*tf.rate();
255  for(int i=nb;i<ne;i++) tf[i-nb] = tf[i];
256  tf.resize(ne-nb);
257  char title[256];
258  sprintf(title,"START : %10.3f (gps) - LENGHT : %g (sec) - RATE : %g (hz)",
259  tf.start()+tStart,tStop-tStart,tf.rate());
260  tf.start(0);
261  if(!options.Contains("SAME")) options=options+TString(" ALP");
262  options.ReplaceAll(" ","");
263  pts->plot(tf, const_cast<char*>(options.Data()), color, 0, tStop-tStart, true, fLow, fHigh, psd, 8);
264  pts->graph[pts->graph.size()-1]->SetTitle(title);
265  if(logx) pts->canvas->SetLogx();
266  if(logy) pts->canvas->SetLogy();
267 
268  return pts;
269 }
270 
271 //______________________________________________________________________________
272 template<class DataType_t>
274 //
275 // Draw waveform spectrogram
276 //
277 // Input:
278 // wavearray - data array
279 // options - draw options (same as TH2D)
280 //
281 // return pointer to CWB::STFT object
282 //
283 
284  return DrawTF(NULL,options);
285 }
286 
287 //______________________________________________________________________________
288 template<class DataType_t>
290 //
291 // Draw waveform spectrogram
292 //
293 // Input:
294 // options - draw options (same as TH2D)
295 //
296 // return pointer to CWB::STFT object
297 //
298 
299  int nfact=4;
300  int nfft=nfact*512;
301  int noverlap=nfft-10;
302  double fparm=nfact*6;
303 
304  double tStart,tStop;
305  if(options.Contains("FULL")) {
306  options.ReplaceAll("FULL","");
307  tStart=this->start();
308  tStop=this->start()+this->size()/this->rate();
309  tStop-=this->start();tStart-=this->start();
310  } else if(options.Contains("CUSTOM")) {
311  options.ReplaceAll("CUSTOM","");
312  tStart=this->start();
313  tStop=this->start()+this->size()/this->rate();
314  tStart = tMin>tStart&&tMin<tStop ? tMin : tStart;
315  tStop = tMax>tStart&&tMax<tStop ? tMax : tStop;
316  tStop-=this->start();tStart-=this->start();
317  } else {
318  GetTimeRange(tStart, tStop);
319  }
320 
321  if(stft!=NULL) delete stft;
322  if(pts!=NULL) delete pts;
323 
325  if(x!=NULL) waveAssign(tf,*x);
326  else waveAssign(tf,*this);
327  tf.start(0);
328  char name[32];sprintf(name,"TF-gID:%d",int(gRandom->Rndm(13)*1.e9));
329  stft = new CWB::STFT(tf,nfft,noverlap,"amplitude","gauss",fparm,name);
330 
331  double fLow = 4.;
332  double fHigh = this->rate()/2.;
333  stft->Draw(tStart,tStop,fLow,fHigh,0,0,1);
334  stft->GetCanvas()->SetLogy(true);
335 
336  return stft;
337 }
338 
339 //______________________________________________________________________________
340 template<class DataType_t>
341 double gwavearray<DataType_t>::GetTimeRange(double& tMin, double& tMax, double efraction) {
342 //
343 // Get mdc time interval
344 //
345 // Input: efraction - energy fraction used to evaluate the range which contains efraction
346 // default = 0.9999999
347 //
348 // Output: tMin - start time interval
349 // tMax - stop time interval
350 //
351 
352  double a,b;
353 
354  double T = GetCentralTime();
355 
356  double E=0.;
357  int size=(int)this->size();
358  double rate=this->rate();
359  for(int j=0;j<size;j++) E += this->data[j]*this->data[j];
360 
361  int OS = int(T*rate);
362  int M = OS<size/2 ? OS : size-OS;
363 
364  int I=0;
365  double sum = this->data[OS]*this->data[OS];
366  for(int j=1; j<int(M); j++) {
367  a = this->data[OS-j];
368  b = this->data[OS+j];
369  sum += a*a+b*b;
370  if(sum/E > efraction) break;
371  I=j;
372  }
373  I++;
374  tMin = tMin>=0 ? double((OS-I)/rate) : 0.;
375  tMax = tMax<size/rate ? double((OS+I)/rate) : size/rate;
376 
377  if(sum/E<efraction) {tMin=0.;tMax=size/rate;}
378 
379  return T;
380 }
381 
382 //______________________________________________________________________________
383 template<class DataType_t>
385 //
386 // Get mdc central time
387 //
388 // Input:
389 //
390 // Retun central time
391 //
392 
393  double a;
394  double E=0.,T=0.;
395  int size=(int)this->size();
396  double rate=this->rate();
397  for(int j=0;j<size;j++) {
398  a = this->data[j];
399  T += a*a*j/rate; // central time
400  E += a*a; // energy
401  }
402  T = E>0 ? T/E : 0.5*size/rate;
403 
404  return T;
405 }
406 
407 //______________________________________________________________________________
408 template<class DataType_t>
410 //
411 // apply time shift
412 //
413 //
414 // Input:
415 // time - time shift (sec)
416 //
417 
418  if(tShift==0) return;
419 
420  // search begin,end of non zero data
421  int ibeg=0; int iend=0;
422  for(int i=0;i<(int)this->size();i++) {
423  if(this->data[i]!=0 && ibeg==0) ibeg=i;
424  if(this->data[i]!=0) iend=i;
425  }
426  int ilen=iend-ibeg+1;
427  // create temporary array for FFTW & add scratch buffer + tShift
428  int isize = 2*ilen+2*fabs(tShift)*this->rate();
429  isize = isize + (isize%4 ? 4 - isize%4 : 0); // force to be multiple of 4
430  wavearray<double> w(isize);
431  w.rate(this->rate()); w=0;
432  // copy this->data data !=0 in the middle of w array & set this->data=0
433  for(int i=0;i<ilen;i++) {w[i+isize/4]=this->data[ibeg+i];this->data[ibeg+i]=0;}
434 
435  double pi = TMath::Pi();
436  // apply time shift to waveform vector
437  w.FFTW(1);
438  TComplex C;
439  double df = w.rate()/w.size();
440  //cout << "tShift : " << tShift << endl;
441  for (int ii=0;ii<(int)w.size()/2;ii++) {
442  TComplex X(w[2*ii],w[2*ii+1]);
443  X=X*C.Exp(TComplex(0.,-2*pi*ii*df*tShift)); // Time Shift
444  w[2*ii]=X.Re();
445  w[2*ii+1]=X.Im();
446  }
447  w.FFTW(-1);
448 
449  // copy shifted data to input this->data array
450  for(int i=0;i<(int)w.size();i++) {
451  int j=ibeg-isize/4+i;
452  if((j>=0)&&(j<(int)this->size())) this->data[j]=w[i];
453  }
454 
455  return;
456 }
457 
458 //______________________________________________________________________________
459 template<class DataType_t>
461 //
462 // apply phase shift
463 //
464 // Input:
465 // pShift - phase shift (degrees)
466 //
467 
468  if(pShift==0) return;
469 
470  // search begin,end of non zero data
471  int ibeg=0; int iend=0;
472  for(int i=0;i<(int)this->size();i++) {
473  if(this->data[i]!=0 && ibeg==0) ibeg=i;
474  if(this->data[i]!=0) iend=i;
475  }
476  int ilen=iend-ibeg+1;
477  // create temporary array for FFTW & add scratch buffer + tShift
478  int isize = 2*ilen;
479  isize = isize + (isize%4 ? 4 - isize%4 : 0); // force to be multiple of 4
480  wavearray<double> w(isize);
481  w.rate(this->rate()); w=0;
482  // copy this->data data !=0 in the middle of w array & set this->data=0
483  for(int i=0;i<ilen;i++) {w[i+isize/4]=this->data[ibeg+i];this->data[ibeg+i]=0;}
484 
485  // apply phase shift to waveform vector
486  w.FFTW(1);
487  TComplex C;
488  //cout << "pShift : " << pShift << endl;
489  pShift*=TMath::Pi()/180.;
490  for (int ii=0;ii<(int)w.size()/2;ii++) {
491  TComplex X(w[2*ii],w[2*ii+1]);
492  X=X*C.Exp(TComplex(0.,-pShift)); // Phase Shift
493  w[2*ii]=X.Re();
494  w[2*ii+1]=X.Im();
495  }
496  w.FFTW(-1);
497 
498  // copy shifted data to input this->data array
499  for(int i=0;i<(int)w.size();i++) {
500  int j=ibeg-isize/4+i;
501  if((j>=0)&&(j<(int)this->size())) this->data[j]=w[i];
502  }
503 
504  return;
505 }
506 
507 template <class DataType_t>
509 {
510  int n=this->size();
511 
512  FILE *fp;
513 
514  if ( (fp = fopen(fname, "w")) == NULL ) {
515  cout << " gwavearray::DumpToFile() error: cannot open file " << fname <<". \n";
516  return;
517  };
518 
519  double dt= this->rate() ? 1./this->rate() : 1.;
520  for (int i = 0; i < n; i++) fprintf( fp,"%e\t%e\n", i*dt, this->data[i]);
521  fclose(fp);
522 }
523 
524 /*
525 //______________________________________________________________________________
526 template <class DataType_t>
527 void gwavearray<DataType_t>::Streamer(TBuffer &R__b)
528 {
529  // Stream an object of class gwavearray<DataType_t>.
530 
531  UInt_t R__s, R__c;
532  if (R__b.IsReading()) {
533  Version_t R__v = R__b.ReadVersion(&R__s, &R__c); if (R__v) { }
534  wavearray<DataType_t>::Streamer(R__b);
535  R__b >> wextern;
536  R__b.CheckByteCount(R__s, R__c, gwavearray<DataType_t>::IsA());
537  } else {
538  R__c = R__b.WriteVersion(gwavearray<DataType_t>::IsA(), kTRUE);
539  wavearray<DataType_t>::Streamer(R__b);
540  R__b << wextern;
541  R__b.SetByteCount(R__c, kTRUE);
542  }
543 }
544 */
545 
546 // instantiations
547 
548 #define CLASS_INSTANTIATION(class_) template class gwavearray< class_ >;
549 
550 CLASS_INSTANTIATION(short)
552 CLASS_INSTANTIATION(unsigned int)
554 CLASS_INSTANTIATION(long long)
555 CLASS_INSTANTIATION(float)
556 CLASS_INSTANTIATION(double)
557 
virtual size_t size() const
Definition: wavearray.hh:127
static const double C
Definition: GNGen.cc:10
int noverlap
Definition: TestDelta.C:20
void PhaseShift(double pShift=0.)
Definition: gwavearray.cc:460
void Init()
Definition: ChirpMass.C:280
void DumpToFile(char *fname)
Definition: gwavearray.cc:508
fprintf(stdout,"start=%f duration=%f rate=%f\n", x.start(), x.size()/x.rate(), x.rate())
double fHigh
virtual void rate(double r)
Definition: wavearray.hh:123
wavearray< double > a(hp.size())
par[0] name
gx Draw(GWAT_TIME)
int n
Definition: cwb_net.C:10
cout<< "skymap size : "<< L<< endl;for(int l=0;l< L;l++) sm.set(l, l);sm > const_cast< char * >("skymap.dat")
TString("c")
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
int nfact
Definition: TestDelta.C:18
wavearray< double > psd(33)
int nfft
Definition: TestDelta.C:19
STL namespace.
Long_t size
#define M
Definition: UniqSLagsList.C:3
virtual void start(double s)
Definition: wavearray.hh:119
int j
Definition: cwb_net.C:10
i drho i
int isize
void waveAssign(wavearray< Tout > &aout, wavearray< Tin > &ain)
Definition: wavearray.hh:342
double GetCentralTime()
Definition: gwavearray.cc:384
double pi
Definition: TestChirp.C:18
void Draw(double t1=0.0, double t2=0.0, double f1=0.0, double f2=0.0, double z1=0.0, double z2=0.0, int dpaletteId=DUMMY_PALETTE_ID, Option_t *option="colfz")
Definition: STFT.cc:76
wavearray< double > w
Definition: Test1.C:27
bool wextern
Definition: gwavearray.hh:90
i() int(T_cor *100))
double Pi
char fname[1024]
GWAT_DRAW
Definition: gwat.hh:10
void TimeShift(double tShift=0.)
Definition: gwavearray.cc:409
double GetTimeRange(double &tMin, double &tMax, double efraction=0.9999999)
Definition: gwavearray.cc:341
int Rate
double GetCentralTime(wavearray< double > x)
Definition: gwat.hh:14
void Init()
Definition: gwavearray.cc:60
double dt
char options[256]
CWB::STFT * DrawTF(TString options="")
Definition: gwavearray.cc:273
char title[256]
Definition: SSeriesExample.C:1
double T
Definition: testWDM_4.C:11
Definition: gwat.hh:12
Definition: gwat.hh:13
WSeries< double > tf
double fLow
double fabs(const Complex &x)
Definition: numpy.cc:37
TCanvas * GetCanvas()
Definition: STFT.hh:52
TString OS
Definition: cwb_rootlogon.C:25
double df
sprintf(tfres,"(1/%g)x(%g) (sec)x(Hz)", 2 *df, df)
DataType_t * data
Definition: wavearray.hh:301
#define CLASS_INSTANTIATION(class_)
Definition: gwavearray.cc:548
fclose(ftrig)
virtual void resize(unsigned int)
Definition: wavearray.cc:445
watplot * DrawTime(TString options="ALP", Color_t color=kBlack)
Definition: gwavearray.cc:106
double fparm
Definition: TestDelta.C:22
void Draw(GWAT_DRAW type=GWAT_TIME, TString options="ALP", Color_t color=kBlack)
Definition: gwavearray.cc:71
watplot * DrawFFT(TString options="ALP", Color_t color=kBlack)
Definition: gwavearray.cc:185
CWB::STFT * stft
Definition: ChirpMass.C:117
exit(0)