Logo Coherent WaveBurst  
Reference Guide
Logo
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
DrawWaveformPE.C
Go to the documentation of this file.
1 // this macro shows how to read the PE results from ROOT and ASCII file
2 // plot the waveforms with tghe same format presented in the PE CED
3 // Author : G.Vedovato
4 
11 
12 std::vector<TString> ReadDataFromROOT(TString ifile, int ifo, wavearray<double> *winj, wavearray<double> *wrec, wavearray<double> *wwht,
15 
19  TString pdir, double P, bool freq=false, bool showerr=true);
20 
21 void DumpWaveform(int ifo, wavearray<double>* wrec, wavearray<double>* wmed,
24 
25 
26 double GetTimeBoundaries(wavearray<double> x, double P, double& bT, double& eT);
29 
30 void FrequencyCut(wavearray<double>* x, double bF, double eF);
31 void PlotSpectrogram(wavearray<double>* x, double tstart=0, double tstop=0, TString title="", TString ofname="", double tchirp=0., double mchirp=0.);
32 
33 
34 #define IN_CWB_ASCII_FILE_L1 "L1_pe_wave.dat"
35 #define IN_CWB_ASCII_FILE_H1 "H1_pe_wave.dat"
36 
37 #define OUT_CWB_ASCII_FILE_L1 "OUT_L1_pe_wave.dat"
38 #define OUT_CWB_ASCII_FILE_H1 "OUT_H1_pe_wave.dat"
39 
40 
41 //#define FREQUENCY_CUT
42 #define FLOW 16
43 #define FHIGH 256
44 
46 TF1* fchirp;
47 
48 void DrawWaveformPE(TString ipath, int gtype=0, int ifo=0, double tshift=0, TString label="PLOT", double P=0.99) {
49 //
50 // ipath- > if(gtype<0) input ced idir
51 // if(gtype>0) input root file name
52 //
53 // gtype -> plot type
54 // < 0 -> plot instantaneous frequency : must be used only when ifo[0]=L1 and ifo[1]=H1
55 // 0 -> plot reconstructed vs whitened
56 // 1 -> plot reconstructed vs median
57 // 2 -> plot rec vs inj (only if injection is present in the input file)
58 // 3 -> plot the difference between med/rec and inj waveforms
59 // 4 -> plot spectrogram ifo
60 // 5 -> plot spectrogram ifo_0 + ifo_1(inv+tshift) : must be used only when ifo[0]=L1 and ifo[1]=H1
61 // 6 -> plot spectrogram ifo_0 - ifo_1(inv+tshift) : must be used only when ifo[0]=L1 and ifo[1]=H1
62 // >=7 -> dump data to ascii file
63 //
64 // ifo -> detector id : 0,1,...
65 //
66 // tshift-> the second detector is shifted by tshift and sign inverted : must be used only when ifo[0]=L1 and ifo[1]=H1
67 // GW150914 -> tshift =-0.0069
68 // GW151226 -> tshift =-0.0011
69 // LVT151012 -> tshift = 0.001
70 //
71 // label -> used for title and output plot files
72 //
73 // P -> [0,1] is used to restrict the plot time interval around the wave energy percentage P
74 //
75 
76  bool exit = label.Contains("EXIT") ? true : false;
77  label.ReplaceAll("EXIT","");
78 
79  fchirp=NULL;
80  stft=NULL;
81 
82  double tchirp=0.;
83  double mchirp=0.;
84 
85  if(label=="GW150914") {tshift=-0.0069;tchirp=8.42;mchirp=30.0;}
86  if(label=="GW151226") {tshift=-0.0011;tchirp=7.65;mchirp=8.9;}
87  if(label=="LVT151012") {tshift= 0.0010;tchirp=7.43;mchirp=24.0;} // ?
88 
90 
98 
105 
106  int nifo=0;
107  std::vector<TString> ifoname;
108 
109  if(gtype<0) {
110  // read data from PE ascii output files
111  ReadDataFromASCII(ipath, ifo, &wrec[ifo], &wmed[ifo], &wl50[ifo], &wu50[ifo], &wl90[ifo], &wu90[ifo],
112  &frec[ifo], &fmed[ifo], &fl50[ifo], &fu50[ifo], &fl90[ifo], &fu90[ifo]);
113  PlotWaveformAsymmErrors("", "", &frec[ifo], &fmed[ifo], &fl50[ifo], &fu50[ifo], &fl90[ifo],
114  &fu90[ifo], &wrec[ifo], ".", P, true);
115  return;
116  } else {
117  // read data from PE root output file
118  for(int n=0;n<NIFO_MAX;n++) {
119  ifoname = ReadDataFromROOT(ipath, n, &winj[n], &wrec[n], &wwht[n], &wmed[n], &wl50[n], &wu50[n], &wl90[n], &wu90[n]);
120  nifo = ifoname.size();
121  if(n>=nifo-1) break;
122  }
123  if(ifo>nifo-1) {cout << "DrawWaveformPE Error : input ifo=" << ifo << " parameter must be < nifo=" << nifo << endl;exit(1);}
124 #ifdef FREQUENCY_CUT
125  // frequency cut for whitened data
126  for(int n=0;n<nifo;n++) FrequencyCut(&wwht[n], FLOW, FHIGH);
127 #endif
128  }
129 
130  // compute med,start,stop times for input wavearray
131  double tmed = wrec[0].start()+(wrec[0].size()/wrec[0].rate())/2.;
132  double tstart = tmed-1;
133  double tstop = tmed+1;
134  GetTimeBoundaries(wrec[0], P, tstart, tstop);
135  //cout << "Time Boundaries : " << tstart-wrec[0].start() << " " << tstop-wrec[0].start() << endl;exit(0);
136 
137  if(gtype==0) {
138  // plot reconstructed vs whitened
139  PlotWaveformAsymmErrors("", "", &wrec[ifo], &wwht[ifo], &wl50[ifo], &wu50[ifo], &wl90[ifo], &wu90[ifo], &wrec[ifo], "", P, false, false);
140  return;
141  }
142  if(gtype==1) {
143  // plot reconstructed vs median
144  PlotWaveformAsymmErrors("", "", &wrec[ifo], &wmed[ifo], &wl50[ifo], &wu50[ifo], &wl90[ifo], &wu90[ifo], &wrec[ifo], "", P);
145  return;
146  }
147  if(gtype==2) {
148  // plot reconstructed vs injected
149  PlotWaveformAsymmErrors("", "", &winj[ifo], &wmed[ifo], &wl50[ifo], &wu50[ifo], &wl90[ifo], &wu90[ifo], &wrec[ifo], "", P);
150  return;
151  }
152  if(gtype==3) {
153  // plot difference between med/rec and inj
155  wdif[ifo] = GetDifWaveform(&wmed[ifo], &winj[ifo]);
156  PlotWaveformAsymmErrors("", "", &wdif[ifo], &wdif[ifo], &wl50[ifo], &wu50[ifo], &wl90[ifo], &wu90[ifo], &wrec[ifo], "", P);
157  return;
158  }
159  if(gtype==4) {
160  // plot spectrogram ifo
161  TString title = TString::Format("%s whitened data %s",label.Data(),ifoname[ifo].Data());
162  TString ofile = TString::Format("%s_whitened_data_%s_Spectrogram.png",label.Data(),ifoname[ifo].Data());
163  PlotSpectrogram(&wwht[ifo], tstart, tstop,title,ofile,tchirp,mchirp);
164  if(exit) exit(0); else return;
165  }
166  if(gtype==5) {
167  // plot spectrogram ifo_0 + ifo_1(inv+tshift)
168  TString title = TString::Format("%s whitened data %s + %s(inv/tshift)",label.Data(),ifoname[0].Data(),ifoname[1].Data());
169  TString ofile = TString::Format("%s_whitened_data_%s_plus_%s_Spectrogram.png",label.Data(),ifoname[0].Data(),ifoname[1].Data());
170  // invert sign of the second ifo
171  wwht[1]*=-1;
172  // shift time of the second ifo
173  CWB::mdc::TimeShift(wwht[1],tshift);
174  // add first and second ifo
175  wwht[0]+=wwht[1];
176  PlotSpectrogram(&wwht[ifo], tstart, tstop,title,ofile,tchirp,mchirp);
177  if(exit) exit(0); else return;
178  }
179  if(gtype==6) {
180  // plot spectrogram ifo_0 - ifo_1(inv+tshift)
181  TString title = TString::Format("%s whitened data %s - %s(inv/tshift)",label.Data(),ifoname[0].Data(),ifoname[1].Data());
182  TString ofile = TString::Format("%s_whitened_data_%s_minus_%s_Spectrogram.png",label.Data(),ifoname[0].Data(),ifoname[1].Data());
183  // invert sign of the second ifo
184  wwht[1]*=-1;
185  // shift time of the second ifo
186  CWB::mdc::TimeShift(wwht[1],tshift);
187  // subctract second from the first ifo
188  wwht[0]-=wwht[1];
189  PlotSpectrogram(&wwht[ifo], tstart, tstop,title,ofile,tchirp,mchirp);
190  if(exit) exit(0); else return;
191  }
192  if(gtype>=7) {
193  // dump data to ascii file
194  DumpWaveform(ifo, &wrec[ifo], &wmed[ifo], &wl50[ifo], &wu50[ifo], &wl90[ifo], &wu90[ifo]);
195  exit(0);
196  }
197 
198  return;
199 }
200 
203  wavearray<double> *wl90, wavearray<double> *wu90) {
204 //
205 // ----------------------------------------------------
206 // ROOT Output PE Parameters
207 // ----------------------------------------------------
208 
209  float pe_erR[11]; // probability distribution of residuals
210  float pe_erF[11]; // probability distribution of frequency residuals
211  int pe_trials; // number of effective trials
212  float pe_pnul[2*NIFO_MAX]; // null pixel statistic, for each detector pnul[0]=avr, pnul=rms
213  float pe_snet[2]; // SNRnet statistic, 0 -> avr, 1 -> rms
214  float pe_ff[2]; // Fitting Factor statistic, 0 -> avr, 1 -> rms
215  float pe_of[2]; // Overlap Factor statistic, 0 -> avr, 1 -> rms
216  float pe_mch[2]; // chirp mass statistic, 0 -> avr, 1 -> rms
217 
218  wavearray<double>* pe_wINJ[NIFO_MAX];
219  wavearray<double>* pe_wREC[NIFO_MAX];
220  wavearray<double>* pe_wWHT[NIFO_MAX];
221  wavearray<double>* pe_wMED[NIFO_MAX];
222  wavearray<double>* pe_wL50[NIFO_MAX];
223  wavearray<double>* pe_wU50[NIFO_MAX];
224  wavearray<double>* pe_wL90[NIFO_MAX];
225  wavearray<double>* pe_wU90[NIFO_MAX];
226  wavearray<double>* pe_wAVR[NIFO_MAX];
227  wavearray<double>* pe_wRMS[NIFO_MAX];
228 
229  TFile* froot = new TFile(ifile,"READ");
230  if(froot==NULL) {cout << "ReadWaveformPE Error : opening input root file" << endl;exit(1);}
231  TTree* gTREE = froot->Get("waveburst");
232  if(gTREE==NULL) {cout << "ReadWaveformPE Error : no waveburst present in the file" << endl;exit(1);}
233 
234  // get detector list
235  TList* list = gTREE->GetUserInfo();
236  int nIFO=list->GetSize();
237  std::vector<TString> ifoname;
238  if (nIFO==0) {cout << "ReadWaveformPE Error : no ifo present in the tree" << endl;exit(1);}
239  for (int n=0;n<list->GetSize();n++) {
240  detector* pDetector;
241  pDetector = (detector*)list->At(n);
242  ifoname.push_back(pDetector->Name);
243  detectorParams dParams = pDetector->getDetectorParams();
244  //pDetector->print();
245  }
246 
247  for(int n=0;n<nIFO;n++) {
248  pe_wINJ[n] = new wavearray<double>;
249  pe_wREC[n] = new wavearray<double>;
250  pe_wWHT[n] = new wavearray<double>;
251  pe_wMED[n] = new wavearray<double>;
252  pe_wL50[n] = new wavearray<double>;
253  pe_wU50[n] = new wavearray<double>;
254  pe_wL90[n] = new wavearray<double>;
255  pe_wU90[n] = new wavearray<double>;
256  pe_wAVR[n] = new wavearray<double>;
257  pe_wRMS[n] = new wavearray<double>;
258  }
259 
260  gTREE->SetBranchAddress("pe_trials",&pe_trials);
261  gTREE->SetBranchAddress("pe_erR",pe_erR);
262  gTREE->SetBranchAddress("pe_erF",pe_erF);
263  gTREE->SetBranchAddress("pe_pnul",pe_pnul);
264  gTREE->SetBranchAddress("pe_snet",pe_snet);
265  gTREE->SetBranchAddress("pe_ff",pe_ff);
266  gTREE->SetBranchAddress("pe_of",pe_of);
267  gTREE->SetBranchAddress("pe_mch",pe_mch);
268 
269  for(int n=0;n<nIFO;n++) {
270  gTREE->SetBranchAddress(TString::Format("pe_wINJ_%d",n).Data(),&pe_wINJ[n]);
271  gTREE->SetBranchAddress(TString::Format("pe_wREC_%d",n).Data(),&pe_wREC[n]);
272  gTREE->SetBranchAddress(TString::Format("pe_wWHT_%d",n).Data(),&pe_wWHT[n]);
273  gTREE->SetBranchAddress(TString::Format("pe_wMED_%d",n).Data(),&pe_wMED[n]);
274  gTREE->SetBranchAddress(TString::Format("pe_wL50_%d",n).Data(),&pe_wL50[n]);
275  gTREE->SetBranchAddress(TString::Format("pe_wU50_%d",n).Data(),&pe_wU50[n]);
276  gTREE->SetBranchAddress(TString::Format("pe_wL90_%d",n).Data(),&pe_wL90[n]);
277  gTREE->SetBranchAddress(TString::Format("pe_wU90_%d",n).Data(),&pe_wU90[n]);
278  gTREE->SetBranchAddress(TString::Format("pe_wAVR_%d",n).Data(),&pe_wAVR[n]);
279  gTREE->SetBranchAddress(TString::Format("pe_wRMS_%d",n).Data(),&pe_wRMS[n]);
280  }
281 
282  gTREE->GetEntry(0);
283 
284  cout.precision(4);
285  cout << "---------------------> size " << pe_wMED[0]->size() << endl;
286  cout << "---------------------> pe_trials " << pe_trials << endl;
287 
288  int n = ifo;
289 
290  *winj = GetAlignedWaveform(pe_wINJ[n],pe_wREC[n]);
291  *wrec = GetAlignedWaveform(pe_wREC[n],pe_wREC[n]);
292  *wwht = GetAlignedWaveform(pe_wWHT[n],pe_wREC[n]);
293  *wmed = GetAlignedWaveform(pe_wMED[n],pe_wREC[n]);
294  *wl50 = GetAlignedWaveform(pe_wL50[n],pe_wREC[n]);
295  *wu50 = GetAlignedWaveform(pe_wU50[n],pe_wREC[n]);
296  *wl90 = GetAlignedWaveform(pe_wL90[n],pe_wREC[n]);
297  *wu90 = GetAlignedWaveform(pe_wU90[n],pe_wREC[n]);
298 
299  for(int n=0;n<nIFO;n++) {
300  delete pe_wINJ[n];
301  delete pe_wREC[n];
302  delete pe_wWHT[n];
303  delete pe_wMED[n];
304  delete pe_wL50[n];
305  delete pe_wU50[n];
306  delete pe_wL90[n];
307  delete pe_wU90[n];
308  delete pe_wAVR[n];
309  delete pe_wRMS[n];
310  }
311 
312  return ifoname;
313 }
314 
315 
316 double GetTimeBoundaries(wavearray<double> x, double P, double& bT, double& eT) {
317 
318  if(P<0) P=0;
319  if(P>1) P=1;
320 
321  int N = x.size();
322 
323  double E = 0; // signal energy
324  double avr = 0; // average
325  for(int i=0;i<N;i++) {avr+=i*x[i]*x[i]; E+=x[i]*x[i];}
326  int M=int(avr/E); // central index
327 
328  // search range which contains percentage P of the total energy E
329  int jB=0;
330  int jE=N-1;
331  double a,b;
332  double sum = ((M>=0)&&(M<N)) ? x[M]*x[M] : 0.;
333  for(int j=1; j<N; j++) {
334  a = ((M-j>=0)&&(M-j<N)) ? x[M-j] : 0.;
335  b = ((M+j>=0)&&(M+j<N)) ? x[M+j] : 0.;
336  if(a) jB=M-j;
337  if(b) jE=M+j;
338  sum += a*a+b*b;
339  if(sum/E > P) break;
340  }
341 
342  bT = x.start()+jB/x.rate();
343  eT = x.start()+jE/x.rate();
344 
345  return eT-bT;
346 }
347 
351  TString pdir, double P, bool freq, bool showerr) {
352 
353  int size = wrec->size();
354 
355  wavearray<double> time(size);
356  wavearray<double> etime(size); etime=0;
357  for (int i=0; i<size; i++) time[i] = i/wrec->rate();
358 
359  double bT, eT;
360  GetTimeBoundaries(*wref, P, bT, eT);
361  bT-=wref->start();
362  eT-=wref->start();
363 
364  // set to 0 the frequency values outside the time range -> fix the y scale autoscale
365  // info : this procedure modify the frequency input data but it is not relevant
366  if(freq) {
367  for(int i=0;i<wrec->size();i++) {
368  if(time[i]>bT && time[i]<eT) continue;
369  wrec->data[i]=0; wmed->data[i]=0; wl50->data[i]=0; wu50->data[i]=0; wl90->data[i]=0; wu90->data[i]=0;
370  }
371  }
372 
373  TGraphAsymmErrors* egr90 = new TGraphAsymmErrors(size,time.data,wmed->data,etime.data,etime.data,wl90->data,wu90->data);
374  egr90->SetLineColor(17);
375  egr90->SetFillStyle(1001);
376  egr90->SetFillColor(17);
377  egr90->GetXaxis()->SetTitle("time (s)");
378  if(freq) egr90->GetYaxis()->SetTitle("frequency (hz)");
379  else egr90->GetYaxis()->SetTitle("magnitude");
380  egr90->SetTitle(title);
381 
382  TGraphAsymmErrors* egr50 = new TGraphAsymmErrors(size,time.data,wmed->data,etime.data,etime.data,wl50->data,wu50->data);
383  egr50->SetLineColor(15);
384  egr50->SetFillStyle(1001);
385  egr50->SetFillColor(15);
386 
387  TGraph* agr = new TGraph(size,time.data,wmed->data);
388  agr->SetLineWidth(1);
389  agr->SetLineColor(kWhite);
390  agr->SetLineStyle(1);
391 
392  TGraph* gr = new TGraph(size,time.data,wrec->data);
393  gr->SetLineWidth(1);
394  gr->SetLineColor(2);
395 
396  TCanvas* canvas = new TCanvas("wrec", "wrec",200,20,800,500);
397  canvas->cd();
398  canvas->SetGridx();
399  canvas->SetGridy();
400 
401  egr90->GetXaxis()->SetRangeUser(bT, eT);
402  if(showerr) {
403  egr90->Draw("A4");
404  egr50->Draw("4same");
405  agr->Draw("Lsame");
406  gr->Draw("Lsame");
407  } else {
408  agr->GetXaxis()->SetTitle("time (s)");
409  if(freq) agr->GetYaxis()->SetTitle("frequency (hz)");
410  else agr->GetYaxis()->SetTitle("magnitude");
411  agr->SetTitle(title);
412  agr->GetXaxis()->SetRangeUser(bT, eT);
413  agr->SetLineColor(16);
414  gr->SetLineWidth(2);
415  agr->Draw("AL");
416  gr->Draw("Lsame");
417  }
418 
419  if(ofname!="") {
420  ofname = TString(pdir)+TString("/")+ofname;
421  canvas->Print(ofname);
422  cout << "write : " << ofname << endl;
423  }
424 
425 /*
426  delete canvas;
427  delete egr50;
428  delete egr90;
429  delete agr;
430  delete gr;
431 */
432  return;
433 }
434 
436 
437  wavearray<double> wf = *wf2;
438  wf=0;
439 
440  if(wf1==NULL) return wf;
441  if(wf1->size()==0) return wf;
442 
443  double R=wf1->rate();
444 
445  double b_wf1 = wf1->start();
446  double e_wf1 = wf1->start()+wf1->size()/R;
447  double b_wf2 = wf2->start();
448  double e_wf2 = wf2->start()+wf2->size()/R;
449 
450  int o_wf1 = b_wf1>b_wf2 ? 0 : int((b_wf2-b_wf1)*R+0.5);
451  int o_wf2 = b_wf1<b_wf2 ? 0 : int((b_wf1-b_wf2)*R+0.5);
452 
453  double startXCOR = b_wf1>b_wf2 ? b_wf1 : b_wf2;
454  double endXCOR = e_wf1<e_wf2 ? e_wf1 : e_wf2;
455  int sizeXCOR = int((endXCOR-startXCOR)*R+0.5);
456 
457  for(int i=0;i<sizeXCOR;i++) wf[i+o_wf2] = wf1->data[i+o_wf1];
458 
459  return wf;
460 }
461 
463 
464  double R=wf1->rate();
465 
466  double b_wf1 = wf1->start();
467  double e_wf1 = wf1->start()+wf1->size()/R;
468  double b_wf2 = wf2->start();
469  double e_wf2 = wf2->start()+wf2->size()/R;
470 
471  int o_wf1 = b_wf1>b_wf2 ? 0 : int((b_wf2-b_wf1)*R+0.5);
472  int o_wf2 = b_wf1<b_wf2 ? 0 : int((b_wf1-b_wf2)*R+0.5);
473 
474  double startXCOR = b_wf1>b_wf2 ? b_wf1 : b_wf2;
475  double endXCOR = e_wf1<e_wf2 ? e_wf1 : e_wf2;
476  int sizeXCOR = int((endXCOR-startXCOR)*R+0.5);
477 
478  wavearray<double> wfdif(sizeXCOR);
479  wfdif=0.;
480  wfdif.rate(R);
481  wfdif.start(b_wf1+double(o_wf1)/R);
482 
483  for(int i=0;i<sizeXCOR;i++) wfdif[i] = wf1->data[i+o_wf1] - wf2->data[i+o_wf2];
484 
485  return wfdif;
486 }
487 
490  wavearray<double>* wl90, wavearray<double>* wu90, wavearray<double>* wref, TString pdir, double P, bool freq) {
491 
492  int size = wrec->size();
493 
494  wavearray<double> time(size);
495  wavearray<double> etime(size); etime=0;
496  for (int i=0; i<size; i++) time[i] = i/wrec->rate();
497 
498  double bT, eT;
499  GetTimeBoundaries(*wref, P, bT, eT);
500  bT-=wref->start();
501  eT-=wref->start();
502 
503  // set to 0 the frequency values outside the time range -> fix the y scale autoscale
504  // info : this procedure modify the frequency input data but it is not relevant
505  if(freq) {
506  for(int i=0;i<wrec->size();i++) {
507  if(time[i]>bT && time[i]<eT) continue;
508  wrec->data[i]=0; wmed->data[i]=0; wl50->data[i]=0; wu50->data[i]=0; wl90->data[i]=0; wu90->data[i]=0;
509  }
510  }
511 
512  TGraphAsymmErrors* egr90 = new TGraphAsymmErrors(size,time.data,wmed->data,etime.data,etime.data,wl90->data,wu90->data);
513  egr90->SetLineColor(17);
514  egr90->SetFillStyle(1001);
515  egr90->SetFillColor(17);
516  egr90->GetXaxis()->SetTitle("time (s)");
517  if(freq) egr90->GetYaxis()->SetTitle("frequency (hz)");
518  else egr90->GetYaxis()->SetTitle("magnitude");
519  egr90->SetTitle(title);
520 
521  TGraphAsymmErrors* egr50 = new TGraphAsymmErrors(size,time.data,wmed->data,etime.data,etime.data,wl50->data,wu50->data);
522  egr50->SetLineColor(15);
523  egr50->SetFillStyle(1001);
524  egr50->SetFillColor(15);
525 
526  TGraph* agr = new TGraph(size,time.data,wmed->data);
527  agr->SetLineWidth(1);
528  agr->SetLineColor(kWhite);
529  agr->SetLineStyle(1);
530 
531  TGraph* gr = new TGraph(size,time.data,wrec->data);
532  gr->SetLineWidth(1);
533  gr->SetLineColor(2);
534 
535  TCanvas* canvas = new TCanvas("wrec", "wrec",200,20,800,500);
536  canvas->cd();
537  canvas->SetGridx();
538  canvas->SetGridy();
539 
540  egr90->GetXaxis()->SetRangeUser(bT, eT);
541  egr90->Draw("A4");
542  egr50->Draw("4same");
543  agr->Draw("Lsame");
544  gr->Draw("Lsame");
545 
546  ofname = TString(pdir)+TString("/")+ofname;
547  if(ofname!="") {
548  canvas->Print(ofname);
549  cout << "write : " << ofname << endl;
550  }
551 
552  delete canvas;
553  delete egr50;
554  delete egr90;
555  delete agr;
556  delete gr;
557 }
558 
559 // Dumps reconstructed waveform/time/errors array in ASCII format.
562  wavearray<double>* wl90, wavearray<double>* wu90) {
563 
564  char ofname[256];
565 
566  if(ifo==0) sprintf(ifname,OUT_CWB_ASCII_FILE_L1);
567  if(ifo==1) sprintf(ifname,OUT_CWB_ASCII_FILE_H1);
568 
569  ofstream out;
570  out.open(ofname,ios::out);
571  if (!out.good()) {cout << "Error Opening Output File : " << ofname << endl;exit(1);}
572  cout << "Create Output File : " << ofname << endl;
573  out.precision(19);
574 
575  // write header
576  out << "#whitened data : time, amp_point, amp_median, amp_lower_50_perc, amp_lower_90_perc, amp_upper_50_perc, amp_upper_90_perc" << endl;
577 
578  // write data
579  int size = wrec->size();
580  double dt=1./wrec->rate();
581  for (int i=0; i<size; i++) {
582  double time = i*dt+wrec->start();
583 
584  double vrec = wrec->data[i];
585  double vmed = wmed->data[i];
586  double vl50 = fabs(wl50->data[i]);
587  double vu50 = fabs(wu50->data[i]);
588  double vl90 = fabs(wl90->data[i]);
589  double vu90 = fabs(wu90->data[i]);
590 
591  double l50 = vmed-vl50;
592  double u50 = vmed+vu50;
593  double l90 = vmed-vl90;
594  double u90 = vmed+vu90;
595 
596  // dump full percentiles
597  out << time << " " << vrec << " " << vmed << " " << l50 << " " << l90 << " " << u50 << " " << u90 << endl;
598  // dump dif percentiles from median
599  //out << time << " " << vrec << " " << vmed << " " << vl50 << " " << vl90 << " " << vu50 << " " << vu90 << endl;
600 
601  }
602 
603  out.close();
604 }
605 
606 // Read reconstructed waveform/time/errors array from ascii file
610  wavearray<double>* frec, wavearray<double>* fmed,
612  wavearray<double>* fl90, wavearray<double>* fu90) {
613 
614 // #whitened data : time, amp_point, amp_mean, amp_rms, amp_median, amp_lower_50_perc, amp_lower_90_perc, amp_upper_50_perc, amp_upper_90_perc,
615 // frq_point, frq_mean, frq_rms, frq_median, frq_lower_50_perc, frq_lower_90_perc, frq_upper_50_perc, frq_upper_90_perc
616 
617 
618  char ifname[256];
619 
620  if(ifo==0) sprintf(ifname,"%s/%s",ipath.Data(),IN_CWB_ASCII_FILE_L1);
621  if(ifo==1) sprintf(ifname,"%s/%s",ipath.Data(),IN_CWB_ASCII_FILE_H1);
622 
623  ifstream in(ifname);
624  if(!in.good()) {cout << "Error Opening File : " << ifname << endl;exit(1);}
625 
626  int size=0;
627  char str[1024];
628  int fpos=0;
629  in.getline(str,1024); // skip first line : header
630  while(true) {
631  in.getline(str,1024);
632  if (!in.good()) break;
633  size++;
634  }
635  cout << "size " << size << endl;
636  in.clear(ios::goodbit);
637  in.seekg(0, ios::beg);
638  if (size==0) {cout << "Error : File " << ifile<< " is empty" << endl;gSystem->Exit(1);}
639 
640  double time;
641  double vrec, vavr, vrms, vmed, vl50, vl90, vu50, vu90;
642  double grec, gavr, grms, gmed, gl50, gl90, gu50, gu90;
643 
644  wrec->rate(2048.);
645  wmed->rate(2048.);
646  wl50->rate(2048.);
647  wu50->rate(2048.);
648  wl90->rate(2048.);
649  wu90->rate(2048.);
650 
651  wrec->resize(size);
652  wmed->resize(size);
653  wl50->resize(size);
654  wu50->resize(size);
655  wl90->resize(size);
656  wu90->resize(size);
657 
658  frec->rate(2048.);
659  fmed->rate(2048.);
660  fl50->rate(2048.);
661  fu50->rate(2048.);
662  fl90->rate(2048.);
663  fu90->rate(2048.);
664 
665  frec->resize(size);
666  fmed->resize(size);
667  fl50->resize(size);
668  fu50->resize(size);
669  fl90->resize(size);
670  fu90->resize(size);
671 
672  // skip header
673  in.getline(str,1024);
674 
675  int i=0;
676  while(1) {
677 
678  in >> time
679  >> vrec >> vavr >> vrms >> vmed >> vl50 >> vl90 >> vu50 >> vu90
680  >> grec >> gavr >> grms >> gmed >> gl50 >> gl90 >> gu50 >> gu90;
681 
682  if (!in.good()) break;
683 
684  wrec->data[i] = vrec;
685  wmed->data[i] = vmed;
686 
687  // read full percentile -> reduce to relative percentiles
688  wl50->data[i] = vmed-vl50;
689  wu50->data[i] = vu50-vmed;
690  wl90->data[i] = vmed-vl90;
691  wu90->data[i] = vu90-vmed;
692 
693  if(i==0) {
694  wrec->start(time);
695  wmed->start(time);
696  wl50->start(time);
697  wu50->start(time);
698  wl90->start(time);
699  wu90->start(time);
700  }
701 
702  frec->data[i] = grec;
703  fmed->data[i] = gmed;
704 
705  // read full percentile -> reduce to relative percentiles
706  fl50->data[i] = gmed-gl50;
707  fu50->data[i] = gu50-gmed;
708  fl90->data[i] = gmed-gl90;
709  fu90->data[i] = gu90-gmed;
710 
711  if(i==0) {
712  frec->start(time);
713  fmed->start(time);
714  fl50->start(time);
715  fu50->start(time);
716  fl90->start(time);
717  fu90->start(time);
718  }
719 
720  i++;
721  }
722 
723  wrec->resize(i);
724  wmed->resize(i);
725  wl50->resize(i);
726  wu50->resize(i);
727  wl90->resize(i);
728  wu90->resize(i);
729 
730  frec->resize(i);
731  fmed->resize(i);
732  fl50->resize(i);
733  fu50->resize(i);
734  fl90->resize(i);
735  fu90->resize(i);
736 
737 }
738 
739 void FrequencyCut(wavearray<double>* x, double bF, double eF) {
740 
741  // cut frequency range bF,eF
742  double F=0.;
743  double dF=(x->rate()/(double)x->size())/2.;
744  x->FFTW(1);
745  for(int j=0;j<x->size()/2;j+=2) {
746  F = j*dF;
747  if(F<bF || F>eF) {x->data[j]=0;x->data[j+1]=0;}
748  }
749  x->FFTW(-1);
750 }
751 
752 void PlotSpectrogram(wavearray<double>* x, double tstart, double tstop, TString title, TString ofname, double tchirp, double mchirp) {
753 
754  if(tstart==0) tstart=x->start();
755  if(tstop==0) tstop=x->stop();
756 
757  // START CHIRP FUNCTION
758  if(mchirp>0) {
759  double Mc = mchirp;
760  double Tc = tchirp;
761 
762  const double G = watconstants::GravitationalConstant();
763  const double SM = watconstants::SolarMass();
764  const double C = watconstants::SpeedOfLightInVacuo();
765  const double Pi = TMath::Pi();
766 
767  double p0 = 256.*Pi/5*pow(G*Mc*SM*Pi/C/C/C, 5./3);
768  double p1 = Tc;
769 
770  double tmin=Tc-0.9;
771  double tmax=Tc;
772 
773  // chirp function -> freq = p0 * (time-p1);
774  fchirp = new TF1("fchirp", "pow([0]*([1]-x),-3./8.)", tmin, tmax);
775  fchirp->SetLineColor(kWhite);
776  fchirp->SetLineWidth(1);
777  fchirp->SetLineStyle(2);
778  fchirp->SetParameter(0, p0);
779  fchirp->SetParameter(1, p1);
780  }
781  // END CHIRP FUNCTION
782 
783  TString xtitle = TString::Format("Time (sec) : GPS OFFSET = %.3f",x->start());
784 
785  int nfact=4;
786  int nfft=nfact*512;
787  int noverlap=nfft-10;
788  double fparm=nfact*6;
789  int ystart = int((tstart-x->start()-1)*x->rate());
790  int ystop = int((tstop-x->start()+1)*x->rate());
791  ystart-=nfft;
792  ystop+=nfft;
793  int ysize=ystop-ystart;
794  wavearray<double> y;y.resize(ysize);y.rate(x->rate());y.start(ystart/x->rate());
795 
796  // stft use dt=y.rate() to normalize data but whitened data are already normalized by dt
797  // so before stft data must be divided by 1./sqrt(dt)
798  for(int i=0;i<(int)y.size();i++) y.data[i]=x->data[i+ystart]/sqrt(y.rate());
799 
800  stft = new CWB::STFT(y,nfft,noverlap,"energy","gauss",fparm);
801 
802  TCanvas* canvas;
803  double tstart = nfft/x->rate()+ystart/x->rate();
804  double tstop = (ysize-nfft)/x->rate()+ystart/x->rate();
805 
806  tstart+=0.9;tstop-=0.9;
807  stft->Draw(tstart,tstop,16,1024,0,0,1);
808  if(fchirp!=NULL) fchirp->Draw("SAME");
809  stft->GetHistogram()->SetTitle(title);
810  stft->GetHistogram()->GetXaxis()->SetTitle(xtitle);
811  canvas = stft->GetCanvas();
812  canvas->SetLogy(true);
813  stft->GetHistogram()->GetXaxis()->SetTitle(xtitle);
814  if(ofname!="") stft->Print(ofname);
815 
816  y.resize(0);
817 }
818 
CWB::config * cfg
Definition: TestCWB_Plugin.C:5
detectorParams getDetectorParams()
Definition: detector.cc:201
virtual size_t size() const
Definition: wavearray.hh:127
static const double C
Definition: GNGen.cc:10
int noverlap
Definition: TestDelta.C:20
char xtitle[1024]
void FrequencyCut(wavearray< double > *x, double bF, double eF)
virtual void rate(double r)
Definition: wavearray.hh:123
void Print(TString pname)
Definition: STFT.cc:297
#define OUT_CWB_ASCII_FILE_H1
std::vector< TString > ReadDataFromROOT(TString ifile, int ifo, wavearray< double > *winj, wavearray< double > *wrec, wavearray< double > *wwht, wavearray< double > *wmed, wavearray< double > *wl50, wavearray< double > *wu50, wavearray< double > *wl90, wavearray< double > *wu90)
wavearray< double > a(hp.size())
int n
Definition: cwb_net.C:10
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
double SolarMass()
Definition: constants.hh:184
int nfact
Definition: TestDelta.C:18
#define FHIGH
int nfft
Definition: TestDelta.C:19
return wmap canvas
TRandom3 P
Definition: compare_bkg.C:296
float mchirp
Definition: cbc_plots.C:436
wavearray< double > GetAlignedWaveform(wavearray< double > *wf, wavearray< double > *wref)
waveform wf
Long_t size
watplot p1(const_cast< char * >("TFMap1"))
#define M
Definition: UniqSLagsList.C:3
void PlotSpectrogram(wavearray< double > *x, double tstart=0, double tstop=0, TString title="", TString ofname="", double tchirp=0., double mchirp=0.)
virtual void start(double s)
Definition: wavearray.hh:119
char ofile[512]
int j
Definition: cwb_net.C:10
i drho i
TList * list
#define N
CWB::STFT * stft
char ifo[NIFO_MAX][8]
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
TTree * gTREE
#define nIFO
ofstream out
Definition: cwb_merge.C:196
double tstart
#define FLOW
char str[1024]
double time[6]
Definition: cbc_plots.C:435
double G
Definition: DrawEBHH.C:12
TGraph * gr
static const double SM
Definition: GNGen.cc:8
i() int(T_cor *100))
#define IN_CWB_ASCII_FILE_L1
TString label
Definition: MergeTrees.C:21
double Pi
TF1 * fchirp
const int NIFO_MAX
Definition: wat.hh:4
TH2D * GetHistogram()
Definition: STFT.hh:53
void PlotWaveformAsymmErrors(TString ofname, TString title, wavearray< double > *wrec, wavearray< double > *wmed, wavearray< double > *wl50, wavearray< double > *wu50, wavearray< double > *wl90, wavearray< double > *wu90, wavearray< double > *wref, TString pdir, double P, bool freq=false, bool showerr=true)
wavearray< double > GetDifWaveform(wavearray< double > *wf1, wavearray< double > *wf2)
double tstop
void ReadDataFromASCII(TString ipath, int ifo, wavearray< double > *wrec, wavearray< double > *wmed, wavearray< double > *wl50, wavearray< double > *wu50, wavearray< double > *wl90, wavearray< double > *wu90, wavearray< double > *frec, wavearray< double > *fmed, wavearray< double > *fl50, wavearray< double > *fu50, wavearray< double > *fl90, wavearray< double > *fu90)
void DumpWaveform(int ifo, wavearray< double > *wrec, wavearray< double > *wmed, wavearray< double > *wl50, wavearray< double > *wu50, wavearray< double > *wl90, wavearray< double > *wu90)
double F
TFile * ifile
#define IN_CWB_ASCII_FILE_H1
TFile * froot
double dt
virtual void FFTW(int=1)
Definition: wavearray.cc:878
double GravitationalConstant()
Definition: constants.hh:113
static void TimeShift(wavearray< double > &x, double tShift=0.)
Definition: mdc.cc:2874
double GetTimeBoundaries(wavearray< double > x, double P, double &bT, double &eT)
int nifo
char title[256]
Definition: SSeriesExample.C:1
ifstream in
char Name[16]
Definition: detector.hh:309
virtual void stop(double s)
Definition: wavearray.hh:121
double fabs(const Complex &x)
Definition: numpy.cc:37
TCanvas * GetCanvas()
Definition: STFT.hh:52
#define OUT_CWB_ASCII_FILE_L1
void DrawWaveformPE(TString ipath, int gtype=0, int ifo=0, double tshift=0, TString label="PLOT", double P=0.99)
sprintf(tfres,"(1/%g)x(%g) (sec)x(Hz)", 2 *df, df)
DataType_t * data
Definition: wavearray.hh:301
virtual void resize(unsigned int)
Definition: wavearray.cc:445
wavearray< double > y
Definition: Test10.C:31
gwavearray< double > * grec
double fparm
Definition: TestDelta.C:22
exit(0)
double SpeedOfLightInVacuo()
Definition: constants.hh:96
double avr