Fork me on GitHub

source: svn/trunk/src/SmearUtil.cc@ 62

Last change on this file since 62 was 62, checked in by severine ovyn, 16 years ago

Read BField info from datacard

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