Logo coherent WaveBurst  
Library Reference Guide
Logo
watplot.hh
Go to the documentation of this file.
1 /*
2 # Copyright (C) 2019 Sergey Klimenko, Gabriele Vedovato
3 #
4 # This program is free software: you can redistribute it and/or modify
5 # it under the terms of the GNU General Public License as published by
6 # the Free Software Foundation, either version 3 of the License, or
7 # (at your option) any later version.
8 #
9 # This program is distributed in the hope that it will be useful,
10 # but WITHOUT ANY WARRANTY; without even the implied warranty of
11 # MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12 # GNU General Public License for more details.
13 #
14 # You should have received a copy of the GNU General Public License
15 # along with this program. If not, see <https://www.gnu.org/licenses/>.
16 */
17 
18 
19 /********************************************************/
20 /* Wavelet Analysis Tool */
21 /* file watplot.hh */
22 /* wat plot routines */
23 /********************************************************/
24 
25 #ifndef WATPLOT_HH
26 #define WATPLOT_HH
27 
28 #include <iostream>
29 #include <vector>
30 #include "TCanvas.h"
31 #include "TStyle.h"
32 #include "TROOT.h"
33 #include "TH1F.h"
34 #include "TH2F.h"
35 #include "TGraph.h"
36 #include "TRandom.h"
37 #include "TString.h"
38 #include "wavearray.hh"
39 #include "wseries.hh"
40 #include "skymap.hh"
41 #include "netcluster.hh"
42 
43 #include "TColor.h"
44 
45 class watplot {
46 public:
47 
48  // Default constructor
49  watplot(char* name=NULL, int=200, int=20, int=800, int=600);
50 
51  // Destructor
52  virtual ~watplot() {
53  clear();
54  if(canvas) delete canvas;
55  }
56 
57  // Clear graph objects
58  inline void clear() {
59  for(size_t i=0; i<hist1D.size(); i++) {
60  if(hist1D[i]) delete hist1D[i];
61  }
62  hist1D.clear();
63  for(size_t i=0; i<graph.size(); i++) {
64  if(graph[i]) delete graph[i];
65  }
66  graph.clear();
67  if(hist2D) { delete hist2D; hist2D = NULL; } // this instruction could produce a crash , why???
68  }
69 
70  // Clear graph objects
71  inline void null() {
72  canvas = NULL;
73  hist2D = NULL;
74  graph.clear();
75  hist1D.clear();
76  }
77 
78 
79  // Plot wavearray<double> time series
80  void plot(wavearray<double>&, char* =NULL, int=1, double=0., double=0.,
81  bool=false, float=0., float=0., bool=false, float=0., bool=false);
82 
83  // Plot wavearray<double>* time series : same parameters as plot(wavearray<double>,...)
84  void plot(wavearray<double>* ts, char* o=NULL, int c=1, double t1=0., double t2=0.,
85  bool fft=false, float f1=0., float f2=0., bool psd=false, float t3=0., bool oneside=false) {
86  plot(*ts, o, c, t1, t2, fft, f1, f2, psd, t3, oneside);
87  }
88 
89  // Plot wavearray<float>* time series : same parameters as plot(wavearray<double>,...)
90  void plot(wavearray<float>* ts, char* o=NULL, int c=1, double t1=0., double t2=0.,
91  bool fft=false, float f1=0., float f2=0., bool psd=false, float t3=0., bool oneside=false) {
93  waveAssign(tf,*ts);
94  plot(tf, o, c, t1, t2, fft, f1, f2, psd, t3, oneside);
95  }
96 
97  // Plot wavearray<float> time series : same parameters as plot(wavearray<double>,...)
98  void plot(wavearray<float> &ts, char* o=NULL, int c=1, double t1=0., double t2=0.,
99  bool fft=false, float f1=0., float f2=0., bool psd=false, float t3=0., bool oneside=false) {
101  waveAssign(tf,ts);
102  plot(tf, o, c, t1, t2, fft, f1, f2, psd, t3, oneside);
103  }
104 
105  // Plot wavearray<int>* time series : same parameters as plot(wavearray<double>,...)
106  void plot(wavearray<int>* ts, char* o=NULL, int c=1, double t1=0., double t2=0.,
107  bool fft=false, float f1=0., float f2=0., bool psd=false, float t3=0., bool oneside=false) {
109  waveAssign(tf,*ts);
110  plot(tf, o, c, t1, t2, fft, f1, f2, psd, t3, oneside);
111  }
112 
113  // Plot wavearray<int> time series : same parameters as plot(wavearray<double>,...)
114  void plot(wavearray<int> &ts, char* o=NULL, int c=1, double t1=0., double t2=0.,
115  bool fft=false, float f1=0., float f2=0., bool psd=false, float t3=0., bool oneside=false) {
117  waveAssign(tf,ts);
118  plot(tf, o, c, t1, t2, fft, f1, f2, psd, t3, oneside);
119  }
120 
121  // wavescan (2D histogram)
122  // void scan(wavearray<double>&, int, int, double=0., double=0., char* =NULL, int=256);
123 
124  // time-frequency series (2D histogram)
125 // void plot(WSeries<double>&, int=0, double=0., double=0., char* =NULL, int=256);
126 // void plot(WSeries<double>* tf, int m=0, double t1=0., double t2=0., char* o=NULL, int p=256) {
127 // plot(*tf, m, t1, t2, o, p);
128 // }
129 
130  void plot(WSeries<double>&, int=0, double=0., double=0., char* =NULL, int=256, int=0);
131  void plot(WSeries<double>* tf, int m=0, double t1=0., double t2=0., char* o=NULL, int p=256, int pid=0) {
132  plot(*tf, m, t1, t2, o, p, pid);
133  }
134 
135  void plsmooth(WSeries<double>&, int=0, double=0., double=0., char* =NULL, int=256, int=0);
136  void plsmooth(WSeries<double>* tf, int sfact=0, double t1=0., double t2=0., char* o=NULL, int p=256, int pid=0) {
137  plsmooth(*tf, sfact, t1, t2, o, p, pid);
138  }
139 
140  // plot skymaps
141  void plot(skymap&, char* = NULL, int=256);
142 
143  // monster event display (multi resolution analysis)
144  // pwc : pointer to netcluster
145  // cid : event cluster id
146  // nifo : number of detectors
147  // type : L/N likelihood/null
148  // irate : rate to be plotted - if 0 then select all rates
149  // opt : draw options
150  // pal : draw palette
151  // wp : wavelet packet
152  void plot(netcluster* pwc, int cid, int nifo, char type='L', int irate=0, char* opt=NULL, int pal=256, bool wp=false);
153 
154  // chirp event display
155  // pwc : pointer to clusterdata
156  // inj_mchirp : chirp mass of injected signal
157  void plot(clusterdata* pcd, double inj_mchirp);
158 
159  // set palette
160  void SetPlotStyle(int paletteId, int NCont=255);
161  // return max val
162  double getmax(WSeries<double> &w, double t1, double t2);
163  // return max value of WSeries* w object
164  double getmax(WSeries<double>* tf, double t1=0., double t2=0.) {
165  return getmax(*tf, t1, t2);
166  }
167 
168  void print(TString fname);
170  void goptions(char* opt=NULL, int col=1, double t1=0., double t2=0.,
171  bool fft=false, float f1=0., float f2=0., bool psd=false, float t3=0., bool oneside=false);
172  void gtitle(TString title="", TString xtitle="", TString ytitle="");
173 
174  // return value according to opt option
175  double procOpt(int opt, double val00, double val90=0.) {
176  switch(opt) {
177  case 1: return val00*val00;
178  case 2: return fabs(val00);
179  case -1: return (val00*val00+val90*val90)/2.;
180  case -2: return sqrt((val00*val00+val90*val90)/2.);
181  default: return val00+val90;
182  }
183  }
184 
185  // return graph pointer
186  TGraph* getGraph(int n){
187  return n<graph.size() ? graph[n] : graph[0];
188  }
189 
190  // data members
191 
192  TCanvas* canvas; // pointer to TCanvas object
193  TH2F* hist2D; // pointer to TH2F object
194  std::vector<TGraph*> graph; // vector of pointers to TGraph objects
195  std::vector<TH1F*> hist1D; // vector of pointers to TH1F objects
196 
198  TString title; // graph title
199  TString xtitle; // x axis name
200  TString ytitle; // y axis name
201  int ncol; // color index of TGraph plots
202  TString opt; // TGraph::Draw options
203  int col; // TGraph line color
204  double t1; // start of time interval in seconds
205  double t2; // end of time interval in seconds
206  bool fft; // true -> plot fft
207  double f1; // set begin frequency (Hz)
208  double f2; // set end frequency (Hz)
209  bool psd; // true -> plot psd using blackmanharris window
210  double t3; // is the chunk length (sec) used to produce the psd
211  bool oneside; //
212 
213  // used by THtml doc
214  ClassDef(watplot,1)
215 };
216 
217 // put operator
219 // get operator
221 // print operator
223 char* operator >> (watplot& graph, char* fname);
224 
225 #endif // WATPLOT_HH
double f1
Definition: watplot.hh:207
void gtitle(TString title="", TString xtitle="", TString ytitle="")
Definition: watplot.cc:1256
double getmax(WSeries< double > *tf, double t1=0., double t2=0.)
Definition: watplot.hh:164
bool oneside
Definition: watplot.hh:211
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:98
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:106
par [0] name
int n
Definition: cwb_net.C:28
TString("c")
double t3
Definition: watplot.hh:210
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<< "\ layers : "<< nLAYERS<< "\ dF(hz) : "<< dF<< "\ 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:198
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:90
netcluster * pwc
Definition: cwb_job_obj.C:38
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:84
std::vector< TGraph * > graph
Definition: watplot.hh:194
double procOpt(int opt, double val00, double val90=0.)
Definition: watplot.hh:175
int m
Definition: cwb_net.C:28
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:150
watplot(char *name=NULL, int=200, int=20, int=800, int=600)
Definition: watplot.cc:86
TGraph * getGraph(int n)
Definition: watplot.hh:186
void waveAssign(wavearray< Tout > &aout, wavearray< Tin > &ain)
Definition: wavearray.hh:360
bool psd
Definition: watplot.hh:209
TH2F * hist2D
Definition: watplot.hh:193
void clear()
Definition: watplot.hh:58
wavearray< double > w
Definition: Test1.C:27
TCanvas * canvas
Definition: watplot.hh:192
bool fft
Definition: watplot.hh:206
double t2
Definition: watplot.hh:205
#define NCont
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:131
double getmax(WSeries< double > &w, double t1, double t2)
Definition: watplot.cc:431
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:1221
TString xtitle
Definition: watplot.hh:199
void plsmooth(WSeries< double > &, int=0, double=0., double=0., char *=NULL, int=256, int=0)
Definition: watplot.cc:284
Definition: skymap.hh:63
int col
Definition: watplot.hh:203
double f2
Definition: watplot.hh:208
wavearray< double > data
Definition: watplot.hh:197
std::vector< TH1F * > hist1D
Definition: watplot.hh:195
int nifo
WSeries< double > tf
wavearray< double > ts(N)
void print(TString fname)
Definition: watplot.cc:1197
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:114
TString ytitle
Definition: watplot.hh:200
double fabs(const Complex &x)
Definition: numpy.cc:55
void blackmanharris(wavearray< double > &w)
Definition: watplot.cc:1171
void null()
Definition: watplot.hh:71
void SetPlotStyle(int paletteId, int NCont=255)
Definition: watplot.cc:1082
virtual ~watplot()
Definition: watplot.hh:52
double t1
Definition: watplot.hh:204
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:136
wavearray< double > & operator>>(watplot &graph, wavearray< double > &x)
Definition: watplot.cc:1276
TString opt
Definition: watplot.hh:202
int ncol
Definition: watplot.hh:201