1 | /*
|
---|
2 | ---- Delphes ----
|
---|
3 | A Fast Simulator for general purpose LHC detector
|
---|
4 | S. Ovyn ~~~~ severine.ovyn@uclouvain.be
|
---|
5 |
|
---|
6 | Center for Particle Physics and Phenomenology (CP3)
|
---|
7 | Universite Catholique de Louvain (UCL)
|
---|
8 | Louvain-la-Neuve, Belgium
|
---|
9 | */
|
---|
10 |
|
---|
11 | /// \file SmearUtil.cc
|
---|
12 | /// \brief RESOLution class, and some generic definitions
|
---|
13 |
|
---|
14 |
|
---|
15 | #include "interface/SmearUtil.h"
|
---|
16 | #include "TRandom.h"
|
---|
17 |
|
---|
18 | #include <iostream>
|
---|
19 | #include <sstream>
|
---|
20 | #include <fstream>
|
---|
21 | #include <iomanip>
|
---|
22 |
|
---|
23 |
|
---|
24 |
|
---|
25 | using namespace std;
|
---|
26 |
|
---|
27 | //------------------------------------------------------------------------------
|
---|
28 |
|
---|
29 | RESOLution::RESOLution() {
|
---|
30 |
|
---|
31 | MAX_TRACKER = 2.5; // tracker coverage
|
---|
32 | MAX_CALO_CEN = 3.0; // central calorimeter coverage
|
---|
33 | MAX_CALO_FWD = 5.0; // forward calorimeter pseudorapidity coverage
|
---|
34 | MAX_MU = 2.4; // muon chambers pseudorapidity coverage
|
---|
35 | MIN_CALO_VFWD= 5.2; // very forward calorimeter (if any), like CASTOR
|
---|
36 | MAX_CALO_VFWD= 6.6; // very forward calorimeter (if any), like CASTOR
|
---|
37 | MIN_ZDC = 8.3; // zero-degree calorimeter, coverage
|
---|
38 |
|
---|
39 | ZDC_S = 140.; // ZDC distance to IP
|
---|
40 | RP220_S = 220; // distance of the RP to the IP, in meters
|
---|
41 | RP220_X = 0.002;// distance of the RP to the beam, in meters
|
---|
42 | FP420_S = 420; // distance of the RP to the IP, in meters
|
---|
43 | FP420_X = 0.004;// distance of the RP to the beam, in meters
|
---|
44 |
|
---|
45 |
|
---|
46 | ELG_Scen = 0.05; // S term for central ECAL
|
---|
47 | ELG_Ncen = 0.25 ; // N term for central ECAL
|
---|
48 | ELG_Ccen = 0.0055 ; // C term for central ECAL
|
---|
49 | ELG_Cfwd = 0.107 ; // S term for forward ECAL
|
---|
50 | ELG_Sfwd = 2.084 ; // C term for forward ECAL
|
---|
51 | ELG_Nfwd = 0.0 ; // N term for central ECAL
|
---|
52 |
|
---|
53 | HAD_Shcal = 1.5 ; // S term for central HCAL // hadronic calorimeter
|
---|
54 | HAD_Nhcal = 0.0 ; // N term for central HCAL
|
---|
55 | HAD_Chcal = 0.05 ; // C term for central HCAL
|
---|
56 | HAD_Shf = 2.7 ; // S term for central HF // forward calorimeter
|
---|
57 | HAD_Nhf = 0.0 ; // N term for central HF
|
---|
58 | HAD_Chf = 0.13 ; // C term for central HF
|
---|
59 |
|
---|
60 | MU_SmearPt = 0.01 ;
|
---|
61 |
|
---|
62 | ELEC_pt = 10.0;
|
---|
63 | MUON_pt = 10.0;
|
---|
64 | JET_pt = 20.0;
|
---|
65 | TAUJET_pt = 10.0;
|
---|
66 |
|
---|
67 |
|
---|
68 | TAU_CONE_ENERGY = 0.15 ; // Delta R = radius of the cone // for "electromagnetic collimation"
|
---|
69 | TAU_EM_COLLIMATION = 0.95;
|
---|
70 | TAU_CONE_TRACKS= 0.4 ; //Delta R for tracker isolation for tau's
|
---|
71 | PT_TRACK_TAU = 2.0 ; // GeV // 6 GeV ????
|
---|
72 |
|
---|
73 |
|
---|
74 | PT_TRACKS_MIN = 0.9 ; // minimal pt needed to reach the calorimeter, in GeV
|
---|
75 | PT_QUARKS_MIN = 2.0 ; // minimal pt needed by quarks to reach the tracker, in GeV (??????)
|
---|
76 | TRACKING_EFF = 90;
|
---|
77 |
|
---|
78 |
|
---|
79 | TAGGING_B = 40;
|
---|
80 | MISTAGGING_C = 10;
|
---|
81 | MISTAGGING_L = 1;
|
---|
82 |
|
---|
83 |
|
---|
84 | CONERADIUS = 0.7; // generic jet radius ; not for tau's !!!
|
---|
85 | JETALGO = 1; // 1 for Cone algorithm, 2 for MidPoint algorithm, 3 for SIScone algorithm, 4 for kt algorithm
|
---|
86 |
|
---|
87 | //General jet parameters
|
---|
88 | SEEDTHRESHOLD = 1.0;
|
---|
89 | OVERLAPTHRESHOLD = 0.75;
|
---|
90 |
|
---|
91 | // Define Cone algorithm.
|
---|
92 | C_ADJACENCYCUT = 2;
|
---|
93 | C_MAXITERATIONS = 100;
|
---|
94 | C_IRATCH = 1;
|
---|
95 |
|
---|
96 | //Define MidPoint algorithm.
|
---|
97 | M_CONEAREAFRACTION = 0.25;
|
---|
98 | M_MAXPAIRSIZE = 2;
|
---|
99 | M_MAXITERATIONS = 100;
|
---|
100 |
|
---|
101 | }
|
---|
102 |
|
---|
103 | //------------------------------------------------------------------------------
|
---|
104 | void RESOLution::ReadDataCard(const string datacard) {
|
---|
105 |
|
---|
106 | string temp_string;
|
---|
107 | istringstream curstring;
|
---|
108 |
|
---|
109 | ifstream fichier_a_lire(datacard.c_str());
|
---|
110 | if(!fichier_a_lire.good()) {
|
---|
111 | cout << datacard << "Datadard " << datacard << " not found, use default values" << endl;
|
---|
112 | return;
|
---|
113 | }
|
---|
114 |
|
---|
115 | while (getline(fichier_a_lire,temp_string)) {
|
---|
116 | curstring.clear(); // needed when using several times istringstream::str(string)
|
---|
117 | curstring.str(temp_string);
|
---|
118 | string varname;
|
---|
119 | float value;
|
---|
120 |
|
---|
121 | if(strstr(temp_string.c_str(),"#")) { }
|
---|
122 | else if(strstr(temp_string.c_str(),"MAX_TRACKER")){curstring >> varname >> value; MAX_TRACKER = value;}
|
---|
123 | else if(strstr(temp_string.c_str(),"MAX_CALO_CEN")){curstring >> varname >> value; MAX_CALO_CEN = value;}
|
---|
124 | else if(strstr(temp_string.c_str(),"MAX_CALO_FWD")){curstring >> varname >> value; MAX_CALO_FWD = value;}
|
---|
125 | else if(strstr(temp_string.c_str(),"MAX_MU")){curstring >> varname >> value; MAX_MU = value;}
|
---|
126 | else if(strstr(temp_string.c_str(),"ELG_Scen")){curstring >> varname >> value; ELG_Scen = value;}
|
---|
127 | else if(strstr(temp_string.c_str(),"ELG_Ncen")){curstring >> varname >> value; ELG_Ncen = value;}
|
---|
128 | else if(strstr(temp_string.c_str(),"ELG_Ccen")){curstring >> varname >> value; ELG_Ccen = value;}
|
---|
129 | else if(strstr(temp_string.c_str(),"ELG_Sfwd")){curstring >> varname >> value; ELG_Sfwd = value;}
|
---|
130 | else if(strstr(temp_string.c_str(),"ELG_Cfwd")){curstring >> varname >> value; ELG_Cfwd = value;}
|
---|
131 | else if(strstr(temp_string.c_str(),"ELG_Nfwd")){curstring >> varname >> value; ELG_Nfwd = value;}
|
---|
132 | else if(strstr(temp_string.c_str(),"HAD_Shcal")){curstring >> varname >> value; HAD_Shcal = value;}
|
---|
133 | else if(strstr(temp_string.c_str(),"HAD_Nhcal")){curstring >> varname >> value; HAD_Nhcal = value;}
|
---|
134 | else if(strstr(temp_string.c_str(),"HAD_Chcal")){curstring >> varname >> value; HAD_Chcal = value;}
|
---|
135 | else if(strstr(temp_string.c_str(),"HAD_Shf")){curstring >> varname >> value; HAD_Shf = value;}
|
---|
136 | else if(strstr(temp_string.c_str(),"HAD_Nhf")){curstring >> varname >> value; HAD_Nhf = value;}
|
---|
137 | else if(strstr(temp_string.c_str(),"HAD_Chf")){curstring >> varname >> value; HAD_Chf = value;}
|
---|
138 | else if(strstr(temp_string.c_str(),"MU_SmearPt")){curstring >> varname >> value; MU_SmearPt = value;}
|
---|
139 | else if(strstr(temp_string.c_str(),"TAU_CONE_ENERGY")){curstring >> varname >> value; TAU_CONE_ENERGY = value;}
|
---|
140 | else if(strstr(temp_string.c_str(),"TAU_CONE_TRACKS")){curstring >> varname >> value; TAU_CONE_TRACKS = value;}
|
---|
141 | else if(strstr(temp_string.c_str(),"PT_TRACK_TAU")){curstring >> varname >> value; PT_TRACK_TAU = value;}
|
---|
142 | else if(strstr(temp_string.c_str(),"PT_TRACKS_MIN")){curstring >> varname >> value; PT_TRACKS_MIN = value;}
|
---|
143 | else if(strstr(temp_string.c_str(),"TAGGING_B")){curstring >> varname >> value; TAGGING_B = (int)value;}
|
---|
144 | else if(strstr(temp_string.c_str(),"MISTAGGING_C")){curstring >> varname >> value; MISTAGGING_C = (int)value;}
|
---|
145 | else if(strstr(temp_string.c_str(),"MISTAGGING_L")){curstring >> varname >> value; MISTAGGING_L = (int)value;}
|
---|
146 | else if(strstr(temp_string.c_str(),"CONERADIUS")){curstring >> varname >> value; CONERADIUS = value;}
|
---|
147 | else if(strstr(temp_string.c_str(),"JETALGO")){curstring >> varname >> value; JETALGO = (int)value;}
|
---|
148 | else if(strstr(temp_string.c_str(),"TRACKING_EFF")){curstring >> varname >> value; TRACKING_EFF = (int)value;}
|
---|
149 | else if(strstr(temp_string.c_str(),"ELEC_pt")){curstring >> varname >> value; ELEC_pt = value;}
|
---|
150 | else if(strstr(temp_string.c_str(),"MUON_pt")){curstring >> varname >> value; MUON_pt = value;}
|
---|
151 | else if(strstr(temp_string.c_str(),"JET_pt")){curstring >> varname >> value; JET_pt = value;}
|
---|
152 | else if(strstr(temp_string.c_str(),"TAUJET_pt")){curstring >> varname >> value; TAUJET_pt = value;}
|
---|
153 |
|
---|
154 | }
|
---|
155 |
|
---|
156 | // General jet variables
|
---|
157 | SEEDTHRESHOLD = 1.0;
|
---|
158 | OVERLAPTHRESHOLD = 0.75;
|
---|
159 |
|
---|
160 | // Define Cone algorithm.
|
---|
161 | C_ADJACENCYCUT = 2;
|
---|
162 | C_MAXITERATIONS = 100;
|
---|
163 | C_IRATCH = 1;
|
---|
164 |
|
---|
165 | //Define MidPoint algorithm.
|
---|
166 | M_CONEAREAFRACTION = 0.25;
|
---|
167 | M_MAXPAIRSIZE = 2;
|
---|
168 | M_MAXITERATIONS = 100;
|
---|
169 |
|
---|
170 | //Define SISCone algorithm.
|
---|
171 | NPASS = 0;
|
---|
172 | PROTOJET_PTMIN = 0.0;
|
---|
173 |
|
---|
174 |
|
---|
175 | }
|
---|
176 |
|
---|
177 | void RESOLution::Logfile(string LogName) {
|
---|
178 | //void RESOLution::Logfile(string outputfilename) {
|
---|
179 |
|
---|
180 | ofstream f_out(LogName.c_str());
|
---|
181 |
|
---|
182 | f_out<<"#*********************************************************************"<<"\n";
|
---|
183 | f_out<<"# *"<<"\n";
|
---|
184 | f_out<<"# ---- DELPHES release 1.0 ---- *"<<"\n";
|
---|
185 | f_out<<"# *"<<"\n";
|
---|
186 | f_out<<"# A Fast Simulator for general purpose LHC detector *"<<"\n";
|
---|
187 | f_out<<"# Written by S. Ovyn and X. Rouby *"<<"\n";
|
---|
188 | f_out<<"# severine.ovyn@uclouvain.be *"<<"\n";
|
---|
189 | f_out<<"# *"<<"\n";
|
---|
190 | f_out<<"# http: *"<<"\n";
|
---|
191 | f_out<<"# *"<<"\n";
|
---|
192 | f_out<<"# Center for Particle Physics and Phenomenology (CP3) *"<<"\n";
|
---|
193 | f_out<<"# Universite Catholique de Louvain (UCL) *"<<"\n";
|
---|
194 | f_out<<"# Louvain-la-Neuve, Belgium *"<<"\n";
|
---|
195 | f_out<<"# *"<<"\n";
|
---|
196 | f_out<<"#....................................................................*"<<"\n";
|
---|
197 | f_out<<"# *"<<"\n";
|
---|
198 | f_out<<"# This package uses: *"<<"\n";
|
---|
199 | f_out<<"# FastJet algorithm: Phys. Lett. B641 (2006) [hep-ph/0512210] *"<<"\n";
|
---|
200 | f_out<<"# Hector: JINST 2:P09005 (2007) [physics.acc-ph:0707.1198v2] *"<<"\n";
|
---|
201 | f_out<<"# ExRootAnalysis *"<<"\n";
|
---|
202 | f_out<<"# *"<<"\n";
|
---|
203 | f_out<<"#....................................................................*"<<"\n";
|
---|
204 | f_out<<"# *"<<"\n";
|
---|
205 | f_out<<"# This file contains all the running parameters (detector and cuts) *"<<"\n";
|
---|
206 | f_out<<"# necessary to reproduce the detector simulation *"<<"\n";
|
---|
207 | f_out<<"# *"<<"\n";
|
---|
208 | f_out<<"#....................................................................*"<<"\n";
|
---|
209 | f_out<<"#>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>"<<"\n";
|
---|
210 | f_out<<"* *"<<"\n";
|
---|
211 | f_out<<"#******************************** *"<<"\n";
|
---|
212 | f_out<<"# Central detector caracteristics *"<<"\n";
|
---|
213 | f_out<<"#******************************** *"<<"\n";
|
---|
214 | f_out<<"* *"<<"\n";
|
---|
215 | f_out << left << setw(30) <<"* Maximum tracking system: "<<""
|
---|
216 | << left << setw(10) <<MAX_TRACKER <<""<< right << setw(15)<<"*"<<"\n";
|
---|
217 | f_out << left << setw(30) <<"* Maximum central calorimeter: "<<""
|
---|
218 | << left << setw(10) <<MAX_CALO_CEN <<""<< right << setw(15)<<"*"<<"\n";
|
---|
219 | f_out << left << setw(30) <<"* Maximum forward calorimeter: "<<""
|
---|
220 | << left << setw(10) <<MAX_CALO_FWD <<""<< right << setw(15)<<"*"<<"\n";
|
---|
221 | f_out << left << setw(30) <<"* Muon chambers coverage: "<<""
|
---|
222 | << left << setw(10) <<MAX_MU <<""<< right << setw(15)<<"*"<<"\n";
|
---|
223 | f_out<<"* *"<<"\n";
|
---|
224 | f_out<<"#************************************* *"<<"\n";
|
---|
225 | f_out<<"# Very forward detector caracteristics *"<<"\n";
|
---|
226 | f_out<<"#************************************* *"<<"\n";
|
---|
227 | f_out<<"* *"<<"\n";
|
---|
228 | f_out << left << setw(55) <<"* Minimum very forward calorimeter: "<<""
|
---|
229 | << left << setw(5) <<MIN_CALO_VFWD <<""<< right << setw(10)<<"*"<<"\n";
|
---|
230 | f_out << left << setw(55) <<"* Maximum very forward calorimeter: "<<""
|
---|
231 | << left << setw(5) <<MAX_CALO_VFWD <<""<< right << setw(10)<<"*"<<"\n";
|
---|
232 | f_out << left << setw(55) <<"* Distance of the ZDC to the IP, in meters: "<<""
|
---|
233 | << left << setw(5) <<ZDC_S <<""<< right << setw(10)<<"*"<<"\n";
|
---|
234 | f_out << left << setw(55) <<"* Distance of the RP to the IP, in meters: "<<""
|
---|
235 | << left << setw(5) <<RP220_S <<""<< right << setw(10)<<"*"<<"\n";
|
---|
236 | f_out << left << setw(55) <<"* Distance of the RP to the beam, in meters: "<<""
|
---|
237 | << left << setw(5) <<RP220_X <<""<< right << setw(10)<<"*"<<"\n";
|
---|
238 | f_out << left << setw(55) <<"* Distance of the RP to the IP, in meters: "<<""
|
---|
239 | << left << setw(5) <<FP420_S <<""<< right << setw(10)<<"*"<<"\n";
|
---|
240 | f_out << left << setw(55) <<"* Distance of the RP to the beam, in meters: "<<""
|
---|
241 | << left << setw(5) <<FP420_X <<""<< right << setw(10)<<"*"<<"\n";
|
---|
242 | f_out<<"* *"<<"\n";
|
---|
243 | f_out<<"#************************************ *"<<"\n";
|
---|
244 | f_out<<"# Electromagnetic smearing parameters *"<<"\n";
|
---|
245 | f_out<<"#************************************ *"<<"\n";
|
---|
246 | f_out<<"* *"<<"\n";
|
---|
247 | //# \sigma/E = C + N/E + S/\sqrt{E}
|
---|
248 | f_out << left << setw(30) <<"* S term for central ECAL: "<<""
|
---|
249 | << left << setw(30) <<ELG_Scen <<""<< right << setw(10)<<"*"<<"\n";
|
---|
250 | f_out << left << setw(30) <<"* N term for central ECAL: "<<""
|
---|
251 | << left << setw(30) <<ELG_Ncen <<""<< right << setw(10)<<"*"<<"\n";
|
---|
252 | f_out << left << setw(30) <<"* C term for central ECAL: "<<""
|
---|
253 | << left << setw(30) <<ELG_Ccen <<""<< right << setw(10)<<"*"<<"\n";
|
---|
254 | f_out << left << setw(30) <<"* S term for forward ECAL: "<<""
|
---|
255 | << left << setw(30) <<ELG_Sfwd <<""<< right << setw(10)<<"*"<<"\n";
|
---|
256 | f_out << left << setw(30) <<"* N term for forward ECAL: "<<""
|
---|
257 | << left << setw(30) <<ELG_Nfwd <<""<< right << setw(10)<<"*"<<"\n";
|
---|
258 | f_out << left << setw(30) <<"* C term for forward ECAL: "<<""
|
---|
259 | << left << setw(30) <<ELG_Cfwd <<""<< right << setw(10)<<"*"<<"\n";
|
---|
260 | f_out<<"* *"<<"\n";
|
---|
261 | f_out<<"#***************************** *"<<"\n";
|
---|
262 | f_out<<"# Hadronic smearing parameters *"<<"\n";
|
---|
263 | f_out<<"#***************************** *"<<"\n";
|
---|
264 | f_out<<"* *"<<"\n";
|
---|
265 | f_out << left << setw(30) <<"* S term for central HCAL: "<<""
|
---|
266 | << left << setw(30) <<HAD_Shcal <<""<< right << setw(10)<<"*"<<"\n";
|
---|
267 | f_out << left << setw(30) <<"* N term for central HCAL: "<<""
|
---|
268 | << left << setw(30) <<HAD_Nhcal <<""<< right << setw(10)<<"*"<<"\n";
|
---|
269 | f_out << left << setw(30) <<"* C term for central HCAL: "<<""
|
---|
270 | << left << setw(30) <<HAD_Chcal <<""<< right << setw(10)<<"*"<<"\n";
|
---|
271 | f_out << left << setw(30) <<"* S term for forward HCAL: "<<""
|
---|
272 | << left << setw(30) <<HAD_Shf <<""<< right << setw(10)<<"*"<<"\n";
|
---|
273 | f_out << left << setw(30) <<"* N term for forward HCAL: "<<""
|
---|
274 | << left << setw(30) <<HAD_Nhf <<""<< right << setw(10)<<"*"<<"\n";
|
---|
275 | f_out << left << setw(30) <<"* C term for forward HCAL: "<<""
|
---|
276 | << left << setw(30) <<HAD_Chf <<""<< right << setw(10)<<"*"<<"\n";
|
---|
277 | f_out<<"* *"<<"\n";
|
---|
278 | f_out<<"#*************************** *"<<"\n";
|
---|
279 | f_out<<"# Tracking system acceptance *"<<"\n";
|
---|
280 | f_out<<"#*************************** *"<<"\n";
|
---|
281 | f_out<<"* *"<<"\n";
|
---|
282 | f_out << left << setw(55) <<"* Minimal pT needed to reach the calorimeter [GeV]: "<<""
|
---|
283 | << left << setw(10) <<PT_TRACKS_MIN <<""<< right << setw(5)<<"*"<<"\n";
|
---|
284 | f_out << left << setw(55) <<"* Efficiency associated to the tracking: "<<""
|
---|
285 | << left << setw(10) <<TRACKING_EFF <<""<< right << setw(5)<<"*"<<"\n";
|
---|
286 | f_out<<"* *"<<"\n";
|
---|
287 | f_out<<"#************************* *"<<"\n";
|
---|
288 | f_out<<"# Muon smearing parameters *"<<"\n";
|
---|
289 | f_out<<"#************************* *"<<"\n";
|
---|
290 | f_out<<"* *"<<"\n";
|
---|
291 | //MU_SmearPt 0.01
|
---|
292 | f_out<<"* *"<<"\n";
|
---|
293 | f_out<<"#****************************** *"<<"\n";
|
---|
294 | f_out<<"# Tau-jet definition parameters *"<<"\n";
|
---|
295 | f_out<<"#****************************** *"<<"\n";
|
---|
296 | f_out<<"* *"<<"\n";
|
---|
297 | f_out << left << setw(45) <<"* Cone radius for calorimeter tagging: "<<""
|
---|
298 | << left << setw(5) <<TAU_CONE_ENERGY <<""<< right << setw(20)<<"*"<<"\n";
|
---|
299 | f_out << left << setw(45) <<"* Fraction of energy in the small cone: "<<""
|
---|
300 | << left << setw(5) <<TAU_EM_COLLIMATION*100 <<""<< right << setw(20)<<"! not in datacard *"<<"\n";
|
---|
301 | f_out << left << setw(45) <<"* Cone radius for tracking tagging: "<<""
|
---|
302 | << left << setw(5) <<TAU_CONE_TRACKS <<""<< right << setw(20)<<"*"<<"\n";
|
---|
303 | f_out << left << setw(45) <<"* Minimum track pT [GeV]: "<<""
|
---|
304 | << left << setw(5) <<PT_TRACK_TAU <<""<< right << setw(20)<<"*"<<"\n";
|
---|
305 | f_out<<"* *"<<"\n";
|
---|
306 | f_out<<"#******************* *"<<"\n";
|
---|
307 | f_out<<"# Minimum pT's [GeV] *"<<"\n";
|
---|
308 | f_out<<"#******************* *"<<"\n";
|
---|
309 | f_out<<"* *"<<"\n";
|
---|
310 | f_out << left << setw(40) <<"* Minimum pT for electrons: "<<""
|
---|
311 | << left << setw(20) <<ELEC_pt <<""<< right << setw(10)<<"*"<<"\n";
|
---|
312 | f_out << left << setw(40) <<"* Minimum pT for muons: "<<""
|
---|
313 | << left << setw(20) <<MUON_pt <<""<< right << setw(10)<<"*"<<"\n";
|
---|
314 | f_out << left << setw(40) <<"* Minimum pT for jets: "<<""
|
---|
315 | << left << setw(20) <<JET_pt <<""<< right << setw(10)<<"*"<<"\n";
|
---|
316 | f_out << left << setw(40) <<"* Minimum pT for Tau-jets: "<<""
|
---|
317 | << left << setw(20) <<TAUJET_pt <<""<< right << setw(10)<<"*"<<"\n";
|
---|
318 | f_out<<"* *"<<"\n";
|
---|
319 | f_out<<"#*************************** *"<<"\n";
|
---|
320 | f_out<<"# B-tagging efficiencies [%] *"<<"\n";
|
---|
321 | f_out<<"#*************************** *"<<"\n";
|
---|
322 | f_out<<"* *"<<"\n";
|
---|
323 | f_out << left << setw(50) <<"* Efficiency to tag a \"b\" as a b-jet: "<<""
|
---|
324 | << left << setw(10) <<TAGGING_B <<""<< right << setw(10)<<"*"<<"\n";
|
---|
325 | f_out << left << setw(50) <<"* Efficiency to mistag a c-jet as a b-jet: "<<""
|
---|
326 | << left << setw(10) <<MISTAGGING_C <<""<< right << setw(10)<<"*"<<"\n";
|
---|
327 | f_out << left << setw(50) <<"* Efficiency to mistag a light jet as a b-jet: "<<""
|
---|
328 | << left << setw(10) <<MISTAGGING_L <<""<< right << setw(10)<<"*"<<"\n";
|
---|
329 | f_out<<"* *"<<"\n";
|
---|
330 | f_out<<"#*************** *"<<"\n";
|
---|
331 | f_out<<"# Jet definition *"<<"\n";
|
---|
332 | f_out<<"#*************** *"<<"\n";
|
---|
333 | f_out<<"* *"<<"\n";
|
---|
334 | f_out<<"* Six algorithms are currently available: *"<<"\n";
|
---|
335 | f_out<<"* - 1) CDF cone algorithm, *"<<"\n";
|
---|
336 | f_out<<"* - 2) CDF MidPoint algorithm, *"<<"\n";
|
---|
337 | f_out<<"* - 3) SIScone algorithm, *"<<"\n";
|
---|
338 | f_out<<"* - 4) kt algorithm, *"<<"\n";
|
---|
339 | f_out<<"* - 5) Cambrigde/Aachen algorithm, *"<<"\n";
|
---|
340 | f_out<<"* - 6) Anti-kt algorithm. *"<<"\n";
|
---|
341 | f_out<<"* *"<<"\n";
|
---|
342 | f_out<<"* You have chosen *"<<"\n";
|
---|
343 | switch(JETALGO) {
|
---|
344 | default:
|
---|
345 | case 1: {
|
---|
346 | f_out<<"* CDF JetClu jet algorithm with parameters: *"<<"\n";
|
---|
347 | f_out << left << setw(40) <<"* - Seed threshold: "<<""
|
---|
348 | << left << setw(10) <<SEEDTHRESHOLD <<""<< right << setw(20)<<"! not in datacard *"<<"\n";
|
---|
349 | f_out << left << setw(40) <<"* - Cone radius: "<<""
|
---|
350 | << left << setw(10) <<CONERADIUS <<""<< right << setw(20)<<"*"<<"\n";
|
---|
351 | f_out << left << setw(40) <<"* - Adjacency cut: "<<""
|
---|
352 | << left << setw(10) <<C_ADJACENCYCUT <<""<< right << setw(20)<<"! not in datacard *"<<"\n";
|
---|
353 | f_out << left << setw(40) <<"* - Max iterations: "<<""
|
---|
354 | << left << setw(10) <<C_MAXITERATIONS <<""<< right << setw(20)<<"! not in datacard *"<<"\n";
|
---|
355 | f_out << left << setw(40) <<"* - Iratch: "<<""
|
---|
356 | << left << setw(10) <<C_IRATCH <<""<< right << setw(20)<<"! not in datacard *"<<"\n";
|
---|
357 | f_out << left << setw(40) <<"* - Overlap threshold: "<<""
|
---|
358 | << left << setw(10) <<OVERLAPTHRESHOLD <<""<< right << setw(20)<<"! not in datacard *"<<"\n";
|
---|
359 | }
|
---|
360 | break;
|
---|
361 | case 2: {
|
---|
362 | f_out<<"* CDF midpoint jet algorithm with parameters: *"<<"\n";
|
---|
363 | f_out << left << setw(40) <<"* - Seed threshold: "<<""
|
---|
364 | << left << setw(20) <<SEEDTHRESHOLD <<""<< right << setw(10)<<"! not in datacard *"<<"\n";
|
---|
365 | f_out << left << setw(40) <<"* - Cone radius: "<<""
|
---|
366 | << left << setw(20) <<CONERADIUS <<""<< right << setw(10)<<"*"<<"\n";
|
---|
367 | f_out << left << setw(40) <<"* - Cone area fraction:"<<""
|
---|
368 | << left << setw(20) <<M_CONEAREAFRACTION <<""<< right << setw(10)<<"! not in datacard *"<<"\n";
|
---|
369 | f_out << left << setw(40) <<"* - Maximum pair size: "<<""
|
---|
370 | << left << setw(20) <<M_MAXPAIRSIZE <<""<< right << setw(10)<<"! not in datacard *"<<"\n";
|
---|
371 | f_out << left << setw(40) <<"* - Max iterations: "<<""
|
---|
372 | << left << setw(20) <<M_MAXITERATIONS <<""<< right << setw(10)<<"! not in datacard *"<<"\n";
|
---|
373 | f_out << left << setw(40) <<"* - Overlap threshold: "<<""
|
---|
374 | << left << setw(20) <<OVERLAPTHRESHOLD <<""<< right << setw(10)<<"! not in datacard *"<<"\n";
|
---|
375 | }
|
---|
376 | break;
|
---|
377 | case 3: {
|
---|
378 | f_out <<"* SISCone jet algorithm with parameters: *"<<"\n";
|
---|
379 | f_out << left << setw(40) <<"* - Cone radius: "<<""
|
---|
380 | << left << setw(20) <<CONERADIUS <<""<< right << setw(10)<<"*"<<"\n";
|
---|
381 | f_out << left << setw(40) <<"* - Overlap threshold: "<<""
|
---|
382 | << left << setw(20) <<OVERLAPTHRESHOLD <<""<< right << setw(10)<<"! not in datacard *"<<"\n";
|
---|
383 | f_out << left << setw(40) <<"* - Number pass max: "<<""
|
---|
384 | << left << setw(20) <<NPASS <<""<< right << setw(10)<<"! not in datacard *"<<"\n";
|
---|
385 | f_out << left << setw(40) <<"* - Minimum pT for protojet: "<<""
|
---|
386 | << left << setw(20) <<PROTOJET_PTMIN <<""<< right << setw(10)<<"! not in datacard *"<<"\n";
|
---|
387 | }
|
---|
388 | break;
|
---|
389 | case 4: {
|
---|
390 | f_out <<"* KT jet algorithm with parameters: *"<<"\n";
|
---|
391 | f_out << left << setw(40) <<"* - Cone radius: "<<""
|
---|
392 | << left << setw(20) <<CONERADIUS <<""<< right << setw(10)<<"*"<<"\n";
|
---|
393 | }
|
---|
394 | break;
|
---|
395 | case 5: {
|
---|
396 | f_out <<"* Cambridge/Aachen jet algorithm with parameters: *"<<"\n";
|
---|
397 | f_out << left << setw(40) <<"* - Cone radius: "<<""
|
---|
398 | << left << setw(20) <<CONERADIUS <<""<< right << setw(10)<<"*"<<"\n";
|
---|
399 | }
|
---|
400 | break;
|
---|
401 | case 6: {
|
---|
402 | f_out <<"* Anti-kt jet algorithm with parameters: *"<<"\n";
|
---|
403 | f_out << left << setw(40) <<"* - Cone radius: "<<""
|
---|
404 | << left << setw(20) <<CONERADIUS <<""<< right << setw(10)<<"*"<<"\n";
|
---|
405 | }
|
---|
406 | break;
|
---|
407 |
|
---|
408 |
|
---|
409 | }
|
---|
410 | f_out<<"* *"<<"\n";
|
---|
411 | f_out<<"#....................................................................*"<<"\n";
|
---|
412 | f_out<<"#>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>"<<"\n";
|
---|
413 |
|
---|
414 | }
|
---|
415 |
|
---|
416 | // **********Provides the smeared TLorentzVector for the electrons********
|
---|
417 | // Smears the electron energy, and changes the 4-momentum accordingly
|
---|
418 | // different smearing if the electron is central (eta < 2.5) or forward
|
---|
419 | void RESOLution::SmearElectron(TLorentzVector &electron) {
|
---|
420 | // the 'electron' variable will be changed by the function
|
---|
421 | float energy = electron.E(); // before smearing
|
---|
422 | float energyS = 0.0; // after smearing // \sigma/E = C + N/E + S/\sqrt{E}
|
---|
423 |
|
---|
424 | if(fabs(electron.Eta()) < MAX_TRACKER) { // if the electron is inside the tracker
|
---|
425 | energyS = gRandom->Gaus(energy, sqrt(
|
---|
426 | pow(ELG_Ncen,2) +
|
---|
427 | pow(ELG_Ccen*energy,2) +
|
---|
428 | pow(ELG_Scen*sqrt(energy),2) ));
|
---|
429 | }
|
---|
430 | if(fabs(electron.Eta()) > MAX_TRACKER && fabs(electron.Eta()) < MAX_CALO_FWD){
|
---|
431 | energyS = gRandom->Gaus(energy, sqrt(
|
---|
432 | pow(ELG_Nfwd,2) +
|
---|
433 | pow(ELG_Cfwd*energy,2) +
|
---|
434 | pow(ELG_Sfwd*sqrt(energy),2) ) );
|
---|
435 | }
|
---|
436 | electron.SetPtEtaPhiE(energyS/cosh(electron.Eta()), electron.Eta(), electron.Phi(), energyS);
|
---|
437 | if(electron.E() < 0)electron.SetPxPyPzE(0,0,0,0); // no negative values after smearing !
|
---|
438 | }
|
---|
439 |
|
---|
440 |
|
---|
441 | // **********Provides the smeared TLorentzVector for the muons********
|
---|
442 | // Smears the muon pT and changes the 4-momentum accordingly
|
---|
443 | void RESOLution::SmearMu(TLorentzVector &muon) {
|
---|
444 | // the 'muon' variable will be changed by the function
|
---|
445 | float pt = muon.Pt(); // before smearing
|
---|
446 | float ptS=pt;
|
---|
447 |
|
---|
448 | if(fabs(muon.Eta()) < MAX_MU )
|
---|
449 | {
|
---|
450 | ptS = gRandom->Gaus(pt, MU_SmearPt*pt ); // after smearing // \sigma/E = C + N/E + S/\sqrt{E}
|
---|
451 | }
|
---|
452 | muon.SetPtEtaPhiE(ptS, muon.Eta(), muon.Phi(), ptS*cosh(muon.Eta()));
|
---|
453 |
|
---|
454 | if(muon.E() < 0)muon.SetPxPyPzE(0,0,0,0); // no negative values after smearing !
|
---|
455 | }
|
---|
456 |
|
---|
457 |
|
---|
458 | // **********Provides the smeared TLorentzVector for the hadrons********
|
---|
459 | // Smears the hadron 4-momentum
|
---|
460 | void RESOLution::SmearHadron(TLorentzVector &hadron, const float frac)
|
---|
461 | // the 'hadron' variable will be changed by the function
|
---|
462 | // the 'frac' variable describes the long-living particles. Should be 0.7 for K0S and Lambda, 1. otherwise
|
---|
463 | {
|
---|
464 | float energy = hadron.E(); // before smearing
|
---|
465 | float energyS = 0.0; // after smearing // \sigma/E = C + N/E + S/\sqrt{E}
|
---|
466 | float energy_ecal = (1.0 - frac)*energy; // electromagnetic calorimeter
|
---|
467 | float energy_hcal = frac*energy; // hadronic calorimeter
|
---|
468 | // frac takes into account the decay of long-living particles, that decay in the calorimeters
|
---|
469 | // some of the particles decay mostly in the ecal, some mostly in the hcal
|
---|
470 |
|
---|
471 | float energyS1,energyS2;
|
---|
472 | if(fabs(hadron.Eta()) < MAX_CALO_CEN) {
|
---|
473 | energyS1 = gRandom->Gaus(energy_hcal, sqrt(
|
---|
474 | pow(HAD_Nhcal,2) +
|
---|
475 | pow(HAD_Chcal*energy_hcal,2) +
|
---|
476 | pow(HAD_Shcal*sqrt(energy_hcal),2) )) ;
|
---|
477 |
|
---|
478 |
|
---|
479 | energyS2 = gRandom->Gaus(energy_ecal, sqrt(
|
---|
480 | pow(ELG_Ncen,2) +
|
---|
481 | pow(ELG_Ccen*energy_ecal,2) +
|
---|
482 | pow(ELG_Scen*sqrt(energy_ecal),2) ) );
|
---|
483 |
|
---|
484 | energyS = ((energyS1>0)?energyS1:0) + ((energyS2>0)?energyS2:0);
|
---|
485 | }
|
---|
486 | if(abs(hadron.Eta()) > MAX_CALO_CEN && fabs(hadron.Eta()) < MAX_CALO_FWD){
|
---|
487 | energyS = gRandom->Gaus(energy, sqrt(
|
---|
488 | pow(HAD_Nhf,2) +
|
---|
489 | pow(HAD_Chf*energy,2) +
|
---|
490 | pow(HAD_Shf*sqrt(energy),2) ));
|
---|
491 | }
|
---|
492 |
|
---|
493 |
|
---|
494 |
|
---|
495 | hadron.SetPtEtaPhiE(energyS/cosh(hadron.Eta()),hadron.Eta(), hadron.Phi(), energyS);
|
---|
496 |
|
---|
497 | if(hadron.E() < 0)hadron.SetPxPyPzE(0,0,0,0);
|
---|
498 | }
|
---|
499 |
|
---|
500 | // **********Provides the energy in the cone of radius TAU_CONE_ENERGY for the tau identification********
|
---|
501 | // to be taken into account, a calo tower should
|
---|
502 | // 1) have a transverse energy \f$ E_T = \sqrt{E_X^2 + E_Y^2} \f$ above a given threshold
|
---|
503 | // 2) be inside a cone with a radius R and the axis defined by (eta,phi)
|
---|
504 | double RESOLution::EnergySmallCone(const vector<PhysicsTower> &towers, const float eta, const float phi) {
|
---|
505 | double Energie=0;
|
---|
506 | for(unsigned int i=0; i < towers.size(); i++) {
|
---|
507 | if(towers[i].fourVector.pt() < SEEDTHRESHOLD) continue;
|
---|
508 | if((DeltaR(phi,eta,towers[i].fourVector.phi(),towers[i].fourVector.eta()) < TAU_CONE_ENERGY)) {
|
---|
509 | Energie += towers[i].fourVector.E;
|
---|
510 | }
|
---|
511 | }
|
---|
512 | return Energie;
|
---|
513 | }
|
---|
514 |
|
---|
515 |
|
---|
516 | // **********Provides the number of tracks in the cone of radius TAU_CONE_TRACKS for the tau identification********
|
---|
517 | // to be taken into account, a track should
|
---|
518 | // 1) avec a transverse momentum \$f p_T \$ above a given threshold
|
---|
519 | // 2) be inside a cone with a radius R and the axis defined by (eta,phi)
|
---|
520 | // IMPORTANT REMARK !!!!!
|
---|
521 | // previously, the argument 'phi' was before the argument 'eta'
|
---|
522 | // this has been changed for consistency with the other functions
|
---|
523 | // double check your running code that uses NumTracks !
|
---|
524 | unsigned int RESOLution::NumTracks(const vector<TLorentzVector> &tracks, const float pt_track, const float eta, const float phi) {
|
---|
525 | unsigned int numtrack=0;
|
---|
526 | for(unsigned int i=0; i < tracks.size(); i++) {
|
---|
527 | if((tracks[i].Pt() < pt_track )||
|
---|
528 | (DeltaR(phi,eta,tracks[i].Phi(),tracks[i].Eta()) > TAU_CONE_TRACKS)
|
---|
529 | )continue;
|
---|
530 | numtrack++;
|
---|
531 | }
|
---|
532 | return numtrack;
|
---|
533 | }
|
---|
534 |
|
---|
535 |
|
---|
536 | //*** Returns the PID of the particle with the highest energy, in a cone with a radius CONERADIUS and an axis (eta,phi) *********
|
---|
537 | //used by Btaggedjet
|
---|
538 | ///// Attention : bug removed => CONERADIUS/2 -> CONERADIUS !!
|
---|
539 | int RESOLution::Bjets(const TSimpleArray<TRootGenParticle> &subarray, const float eta, const float phi) {
|
---|
540 | float emax=0;
|
---|
541 | int Ppid=0;
|
---|
542 | if(subarray.GetEntries()>0) {
|
---|
543 | for(int i=0; i < subarray.GetEntries();i++) { // should have pt>PT_JETMIN and a small cone radius (r<CONE_JET)
|
---|
544 | float genDeltaR = DeltaR(subarray[i]->Phi,subarray[i]->Eta,phi,eta);
|
---|
545 | if(genDeltaR < CONERADIUS && subarray[i]->E > emax) {
|
---|
546 | emax=subarray[i]->E;
|
---|
547 | Ppid=abs(subarray[i]->PID);
|
---|
548 | }
|
---|
549 | }
|
---|
550 | }
|
---|
551 | return Ppid;
|
---|
552 | }
|
---|
553 |
|
---|
554 |
|
---|
555 | //******************** Simulates the b-tagging efficiency for real bjet, or the misendentification for other jets****************
|
---|
556 | bool RESOLution::Btaggedjet(const TLorentzVector &JET, const TSimpleArray<TRootGenParticle> &subarray) {
|
---|
557 | if( rand()%100 < (TAGGING_B+1) && Bjets(subarray,JET.Eta(),JET.Phi())==pB ) return true; // b-tag of b-jets is 40%
|
---|
558 | else if( rand()%100 < (MISTAGGING_C+1) && Bjets(subarray,JET.Eta(),JET.Phi())==pC ) return true; // b-tag of c-jets is 10%
|
---|
559 | else if( rand()%100 < (MISTAGGING_L+1) && Bjets(subarray,JET.Eta(),JET.Phi())!=0) return true; // b-tag of light jets is 1%
|
---|
560 | return false;
|
---|
561 | }
|
---|
562 |
|
---|
563 | //***********************Isolation criteria***********************
|
---|
564 | //****************************************************************
|
---|
565 | bool RESOLution::Isolation(Float_t phi,Float_t eta,const vector<TLorentzVector> &tracks,float PT_TRACK2)
|
---|
566 | {
|
---|
567 | bool isolated = false;
|
---|
568 | Float_t deltar=5000.; // Initial value; should be high; no further repercussion
|
---|
569 | // loop on all final charged particles, with p_t >2, close enough from the electron
|
---|
570 | for(unsigned int i=0; i < tracks.size(); i++)
|
---|
571 | {
|
---|
572 | if(tracks[i].Pt() < PT_TRACK2)continue;
|
---|
573 | Float_t genDeltaR = DeltaR(phi,eta,tracks[i].Phi(),tracks[i].Eta()); // slower to evaluate
|
---|
574 | if(
|
---|
575 | (genDeltaR > deltar) ||
|
---|
576 | (genDeltaR==0)
|
---|
577 | ) continue ;
|
---|
578 | deltar=genDeltaR;
|
---|
579 | }
|
---|
580 | if(deltar > 0.5)isolated = true; // returns the closest distance
|
---|
581 | return isolated;
|
---|
582 | }
|
---|
583 |
|
---|
584 |
|
---|
585 | //**************************** Returns the delta Phi ****************************
|
---|
586 | float DeltaPhi(const float phi1, const float phi2) {
|
---|
587 | float deltaphi=phi1-phi2; // in here, -PI < phi < PI
|
---|
588 | if(fabs(deltaphi) > PI) deltaphi=2.*PI-fabs(deltaphi);// put deltaphi between 0 and PI
|
---|
589 | else deltaphi=fabs(deltaphi);
|
---|
590 |
|
---|
591 | return deltaphi;
|
---|
592 | }
|
---|
593 |
|
---|
594 | //**************************** Returns the delta R****************************
|
---|
595 | float DeltaR(const float phi1, const float eta1, const float phi2, const float eta2) {
|
---|
596 | return sqrt(pow(DeltaPhi(phi1,phi2),2) + pow(eta1-eta2,2));
|
---|
597 | }
|
---|
598 |
|
---|
599 | int sign(const int myint) {
|
---|
600 | if (myint >0) return 1;
|
---|
601 | else if (myint <0) return -1;
|
---|
602 | else return 0;
|
---|
603 | }
|
---|
604 |
|
---|
605 | int sign(const float myfloat) {
|
---|
606 | if (myfloat >0) return 1;
|
---|
607 | else if (myfloat <0) return -1;
|
---|
608 | else return 0;
|
---|
609 | }
|
---|
610 |
|
---|
611 | int Charge(int pid)
|
---|
612 | {
|
---|
613 | int charge;
|
---|
614 | if(
|
---|
615 | (pid == pGAMMA) ||
|
---|
616 | (pid == pPI0) ||
|
---|
617 | (pid == pK0L) ||
|
---|
618 | (pid == pN) ||
|
---|
619 | (pid == pSIGMA0) ||
|
---|
620 | (pid == pDELTA0) ||
|
---|
621 | (pid == pK0S) // not charged particles : invisible by tracker
|
---|
622 | )
|
---|
623 | charge = 0;
|
---|
624 | else charge = (sign(pid));
|
---|
625 | return charge;
|
---|
626 |
|
---|
627 | }
|
---|