Logo coherent WaveBurst  
Library Reference Guide
Logo
CBCTool.hh
Go to the documentation of this file.
1 /*
2 # Copyright (C) 2019 Francesco Salemi
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  * Package: CBCTool Class Library
20  * File name: CBCTool.hh
21  * Author: Francesco Salemi (francesco.salemi@aei.mpg.de)
22  **********************************************************/
23 
24 #ifndef CBCTOOL_HH
25 #define CBCTOOL_HH
26 
27 #include "Math/BrentRootFinder.h"
28 #include "Math/WrappedTF1.h"
29 #include "TApplication.h"
30 #include "TCanvas.h"
31 #include "TChain.h"
32 #include "TF1.h"
33 #include "TFile.h"
34 #include "TGaxis.h"
35 #include "TGraphSmooth.h"
36 #include "TH1D.h"
37 #include "TH1I.h"
38 #include "TH2F.h"
39 #include "TLegend.h"
40 #include "TMath.h"
41 #include "TMultiGraph.h"
42 #include "TObjArray.h"
43 #include "TObjString.h"
44 #include "TROOT.h"
45 #include "TRandom.h"
46 #include "TRandom3.h"
47 #include "TRatioPlot.h"
48 #include "TString.h"
49 #include "TStyle.h"
50 #include "TSystem.h"
51 #include "TSystemDirectory.h"
52 #include "TTree.h"
53 #include "TTreeIndex.h"
54 
55 #include "TGraph2DErrors.h"
56 #include "TPaveText.h"
57 #include "TText.h"
58 
59 #include <ctype.h>
60 #include <math.h>
61 #include <stdlib.h>
62 #include <string.h>
63 #include <fstream>
64 #include <iostream>
65 #include <string>
66 #include <vector>
67 
68 #include "Meyer.hh"
69 #include "netevent.hh"
70 #include "network.hh"
71 #include "wavearray.hh"
72 #include "wavecomplex.hh"
73 
74 #include "FrameL.h"
75 
76 #include "History.hh"
77 #include "Toolbox.hh"
78 #include "Toolfun.hh"
79 
80 #define LST_TREE_NAME "frl"
81 
82 using namespace std;
83 
84 namespace CWB {
85 
86 class CBCTool {
87  public:
88  /* ******************************** */
89  /* * post processing methods * */
90  /* ******************************** */
91 
92  static Long64_t GetTreeIndex(TTree* tree, const char* selection);
93 
94  static void doROCPlot(int bkg_entries,
95  double* rho_bkg,
96  int* index,
97  float RHO_BIN,
98  double liveTot,
99  TF1* AverageRad,
100  TCanvas* c1,
101  TString odir,
102  bool write_ascii);
103 
104  static void doROCPlot(int bkg_entries,
105  double* rho_bkg,
106  int* index,
107  float RHO_BIN,
108  float RHO_NBINS,
109  float RHO_MIN,
110  double liveTot,
111  double* Rrho,
112  double* eRrho,
113  double* Trho,
114  TCanvas* c1,
115  TString odir,
116  bool write_ascii);
117 
118  static TF1* doRangePlot(int RHO_NBINS,
119  double* Trho,
120  double* Rrho,
121  double* eRrho,
122  float RHO_MIN,
123  float T_out,
124  TCanvas* c1,
125  TString networkname,
126  TString odir,
127  bool write_ascii);
128 
129  static double getLiveTime(int nIFO,
130  TChain& liv,
131  int lag_number,
132  int slag_number,
133  int dummy);
134 
135  static double getLiveTime2(TChain& liv);
136 
137  static double getZeroLiveTime(int nIFO, TChain& liv);
138 
139  static double getZeroLiveTime2(int nIFO, TChain& liv);
140 
141  static double getFAR(float rho, TH1* hc, double liveTot);
142 
143  static TH2F* FillSLagHist(int NIFO_MAX,
144  TChain& live,
145  int NSlag,
146  double SlagMin,
147  double SlagMax);
148 
149  static void doChirpFARPlot(int sel_events,
150  double* recMchirp,
151  double* injMchirp,
152  double* far,
153  TCanvas* c1,
154  TString odir);
155 
156  static double calc_isco_radius(double a);
157 
158  static double calc_isco_freq(double a);
159 
160  static double _final_spin_diff(double a_f,
161  double eta,
162  double delta_m,
163  double S,
164  double Delta);
165 
166  static double _final_spin_diff(Double_t* x, Double_t* par);
167 
168  static double bbh_final_mass_and_spin_non_precessing(double m1,
169  double m2,
170  double chi1,
171  double chi2);
172 
173  static double chip(double m1,
174  double m2,
175  double s1x,
176  double s1y,
177  double s1z,
178  double s2x,
179  double s2y,
180  double s2z);
181 
182  static void AddChip(TString filein, TString treename);
183 
184  static void CreateDistanceParplots(char* sim_file_name,
185  char* mdc_file_name,
186  char* netdir,
187  TString opt,
188  double MINX,
189  double MAXX,
190  double MAXDISTANCE,
191  int NBIN_DIST,
192  float T_ifar,
193  float T_win,
194  int nIFO);
195 
196  static TGraphErrors* CreateGraphRadiusIFAR(char* sim_file_name,
197  char* mdc_file_name,
198  TString SEL,
199  float shell_volume,
200  Color_t color,
201  TString opt,
202  double liveTot,
203  float T_ifar,
204  float T_win,
205  int TRIALS,
206  int nIFO,
207  float VT,
208  float Tscale);
209 
210  private:
211  ClassDef(CBCTool, 1)
212 };
213 
214 } // namespace CWB
215 
216 #endif
TTree * tree
Definition: TimeSortTree.C:20
double rho
float Tscale
double T_ifar
TString live
Definition: ced.hh:42
TString opt
double m1
wavearray< double > a(hp.size())
void AddChip(TString filein, TString treename)
Definition: AddChip.C:19
#define RHO_NBINS
TString("c")
#define RHO_BIN
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
char odir[1024]
STL namespace.
void CreateDistanceParplots(char *sim_file_name, char *mdc_file_name, char *netdir, TString opt="", double MINX=0.0, double MAXX=1.0, double MAXDISTANCE=5000., int NBIN_DIST=10, float T_ifar=0.0, float T_win=0.2, int nIFO=2)
TGraphErrors * CreateGraphRadiusIFAR(char *sim_file_name, char *mdc_file_name, TString SEL, float shell_volume, Color_t color=kBlue, TString opt="default", double liveTot=1e6, float T_ifar=0.0, float T_win=0.2, int TRIALS=1, int nIFO=2)
#define nIFO
TCanvas * c1
char netdir[1024]
const int NIFO_MAX
Definition: wat.hh:22
double liveTot
vector< mdcpar > par
double VT
#define RHO_MIN
char sim_file_name[1024]
wavearray< int > index
double T_win
double m2
Meyer< double > S(1024, 2)
double T_out
bool write_ascii
float * shell_volume
char mdc_file_name[1024]
static double Delta
Definition: geodesics.cc:26
int TRIALS