Logo Coherent WaveBurst  
Reference Guide
Logo
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
watplot.hh
Go to the documentation of this file.
1 /********************************************************/
2 /* Wavelet Analysis Tool */
3 /* file watplot.hh */
4 /* wat plot routines */
5 /********************************************************/
6 
7 #ifndef WATPLOT_HH
8 #define WATPLOT_HH
9 
10 #include <iostream>
11 #include <vector>
12 #include "TCanvas.h"
13 #include "TStyle.h"
14 #include "TROOT.h"
15 #include "TH1F.h"
16 #include "TH2F.h"
17 #include "TGraph.h"
18 #include "TRandom.h"
19 #include "TString.h"
20 #include "wavearray.hh"
21 #include "wseries.hh"
22 #include "skymap.hh"
23 #include "netcluster.hh"
24 
25 #include "TColor.h"
26 
27 class watplot {
28 public:
29 
30  // Default constructor
31  watplot(char* name=NULL, int=200, int=20, int=800, int=600);
32 
33  // Destructor
34  virtual ~watplot() {
35  clear();
36  if(canvas) delete canvas;
37  }
38 
39  // Clear graph objects
40  inline void clear() {
41  for(size_t i=0; i<hist1D.size(); i++) {
42  if(hist1D[i]) delete hist1D[i];
43  }
44  hist1D.clear();
45  for(size_t i=0; i<graph.size(); i++) {
46  if(graph[i]) delete graph[i];
47  }
48  graph.clear();
49  if(hist2D) { delete hist2D; hist2D = NULL; } // this instruction could produce a crash , why???
50  }
51 
52  // Clear graph objects
53  inline void null() {
54  canvas = NULL;
55  hist2D = NULL;
56  graph.clear();
57  hist1D.clear();
58  }
59 
60 
61  // Plot wavearray<double> time series
62  void plot(wavearray<double>&, char* =NULL, int=1, double=0., double=0.,
63  bool=false, float=0., float=0., bool=false, float=0., bool=false);
64 
65  // Plot wavearray<double>* time series : same parameters as plot(wavearray<double>,...)
66  void plot(wavearray<double>* ts, char* o=NULL, int c=1, double t1=0., double t2=0.,
67  bool fft=false, float f1=0., float f2=0., bool psd=false, float t3=0., bool oneside=false) {
68  plot(*ts, o, c, t1, t2, fft, f1, f2, psd, t3, oneside);
69  }
70 
71  // Plot wavearray<float>* time series : same parameters as plot(wavearray<double>,...)
72  void plot(wavearray<float>* ts, char* o=NULL, int c=1, double t1=0., double t2=0.,
73  bool fft=false, float f1=0., float f2=0., bool psd=false, float t3=0., bool oneside=false) {
75  waveAssign(tf,*ts);
76  plot(tf, o, c, t1, t2, fft, f1, f2, psd, t3, oneside);
77  }
78 
79  // Plot wavearray<float> time series : same parameters as plot(wavearray<double>,...)
80  void plot(wavearray<float> &ts, char* o=NULL, int c=1, double t1=0., double t2=0.,
81  bool fft=false, float f1=0., float f2=0., bool psd=false, float t3=0., bool oneside=false) {
83  waveAssign(tf,ts);
84  plot(tf, o, c, t1, t2, fft, f1, f2, psd, t3, oneside);
85  }
86 
87  // Plot wavearray<int>* time series : same parameters as plot(wavearray<double>,...)
88  void plot(wavearray<int>* ts, char* o=NULL, int c=1, double t1=0., double t2=0.,
89  bool fft=false, float f1=0., float f2=0., bool psd=false, float t3=0., bool oneside=false) {
91  waveAssign(tf,*ts);
92  plot(tf, o, c, t1, t2, fft, f1, f2, psd, t3, oneside);
93  }
94 
95  // Plot wavearray<int> time series : same parameters as plot(wavearray<double>,...)
96  void plot(wavearray<int> &ts, char* o=NULL, int c=1, double t1=0., double t2=0.,
97  bool fft=false, float f1=0., float f2=0., bool psd=false, float t3=0., bool oneside=false) {
99  waveAssign(tf,ts);
100  plot(tf, o, c, t1, t2, fft, f1, f2, psd, t3, oneside);
101  }
102 
103  // wavescan (2D histogram)
104  // void scan(wavearray<double>&, int, int, double=0., double=0., char* =NULL, int=256);
105 
106  // time-frequency series (2D histogram)
107 // void plot(WSeries<double>&, int=0, double=0., double=0., char* =NULL, int=256);
108 // void plot(WSeries<double>* tf, int m=0, double t1=0., double t2=0., char* o=NULL, int p=256) {
109 // plot(*tf, m, t1, t2, o, p);
110 // }
111 
112  void plot(WSeries<double>&, int=0, double=0., double=0., char* =NULL, int=256, int=0);
113  void plot(WSeries<double>* tf, int m=0, double t1=0., double t2=0., char* o=NULL, int p=256, int pid=0) {
114  plot(*tf, m, t1, t2, o, p, pid);
115  }
116 
117  void plsmooth(WSeries<double>&, int=0, double=0., double=0., char* =NULL, int=256, int=0);
118  void plsmooth(WSeries<double>* tf, int sfact=0, double t1=0., double t2=0., char* o=NULL, int p=256, int pid=0) {
119  plsmooth(*tf, sfact, t1, t2, o, p, pid);
120  }
121 
122  // plot skymaps
123  void plot(skymap&, char* = NULL, int=256);
124 
125  // monster event display (multi resolution analysis)
126  // pwc : pointer to netcluster
127  // cid : event cluster id
128  // nifo : number of detectors
129  // type : L/N likelihood/null
130  // irate : rate to be plotted - if 0 then select all rates
131  // opt : draw options
132  // pal : draw palette
133  // wp : wavelet packet
134  void plot(netcluster* pwc, int cid, int nifo, char type='L', int irate=0, char* opt=NULL, int pal=256, bool wp=false);
135 
136  // chirp event display
137  // pwc : pointer to clusterdata
138  // inj_mchirp : chirp mass of injected signal
139  void plot(clusterdata* pcd, double inj_mchirp);
140 
141  // set palette
142  void SetPlotStyle(int paletteId, int NCont=255);
143  // return max val
144  double getmax(WSeries<double> &w, double t1, double t2);
145  // return max value of WSeries* w object
146  double getmax(WSeries<double>* tf, double t1=0., double t2=0.) {
147  return getmax(*tf, t1, t2);
148  }
149 
150  void print(TString fname);
152  void goptions(char* opt=NULL, int col=1, double t1=0., double t2=0.,
153  bool fft=false, float f1=0., float f2=0., bool psd=false, float t3=0., bool oneside=false);
154  void gtitle(TString title="", TString xtitle="", TString ytitle="");
155 
156  // return value according to opt option
157  double procOpt(int opt, double val00, double val90=0.) {
158  switch(opt) {
159  case 1: return val00*val00;
160  case 2: return fabs(val00);
161  case -1: return (val00*val00+val90*val90)/2.;
162  case -2: return sqrt((val00*val00+val90*val90)/2.);
163  default: return val00+val90;
164  }
165  }
166 
167  // return graph pointer
168  TGraph* getGraph(int n){
169  return n<graph.size() ? graph[n] : graph[0];
170  }
171 
172  // data members
173 
174  TCanvas* canvas; // pointer to TCanvas object
175  TH2F* hist2D; // pointer to TH2F object
176  std::vector<TGraph*> graph; // vector of pointers to TGraph objects
177  std::vector<TH1F*> hist1D; // vector of pointers to TH1F objects
178 
180  TString title; // graph title
181  TString xtitle; // x axis name
182  TString ytitle; // y axis name
183  int ncol; // color index of TGraph plots
184  TString opt; // TGraph::Draw options
185  int col; // TGraph line color
186  double t1; // start of time interval in seconds
187  double t2; // end of time interval in seconds
188  bool fft; // true -> plot fft
189  double f1; // set begin frequency (Hz)
190  double f2; // set end frequency (Hz)
191  bool psd; // true -> plot psd using blackmanharris window
192  double t3; // is the chunk length (sec) used to produce the psd
193  bool oneside; //
194 
195  // used by THtml doc
196  ClassDef(watplot,1)
197 };
198 
199 // put operator
201 // get operator
203 // print operator
205 char* operator >> (watplot& graph, char* fname);
206 
207 #endif // WATPLOT_HH
void graph(TString ifile)
Definition: Average.C:543
double f1
Definition: watplot.hh:189
void gtitle(TString title="", TString xtitle="", TString ytitle="")
Definition: watplot.cc:1237
double getmax(WSeries< double > *tf, double t1=0., double t2=0.)
Definition: watplot.hh:146
bool oneside
Definition: watplot.hh:193
void plot(wavearray< float > &ts, char *o=NULL, int c=1, double t1=0., double t2=0., bool fft=false, float f1=0., float f2=0., bool psd=false, float t3=0., bool oneside=false)
Definition: watplot.hh:80
int irate
void plot(wavearray< int > *ts, char *o=NULL, int c=1, double t1=0., double t2=0., bool fft=false, float f1=0., float f2=0., bool psd=false, float t3=0., bool oneside=false)
Definition: watplot.hh:88
par[0] name
int n
Definition: cwb_net.C:10
TString("c")
double t3
Definition: watplot.hh:192
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
TString title
Definition: watplot.hh:180
void plot(wavearray< float > *ts, char *o=NULL, int c=1, double t1=0., double t2=0., bool fft=false, float f1=0., float f2=0., bool psd=false, float t3=0., bool oneside=false)
Definition: watplot.hh:72
netcluster * pwc
Definition: cwb_job_obj.C:20
void plot(wavearray< double > *ts, char *o=NULL, int c=1, double t1=0., double t2=0., bool fft=false, float f1=0., float f2=0., bool psd=false, float t3=0., bool oneside=false)
Definition: watplot.hh:66
std::vector< TGraph * > graph
Definition: watplot.hh:176
double procOpt(int opt, double val00, double val90=0.)
Definition: watplot.hh:157
int m
Definition: cwb_net.C:10
i drho i
void plot(wavearray< double > &, char *=NULL, int=1, double=0., double=0., bool=false, float=0., float=0., bool=false, float=0., bool=false)
Definition: watplot.cc:132
watplot(char *name=NULL, int=200, int=20, int=800, int=600)
Definition: watplot.cc:68
TGraph * getGraph(int n)
Definition: watplot.hh:168
void waveAssign(wavearray< Tout > &aout, wavearray< Tin > &ain)
Definition: wavearray.hh:342
bool psd
Definition: watplot.hh:191
TH2F * hist2D
Definition: watplot.hh:175
#define NCont
void clear()
Definition: watplot.hh:40
wavearray< double > w
Definition: Test1.C:27
TCanvas * canvas
Definition: watplot.hh:174
bool fft
Definition: watplot.hh:188
double t2
Definition: watplot.hh:187
char fname[1024]
void plot(WSeries< double > *tf, int m=0, double t1=0., double t2=0., char *o=NULL, int p=256, int pid=0)
Definition: watplot.hh:113
double getmax(WSeries< double > &w, double t1, double t2)
Definition: watplot.cc:413
void goptions(char *opt=NULL, int col=1, double t1=0., double t2=0., bool fft=false, float f1=0., float f2=0., bool psd=false, float t3=0., bool oneside=false)
Definition: watplot.cc:1202
TString xtitle
Definition: watplot.hh:181
void plsmooth(WSeries< double > &, int=0, double=0., double=0., char *=NULL, int=256, int=0)
Definition: watplot.cc:266
Definition: skymap.hh:45
int col
Definition: watplot.hh:185
double f2
Definition: watplot.hh:190
wavearray< double > data
Definition: watplot.hh:179
std::vector< TH1F * > hist1D
Definition: watplot.hh:177
int nifo
WSeries< double > tf
wavearray< double > ts(N)
void print(TString fname)
Definition: watplot.cc:1178
void plot(wavearray< int > &ts, char *o=NULL, int c=1, double t1=0., double t2=0., bool fft=false, float f1=0., float f2=0., bool psd=false, float t3=0., bool oneside=false)
Definition: watplot.hh:96
TString ytitle
Definition: watplot.hh:182
double fabs(const Complex &x)
Definition: numpy.cc:37
void blackmanharris(wavearray< double > &w)
Definition: watplot.cc:1152
void null()
Definition: watplot.hh:53
void SetPlotStyle(int paletteId, int NCont=255)
Definition: watplot.cc:1063
virtual ~watplot()
Definition: watplot.hh:34
double t1
Definition: watplot.hh:186
void plsmooth(WSeries< double > *tf, int sfact=0, double t1=0., double t2=0., char *o=NULL, int p=256, int pid=0)
Definition: watplot.hh:118
wavearray< double > & operator>>(watplot &graph, wavearray< double > &x)
Definition: watplot.cc:1257
TString opt
Definition: watplot.hh:184
int ncol
Definition: watplot.hh:183