Logo coherent WaveBurst  
Library Reference Guide
Logo
DrawCosOmegaPRC.C
Go to the documentation of this file.
1 /*
2 # Copyright (C) 2019 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 // Draw PRC Cos Omega distribution
21 // Note : this macro is used to generate the PRC report (cwb_report merge_label prc)
22 // Author : Gabriele Vedovato
23 
24 #include "TROOT.h"
25 #include "TSystem.h"
26 #include "TFile.h"
27 #include "TTree.h"
28 #include <fstream>
29 #include <iostream>
30 #include "TGraph.h"
31 #include "TMultiGraph.h"
32 #include "TCanvas.h"
33 #include "TH1F.h"
34 #include "TMath.h"
35 #include "TH1F.h"
36 #include "TH2F.h"
37 #include <TComplex.h>
38 #include <TStyle.h>
39 #include <TRandom.h>
40 #include "TVector3.h"
41 #include "TRotation.h"
42 #include "Math/Vector3Dfwd.h"
43 #include "Math/Rotation3D.h"
44 
45 TH1F* co_hist;
46 TCanvas* co_canvas;
47 
48 void
50  int pp_inetcc, float T_cor, int pp_irho, float T_cut, float T_vED, float T_pen, float T_ifar) {
51 
52  if (!gROOT->GetClass("Polar3DVector")) gSystem->Load("libMathCore");
53 
54  using namespace ROOT::Math;
55 
56  char wave_file_name[1024];
57  sprintf(wave_file_name,"merge/wave_%s.%s.root",data_label.Data(),merge_label.Data());
58 
59  TFile* ifile = TFile::Open(wave_file_name);
60  TTree* itree = (TTree *) gROOT->FindObject("waveburst");
61 
62  char cut[1024];
63  char tmp[1024];
64  sprintf(cut,"abs(time[0]-time[%d])<%f && netcc[%d]>%f && rho[%d]>%f",
65  nIFO,T_win,pp_inetcc,T_cor,pp_irho,T_cut);
66  if(T_vED>0) {strcpy(tmp,cut);sprintf(cut,"%s && neted[0]/ecor<%f",tmp,T_vED);}
67  if(T_pen>0) {strcpy(tmp,cut);sprintf(cut,"%s && penalty>%f",tmp,T_pen);}
68  if(T_ifar>0) {strcpy(tmp,cut);sprintf(cut,"%s && ifar>(24.*3600.*365.)*%f",tmp,T_ifar);}
69 
70  itree->Draw("theta[0]:phi[0]:theta[1]:phi[1]",cut,"goff");
71  int size = itree->GetSelectedRows();
72  cout << "detected events : " << size << endl;
73 
74  double* th0 = itree->GetV1();
75  double* ph0 = itree->GetV2();
76  double* th1 = itree->GetV3();
77  double* ph1 = itree->GetV4();
78 
79  co_hist = new TH1F("co_hist","co_hist",24,-1,1);
80  //co_hist = new TH1F("co_hist","co_hist",100,0,180);
81 
82  double deg2rad = TMath::Pi()/180.;
83 
84  for(int i=0;i<size;i++) {
85 
86  //cout << i << " " << th0[i] << " " << ph0[i] << endl;
87  //cout << i << " " << th1[i] << " " << ph1[i] << endl;
88 
89  Polar3DVector v0(1, th0[i]*deg2rad, ph0[i]*deg2rad);
90  Polar3DVector v1(1, th1[i]*deg2rad, ph1[i]*deg2rad);
91 
92  double Dot = v0.Dot(v1);
93  double dOmega = 180.*TMath::ACos(Dot)/TMath::Pi();
94 
95 /*
96  double inj_phi = ph1[i]*deg2rad;
97  double est_phi = ph0[i]*deg2rad;
98  double inj_theta = th1[i]*deg2rad;
99  double est_theta = th0[i]*deg2rad;
100  cos_delta_theta = cos(inj_theta)*cos(est_theta) + sin(inj_theta)*sin(est_theta)*cos(inj_phi-est_phi);
101  //co_hist->Fill(cos_delta_theta);
102 */
103 
104  //cout << i << "Dot : " << Dot << " dOmega : " << dOmega << endl;
105  co_hist->Fill(Dot);
106  //co_hist->Fill(dOmega);
107  }
108 
109  // Normalization
110  double integral = 0;
111  for (int i=0;i<=co_hist->GetNbinsX();i++) integral+=co_hist->GetBinContent(i);
112  for (int i=0;i<=co_hist->GetNbinsX();i++) co_hist->SetBinContent(i,co_hist->GetBinContent(i)/integral);
113 
114  char title[256];
115  sprintf(title,"angle offset : #events = %d",size);
116  co_hist->SetTitle(title);
117  co_hist->GetXaxis()->SetTitle("cos(angle offset)");
118  co_hist->GetYaxis()->SetTitle("probability density");
119  co_hist->GetXaxis()->SetRangeUser(-1,1);
120  co_hist->GetYaxis()->SetRangeUser(0.001,1);
121 // co_hist->GetYaxis()->SetRangeUser(1e-3,1);
122  co_hist->GetXaxis()->SetTitleOffset(1.35);
123  co_hist->GetYaxis()->SetTitleOffset(1.35);
124  co_hist->GetXaxis()->CenterTitle(true);
125  co_hist->GetYaxis()->CenterTitle(true);
126  co_hist->GetXaxis()->SetLabelFont(42);
127  co_hist->GetXaxis()->SetTitleFont(42);
128  co_hist->GetYaxis()->SetLabelFont(42);
129  co_hist->GetYaxis()->SetTitleFont(42);
130  co_hist->SetLineColor(kRed);
131  co_hist->SetFillColor(kRed);
132 
133  co_hist->SetStats(kFALSE);
134  co_canvas = new TCanvas("fom", "PRC", 300,40, 600, 500);
135  co_canvas->SetGridx();
136  co_canvas->SetGridy();
137  co_canvas->SetLogy();
138  co_hist->Draw("HIST");
139 
140  char ofname[1024];
141  sprintf(ofname,"%s/cos_theta.gif",odir.Data());
142  co_canvas->Print(ofname);
143  char cmd[1024];
144  TString pfname(ofname);
145  pfname.ReplaceAll(".gif",".png");
146  sprintf(cmd,"convert %s %s",ofname,pfname.Data());
147  cout << cmd << endl;
148  gSystem->Exec(cmd);
149  sprintf(cmd,"rm %s",ofname);
150  cout << cmd << endl;
151  gSystem->Exec(cmd);
152  //exit(0);
153 
154 }
155 
double T_ifar
void DrawCosOmegaPRC(TString data_label, TString odir, TString merge_label, int nIFO, float T_win, int pp_inetcc, float T_cor, int pp_irho, float T_cut, float T_vED, float T_pen, float T_ifar)
double T_pen
TString("c")
double T_cor
char odir[1024]
i pp_inetcc
Long_t size
i drho i
#define nIFO
char data_label[512]
Definition: test_config1.C:160
double deg2rad
double Pi
double * tmp
Definition: testWDM_5.C:31
char cut[512]
TString merge_label
TFile * ifile
char title[256]
Definition: SSeriesExample.C:1
TH1F * co_hist
double T_win
char cmd[1024]
strcpy(RunLabel, RUN_LABEL)
sprintf(tfres,"(1/%g)x(%g) (sec)x(Hz)", 2 *df, df)
double T_cut
TTree * itree
double T_vED
i drho pp_irho
TCanvas * co_canvas