Logo Coherent WaveBurst  
Reference Guide
Logo
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
gwavearray.hh
Go to the documentation of this file.
1 /**********************************************************
2  * Package: graphic wat Class Library
3  * File name: gwavearray.hh
4  * Author: Gabriele Vedovato (vedovato@lnl.infn.it)
5  **********************************************************/
6 
7 #ifndef GWAVEARRAY_HH
8 #define GWAVEARRAY_HH
9 
10 #ifndef WAVEARRAY_HH
11 #include "wavearray.hh"
12 #endif
13 
14 #include "gwat.hh"
15 #include "watplot.hh"
16 #include "time.hh"
17 #include "STFT.hh"
18 #include "TComplex.h"
19 
20 
21 template<class DataType_t>
22 class gwavearray : public wavearray<DataType_t>
23 {
24 public:
25 
26  gwavearray() : wavearray<DataType_t>() {Init();wextern=false;} // Default constructor
27  gwavearray(const wavearray<DataType_t>&); // copy Constructor
29 
30  virtual ~gwavearray(); // Destructor
31 
32  void Draw(GWAT_DRAW type=GWAT_TIME,
33  TString options="ALP", Color_t color=kBlack);
34 
35  void PlotTime(double edge=0) {
36  SetTimeRange(this->start()+edge,this->start()+this->size()/this->rate()-edge);
37  DrawTime("CUSTOM ALP");} // *MENU*
38  void PlotFFT(double edge=0) {
39  SetTimeRange(this->start()+edge,this->start()+this->size()/this->rate()-edge);
40  DrawFFT("CUSTOM ALP");} // *MENU*
41  void PlotPSD(double edge=0) {
42  SetTimeRange(this->start()+edge,this->start()+this->size()/this->rate()-edge);
43  DrawFFT("CUSTOM ALP PSD");} // *MENU*
44  void PlotTF(double edge=0) {DrawTF();;}
45 
46  watplot* DrawTime(TString options="ALP", Color_t color=kBlack);
47  watplot* DrawFFT(TString options="ALP", Color_t color=kBlack);
49 
51  TString options="ALP", Color_t color=kBlack);
52 
53  watplot* DrawTime(wavearray<DataType_t>* x, TString options="ALP", Color_t color=kBlack);
54  watplot* DrawFFT(wavearray<DataType_t>* x, TString options="ALP", Color_t color=kBlack);
56 
57  void SetTimeRange(double tMin, double tMax) {
58  if(tMin>=tMax) return;
59  double tStart=this->start();
60  double tStop=this->start()+this->size()/this->rate();
61  this->tMin = tMin<tStart||tMin>tStop ? tStart : tMin;
62  this->tMax = tMax<tStart||tMax>tStop ? tStop : tMax;
63  }
64  double GetTimeRange(double& tMin, double& tMax, double efraction = 0.9999999);
65  double GetCentralTime();
66  void TimeShift(double tShift=0.);
67  void PhaseShift(double pShift=0.);
68 
69  CWB::STFT* GetSTFT() {return this->stft;}
70  watplot* GetWATPLOT() {return this->pts;}
71 
72  // print wseries parameters
74  virtual void Browse(TBrowser *b) {PlotTime();}
75  void SetWATPLOT(watplot* pts) {this->pts=pts;}
76 
77  void SetComment(TString comment) {this->comment=comment;}
79  void PrintComment() {cout << comment << endl;} // *MENU*
80 
81  void DumpToFile(char* fname); // *MENU*
82 
83 private:
84 
86  watplot* pts; //!
87 
88  void Init();
89 
90  bool wextern; // external wavearray
91 
92  TString comment; // store comment
93 
94  double tMin;
95  double tMax;
96  double tSave; // start time : used to manage DrawTime with option "SAME"
97 
98  ClassDef(gwavearray,3)
99 };
100 
101 #endif
virtual size_t size() const
Definition: wavearray.hh:127
void PlotTime(double edge=0)
Definition: gwavearray.hh:35
void PhaseShift(double pShift=0.)
Definition: gwavearray.cc:460
void PrintComment()
Definition: gwavearray.hh:79
void DumpToFile(char *fname)
Definition: gwavearray.cc:508
void PlotPSD(double edge=0)
Definition: gwavearray.hh:41
double tSave
Definition: gwavearray.hh:96
virtual double start() const
Definition: wavearray.hh:120
void print()
Definition: wavearray.cc:2203
TString comment
Definition: gwavearray.hh:92
TString("c")
virtual double edge() const
Definition: wavearray.hh:126
void PlotTF(double edge=0)
Definition: gwavearray.hh:44
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
CWB::STFT * GetSTFT()
Definition: gwavearray.hh:69
void PlotFFT(double edge=0)
Definition: gwavearray.hh:38
double GetCentralTime()
Definition: gwavearray.cc:384
CWB::STFT * stft
Definition: gwavearray.hh:85
bool wextern
Definition: gwavearray.hh:90
double tMin
Definition: gwavearray.hh:94
watplot * pts
Definition: gwavearray.hh:86
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
void SetComment(TString comment)
Definition: gwavearray.hh:77
void SetTimeRange(double tMin, double tMax)
Definition: gwavearray.hh:57
virtual ~gwavearray()
Definition: gwavearray.cc:30
void SetWATPLOT(watplot *pts)
Definition: gwavearray.hh:75
virtual void Browse(TBrowser *b)
Definition: gwavearray.hh:74
void Init()
Definition: gwavearray.cc:60
TString GetComment()
Definition: gwavearray.hh:78
char options[256]
CWB::STFT * DrawTF(TString options="")
Definition: gwavearray.cc:273
virtual double rate() const
Definition: wavearray.hh:124
watplot * GetWATPLOT()
Definition: gwavearray.hh:70
watplot * DrawTime(TString options="ALP", Color_t color=kBlack)
Definition: gwavearray.cc:106
double tMax
Definition: gwavearray.hh:95
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
void PrintParameters()
Definition: gwavearray.hh:73