Fork me on GitHub

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

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

add parameters for RP

File size: 52.9 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 "SmearUtil.h"
16#include "TRandom.h"
17
18#include <iostream>
19#include <fstream>
20#include <sstream>
21#include <iomanip>
22using namespace std;
23
24
25ParticleUtil::ParticleUtil(const TLorentzVector &genMomentum, int pid) {
26 _pid=pid;
27 _e = genMomentum.E();
28 _px = genMomentum.Px();
29 _py = genMomentum.Py();
30 _pz = genMomentum.Pz();
31 _pt = genMomentum.Pt();
32
33 //_e, _px, _py, _pz, _pt;
34 //float _eta, _etaCalo, _phi, _phiCalo;
35 //int _pid;
36}
37
38//------------------------------------------------------------------------------
39
40RESOLution::RESOLution() {
41
42 // Detector characteristics
43 CEN_max_tracker = 2.5; // Maximum tracker coverage
44 CEN_max_calo_cen = 3.0; // central calorimeter coverage
45 CEN_max_calo_fwd = 5.0; // forward calorimeter pseudorapidity coverage
46 CEN_max_mu = 2.4; // muon chambers pseudorapidity coverage
47
48 // Energy resolution for electron/photon
49 // \sigma/E = C + N/E + S/\sqrt{E}
50 ELG_Scen = 0.05; // S term for central ECAL
51 ELG_Ncen = 0.25; // N term for central ECAL
52 ELG_Ccen = 0.005; // C term for central ECAL
53 ELG_Cfwd = 0.107; // S term for forward ECAL
54 ELG_Sfwd = 2.084; // C term for forward ECAL
55 ELG_Nfwd = 0.0; // N term for central ECAL
56
57 // Energy resolution for hadrons in ecal/hcal/hf
58 // \sigma/E = C + N/E + S/\sqrt{E}
59 HAD_Shcal = 1.5; // S term for central HCAL // hadronic calorimeter
60 HAD_Nhcal = 0.; // N term for central HCAL
61 HAD_Chcal = 0.05; // C term for central HCAL
62 HAD_Shf = 2.7; // S term for HF // forward calorimeter
63 HAD_Nhf = 0.; // N term for HF
64 HAD_Chf = 0.13; // C term for HF
65
66 // Muon smearing
67 MU_SmearPt = 0.01;
68
69 // Tracking efficiencies
70 TRACK_ptmin = 0.9; // minimal pt needed to reach the calorimeter in GeV
71 TRACK_eff = 100; // efficiency associated to the tracking
72
73 // Calorimetric towers
74 TOWER_number = 40;
75 const float tower_eta_edges[41] = {
76 0., 0.087, 0.174, 0.261, 0.348, 0.435, 0.522, 0.609, 0.696, 0.783, 0.870, 0.957, 1.044, 1.131, 1.218, 1.305, 1.392, 1.479, 1.566,
77 1.653, 1.740, 1.830, 1.930, 2.043, 2.172, 2.322, 2.500, 2.650, 2.868, 2.950, 3.125, 3.300, 3.475, 3.650, 3.825, 4.000, 4.175,
78 4.350, 4.525, 4.700, 5.000}; // temporary object
79 TOWER_eta_edges = new float[TOWER_number+1];
80 for(unsigned int i=0; i<TOWER_number +1; i++) TOWER_eta_edges[i] = tower_eta_edges[i];
81
82 const float tower_dphi[40] = {
83 5, 5, 5, 5, 5, 5, 5, 5, 5, 5, 5, 5, 5, 5, 5, 5, 5, 5, 5, 10,
84 10,10,10,10,10, 10,10,10,10,10, 10,10,10,10,10, 10,10,10,20, 20 }; // temporary object
85 TOWER_dphi = new float[TOWER_number];
86 for(unsigned int i=0; i<TOWER_number; i++) TOWER_dphi[i] = tower_dphi[i];
87
88
89 // Thresholds for reconstructed objetcs
90 PTCUT_elec = 10.0;
91 PTCUT_muon = 10.0;
92 PTCUT_jet = 20.0;
93 PTCUT_gamma = 10.0;
94 PTCUT_taujet = 10.0;
95
96 // General jet variable
97 JET_coneradius = 0.7; // generic jet radius ; not for tau's !!!
98 JET_jetalgo = 1; // 1 for Cone algorithm, 2 for MidPoint algorithm, 3 for SIScone algorithm, 4 for kt algorithm
99 JET_seed = 1.0; // minimum seed to start jet reconstruction
100
101 // Tagging definition
102 BTAG_b = 40;
103 BTAG_mistag_c = 10;
104 BTAG_mistag_l = 1;
105
106 // FLAGS
107 FLAG_bfield = 1; //1 to run the bfield propagation else 0
108 FLAG_vfd = 1; //1 to run the very forward detectors else 0
109 FLAG_trigger = 1; //1 to run the trigger selection else 0
110 FLAG_frog = 1; //1 to run the FROG event display
111
112 // In case BField propagation allowed
113 TRACK_radius = 129; //radius of the BField coverage
114 TRACK_length = 300; //length of the BField coverage
115 TRACK_bfield_x = 0; //X composant of the BField
116 TRACK_bfield_y = 0; //Y composant of the BField
117 TRACK_bfield_z = 3.8; //Z composant of the BField
118
119 // In case Very forward detectors allowed
120 VFD_min_calo_vfd = 5.2; // very forward calorimeter (if any) like CASTOR
121 VFD_max_calo_vfd = 6.6;
122 VFD_min_zdc = 8.3;
123 VFD_s_zdc = 140; // distance of the Zero Degree Calorimeter, from the Interaction poin, in [m]
124
125 RP_220_s = 220; // distance of the RP to the IP, in meters
126 RP_220_x = 0.002; // distance of the RP to the beam, in meters
127 RP_420_s = 420; // distance of the RP to the IP, in meters
128 RP_420_x = 0.004; // distance of the RP to the beam, in meters
129 RP_beam1Card = "data/LHCB1IR5_v6.500.tfs";
130 RP_beam2Card = "data/LHCB1IR5_v6.500.tfs";
131
132 // In case FROG event display allowed
133 NEvents_Frog = 10;
134
135 //********************************************
136 //jet stuffs not defined in the input datacard
137 //********************************************
138
139 JET_overlap = 0.75;
140 // MidPoint algorithm definition
141 JET_M_coneareafraction = 0.25;
142 JET_M_maxpairsize = 2;
143 JET_M_maxiterations = 100;
144 // Define Cone algorithm.
145 JET_C_adjacencycut = 2;
146 JET_C_maxiterations = 100;
147 JET_C_iratch = 1;
148 //Define SISCone algorithm.
149 JET_S_npass = 0;
150 JET_S_protojet_ptmin= 0.0;
151
152 //For Tau-jet definition
153 TAU_energy_scone = 0.15; // radius R of the cone for tau definition, based on energy threshold
154 TAU_track_scone = 0.4; // radius R of the cone for tau definition, based on track number
155 TAU_track_pt = 2; // minimal pt [GeV] for tracks to be considered in tau definition
156 TAU_energy_frac = 0.95; // fraction of energy required in the central part of the cone, for tau jets
157
158 PT_QUARKS_MIN = 2.0 ; // minimal pt needed by quarks to do b-tag
159
160 //for very forward detectors
161 RP_offsetEl_s = 120;
162 RP_offsetEl_x = 0.097;
163
164
165}
166
167
168RESOLution::RESOLution(const RESOLution & DET) {
169 // Detector characteristics
170 CEN_max_tracker = DET.CEN_max_tracker;
171 CEN_max_calo_cen = DET.CEN_max_calo_cen;
172 CEN_max_calo_fwd = DET.CEN_max_calo_fwd;
173 CEN_max_mu = DET.CEN_max_mu;
174
175 // Energy resolution for electron/photon
176 ELG_Scen = DET.ELG_Scen;
177 ELG_Ncen = DET.ELG_Ncen;
178 ELG_Ccen = DET.ELG_Ccen;
179 ELG_Cfwd = DET.ELG_Cfwd;
180 ELG_Sfwd = DET.ELG_Sfwd;
181 ELG_Nfwd = DET.ELG_Nfwd;
182
183 // Energy resolution for hadrons in ecal/hcal/hf
184 HAD_Shcal = DET.HAD_Shcal;
185 HAD_Nhcal = DET.HAD_Nhcal;
186 HAD_Chcal = DET.HAD_Chcal;
187 HAD_Shf = DET.HAD_Shf;
188 HAD_Nhf = DET.HAD_Nhf;
189 HAD_Chf = DET.HAD_Chf;
190
191 // Muon smearing
192 MU_SmearPt = DET.MU_SmearPt;
193
194 // Tracking efficiencies
195 TRACK_ptmin = DET.TRACK_ptmin;
196 TRACK_eff = DET.TRACK_eff;
197
198 // Calorimetric towers
199 TOWER_number = DET.TOWER_number;
200 TOWER_eta_edges = new float[TOWER_number+1];
201 for(unsigned int i=0; i<TOWER_number +1; i++) TOWER_eta_edges[i] = DET.TOWER_eta_edges[i];
202
203 TOWER_dphi = new float[TOWER_number];
204 for(unsigned int i=0; i<TOWER_number; i++) TOWER_dphi[i] = DET.TOWER_dphi[i];
205
206 // Thresholds for reconstructed objetcs
207 PTCUT_elec = DET.PTCUT_elec;
208 PTCUT_muon = DET.PTCUT_muon;
209 PTCUT_jet = DET.PTCUT_jet;
210 PTCUT_gamma = DET.PTCUT_gamma;
211 PTCUT_taujet = DET.PTCUT_taujet;
212
213 // General jet variable
214 JET_coneradius = DET.JET_coneradius;
215 JET_jetalgo = DET.JET_jetalgo;
216 JET_seed = DET.JET_seed;
217
218 // Tagging definition
219 BTAG_b = DET.BTAG_b;
220 BTAG_mistag_c = DET.BTAG_mistag_c;
221 BTAG_mistag_l = DET.BTAG_mistag_l;
222
223 // FLAGS
224 FLAG_bfield = DET.FLAG_bfield;
225 FLAG_vfd = DET.FLAG_vfd;
226 FLAG_trigger = DET.FLAG_trigger;
227 FLAG_frog = DET.FLAG_frog;
228
229 // In case BField propagation allowed
230 TRACK_radius = DET.TRACK_radius;
231 TRACK_length = DET.TRACK_length;
232 TRACK_bfield_x = DET.TRACK_bfield_x;
233 TRACK_bfield_y = DET.TRACK_bfield_y;
234 TRACK_bfield_z = DET.TRACK_bfield_z;
235
236 // In case Very forward detectors allowed
237 VFD_min_calo_vfd = DET.VFD_min_calo_vfd;
238 VFD_max_calo_vfd = DET.VFD_max_calo_vfd;
239 VFD_min_zdc = DET.VFD_min_zdc;
240 VFD_s_zdc = DET.VFD_s_zdc;
241
242 RP_220_s = DET.RP_220_s;
243 RP_220_x = DET.RP_220_x;
244 RP_420_s = DET.RP_420_s;
245 RP_420_x = DET.RP_420_x;
246 RP_beam1Card = DET.RP_beam1Card;
247 RP_beam2Card = DET.RP_beam2Card;
248 RP_offsetEl_s = DET.RP_offsetEl_s;
249 RP_offsetEl_x = DET.RP_offsetEl_x;
250
251 // In case FROG event display allowed
252 NEvents_Frog = DET.NEvents_Frog;
253
254 JET_overlap = DET.JET_overlap;
255 // MidPoint algorithm definition
256 JET_M_coneareafraction = DET.JET_M_coneareafraction;
257 JET_M_maxpairsize = DET.JET_M_maxpairsize;
258 JET_M_maxiterations = DET.JET_M_maxiterations;
259 // Define Cone algorithm.
260 JET_C_adjacencycut = DET.JET_C_adjacencycut;
261 JET_C_maxiterations = DET.JET_C_maxiterations;
262 JET_C_iratch = DET.JET_C_iratch;
263 //Define SISCone algorithm.
264 JET_S_npass = DET.JET_S_npass;
265 JET_S_protojet_ptmin = DET.JET_S_protojet_ptmin;
266
267 //For Tau-jet definition
268 TAU_energy_scone = DET.TAU_energy_scone;
269 TAU_track_scone = DET.TAU_track_scone;
270 TAU_track_pt = DET.TAU_track_pt;
271 TAU_energy_frac = DET.TAU_energy_frac;
272
273 PT_QUARKS_MIN = DET.PT_QUARKS_MIN;
274}
275
276RESOLution& RESOLution::operator=(const RESOLution& DET) {
277 if(this==&DET) return *this;
278 // Detector characteristics
279 CEN_max_tracker = DET.CEN_max_tracker;
280 CEN_max_calo_cen = DET.CEN_max_calo_cen;
281 CEN_max_calo_fwd = DET.CEN_max_calo_fwd;
282 CEN_max_mu = DET.CEN_max_mu;
283
284 // Energy resolution for electron/photon
285 ELG_Scen = DET.ELG_Scen;
286 ELG_Ncen = DET.ELG_Ncen;
287 ELG_Ccen = DET.ELG_Ccen;
288 ELG_Cfwd = DET.ELG_Cfwd;
289 ELG_Sfwd = DET.ELG_Sfwd;
290 ELG_Nfwd = DET.ELG_Nfwd;
291
292 // Energy resolution for hadrons in ecal/hcal/hf
293 HAD_Shcal = DET.HAD_Shcal;
294 HAD_Nhcal = DET.HAD_Nhcal;
295 HAD_Chcal = DET.HAD_Chcal;
296 HAD_Shf = DET.HAD_Shf;
297 HAD_Nhf = DET.HAD_Nhf;
298 HAD_Chf = DET.HAD_Chf;
299
300 // Muon smearing
301 MU_SmearPt = DET.MU_SmearPt;
302
303 // Tracking efficiencies
304 TRACK_ptmin = DET.TRACK_ptmin;
305 TRACK_eff = DET.TRACK_eff;
306
307 // Calorimetric towers
308 TOWER_number = DET.TOWER_number;
309 TOWER_eta_edges = new float[TOWER_number+1];
310 for(unsigned int i=0; i<TOWER_number +1; i++) TOWER_eta_edges[i] = DET.TOWER_eta_edges[i];
311
312 TOWER_dphi = new float[TOWER_number];
313 for(unsigned int i=0; i<TOWER_number; i++) TOWER_dphi[i] = DET.TOWER_dphi[i];
314
315 // Thresholds for reconstructed objetcs
316 PTCUT_elec = DET.PTCUT_elec;
317 PTCUT_muon = DET.PTCUT_muon;
318 PTCUT_jet = DET.PTCUT_jet;
319 PTCUT_gamma = DET.PTCUT_gamma;
320 PTCUT_taujet = DET.PTCUT_taujet;
321
322 // General jet variable
323 JET_coneradius = DET.JET_coneradius;
324 JET_jetalgo = DET.JET_jetalgo;
325 JET_seed = DET.JET_seed;
326
327 // Tagging definition
328 BTAG_b = DET.BTAG_b;
329 BTAG_mistag_c = DET.BTAG_mistag_c;
330 BTAG_mistag_l = DET.BTAG_mistag_l;
331
332 // FLAGS
333 FLAG_bfield = DET.FLAG_bfield;
334 FLAG_vfd = DET.FLAG_vfd;
335 FLAG_trigger = DET.FLAG_trigger;
336 FLAG_frog = DET.FLAG_frog;
337
338 // In case BField propagation allowed
339 TRACK_radius = DET.TRACK_radius;
340 TRACK_length = DET.TRACK_length;
341 TRACK_bfield_x = DET.TRACK_bfield_x;
342 TRACK_bfield_y = DET.TRACK_bfield_y;
343 TRACK_bfield_z = DET.TRACK_bfield_z;
344
345 // In case Very forward detectors allowed
346 VFD_min_calo_vfd = DET.VFD_min_calo_vfd;
347 VFD_max_calo_vfd = DET.VFD_max_calo_vfd;
348 VFD_min_zdc = DET.VFD_min_zdc;
349 VFD_s_zdc = DET.VFD_s_zdc;
350
351 RP_220_s = DET.RP_220_s;
352 RP_220_x = DET.RP_220_x;
353 RP_420_s = DET.RP_420_s;
354 RP_420_x = DET.RP_420_x;
355 RP_offsetEl_s = DET.RP_offsetEl_s;
356 RP_offsetEl_x = DET.RP_offsetEl_x;
357 RP_beam1Card = DET.RP_beam1Card;
358 RP_beam2Card = DET.RP_beam2Card;
359
360
361 // In case FROG event display allowed
362 NEvents_Frog = DET.NEvents_Frog;
363
364 JET_overlap = DET.JET_overlap;
365 // MidPoint algorithm definition
366 JET_M_coneareafraction = DET.JET_M_coneareafraction;
367 JET_M_maxpairsize = DET.JET_M_maxpairsize;
368 JET_M_maxiterations = DET.JET_M_maxiterations;
369 // Define Cone algorithm.
370 JET_C_adjacencycut = DET.JET_C_adjacencycut;
371 JET_C_maxiterations = DET.JET_C_maxiterations;
372 JET_C_iratch = DET.JET_C_iratch;
373 //Define SISCone algorithm.
374 JET_S_npass = DET.JET_S_npass;
375 JET_S_protojet_ptmin = DET.JET_S_protojet_ptmin;
376
377 //For Tau-jet definition
378 TAU_energy_scone = DET.TAU_energy_scone;
379 TAU_track_scone = DET.TAU_track_scone;
380 TAU_track_pt = DET.TAU_track_pt;
381 TAU_energy_frac = DET.TAU_energy_frac;
382
383 PT_QUARKS_MIN = DET.PT_QUARKS_MIN;
384 return *this;
385}
386
387
388
389
390//------------------------------------------------------------------------------
391void RESOLution::ReadDataCard(const string datacard) {
392
393 string temp_string;
394 istringstream curstring;
395
396 ifstream fichier_a_lire(datacard.c_str());
397 if(!fichier_a_lire.good()) {
398 cout <<"** WARNING: Datadard not found, use default values **" << endl;
399 return;
400 }
401
402 while (getline(fichier_a_lire,temp_string)) {
403 curstring.clear(); // needed when using several times istringstream::str(string)
404 curstring.str(temp_string);
405 string varname;
406 float value; int ivalue; string svalue;
407
408 if(strstr(temp_string.c_str(),"#")) { }
409 else if(strstr(temp_string.c_str(),"CEN_max_tracker")) {curstring >> varname >> value; CEN_max_tracker = value;}
410 else if(strstr(temp_string.c_str(),"CEN_max_calo_cen")) {curstring >> varname >> value; CEN_max_calo_cen = value;}
411 else if(strstr(temp_string.c_str(),"CEN_max_calo_fwd")) {curstring >> varname >> value; CEN_max_calo_fwd = value;}
412 else if(strstr(temp_string.c_str(),"CEN_max_mu")) {curstring >> varname >> value; CEN_max_mu = value;}
413
414 else if(strstr(temp_string.c_str(),"VFD_min_calo_vfd")) {curstring >> varname >> value; VFD_min_calo_vfd = value;}
415 else if(strstr(temp_string.c_str(),"VFD_max_calo_vfd")) {curstring >> varname >> value; VFD_max_calo_vfd = value;}
416 else if(strstr(temp_string.c_str(),"VFD_min_zdc")) {curstring >> varname >> value; VFD_min_zdc = value;}
417 else if(strstr(temp_string.c_str(),"VFD_s_zdc")) {curstring >> varname >> value; VFD_s_zdc = value;}
418
419 else if(strstr(temp_string.c_str(),"RP_220_s")) {curstring >> varname >> value; RP_220_s = value;}
420 else if(strstr(temp_string.c_str(),"RP_220_x")) {curstring >> varname >> value; RP_220_x = value;}
421 else if(strstr(temp_string.c_str(),"RP_420_s")) {curstring >> varname >> value; RP_420_s = value;}
422 else if(strstr(temp_string.c_str(),"RP_420_x")) {curstring >> varname >> value; RP_420_x = value;}
423 //else if(strstr(temp_string.c_str(),"Beam1Card")) {curstring >> varname >> svalue;Beam1Card = svalue;}
424 //else if(strstr(temp_string.c_str(),"Beam2Card")) {curstring >> varname >> svalue;Beam2Card = svalue;}
425
426 else if(strstr(temp_string.c_str(),"ELG_Scen")) {curstring >> varname >> value; ELG_Scen = value;}
427 else if(strstr(temp_string.c_str(),"ELG_Ncen")) {curstring >> varname >> value; ELG_Ncen = value;}
428 else if(strstr(temp_string.c_str(),"ELG_Ccen")) {curstring >> varname >> value; ELG_Ccen = value;}
429 else if(strstr(temp_string.c_str(),"ELG_Sfwd")) {curstring >> varname >> value; ELG_Sfwd = value;}
430 else if(strstr(temp_string.c_str(),"ELG_Cfwd")) {curstring >> varname >> value; ELG_Cfwd = value;}
431 else if(strstr(temp_string.c_str(),"ELG_Nfwd")) {curstring >> varname >> value; ELG_Nfwd = value;}
432 else if(strstr(temp_string.c_str(),"HAD_Shcal")) {curstring >> varname >> value; HAD_Shcal = value;}
433 else if(strstr(temp_string.c_str(),"HAD_Nhcal")) {curstring >> varname >> value; HAD_Nhcal = value;}
434 else if(strstr(temp_string.c_str(),"HAD_Chcal")) {curstring >> varname >> value; HAD_Chcal = value;}
435 else if(strstr(temp_string.c_str(),"HAD_Shf")) {curstring >> varname >> value; HAD_Shf = value;}
436 else if(strstr(temp_string.c_str(),"HAD_Nhf")) {curstring >> varname >> value; HAD_Nhf = value;}
437 else if(strstr(temp_string.c_str(),"HAD_Chf")) {curstring >> varname >> value; HAD_Chf = value;}
438 else if(strstr(temp_string.c_str(),"MU_SmearPt")) {curstring >> varname >> value; MU_SmearPt = value;}
439
440 else if(strstr(temp_string.c_str(),"TRACK_radius")) {curstring >> varname >> ivalue;TRACK_radius = ivalue;}
441 else if(strstr(temp_string.c_str(),"TRACK_length")) {curstring >> varname >> ivalue;TRACK_length = ivalue;}
442 else if(strstr(temp_string.c_str(),"TRACK_bfield_x")) {curstring >> varname >> value; TRACK_bfield_x = value;}
443 else if(strstr(temp_string.c_str(),"TRACK_bfield_y")) {curstring >> varname >> value; TRACK_bfield_y = value;}
444 else if(strstr(temp_string.c_str(),"TRACK_bfield_z")) {curstring >> varname >> value; TRACK_bfield_z = value;}
445 else if(strstr(temp_string.c_str(),"FLAG_bfield")) {curstring >> varname >> ivalue; FLAG_bfield = ivalue;}
446 else if(strstr(temp_string.c_str(),"TRACK_ptmin")) {curstring >> varname >> value; TRACK_ptmin = value;}
447 else if(strstr(temp_string.c_str(),"TRACK_eff")) {curstring >> varname >> ivalue;TRACK_eff = ivalue;}
448
449 else if(strstr(temp_string.c_str(),"TOWER_number")) {curstring >> varname >> ivalue;TOWER_number = ivalue;}
450 else if(strstr(temp_string.c_str(),"TOWER_eta_edges")){
451 curstring >> varname; for(unsigned int i=0; i<TOWER_number+1; i++) {curstring >> value; TOWER_eta_edges[i] = value;} }
452 else if(strstr(temp_string.c_str(),"TOWER_dphi")){
453 curstring >> varname; for(unsigned int i=0; i<TOWER_number; i++) {curstring >> value; TOWER_dphi[i] = value;} }
454
455 else if(strstr(temp_string.c_str(),"PTCUT_elec")) {curstring >> varname >> value; PTCUT_elec = value;}
456 else if(strstr(temp_string.c_str(),"PTCUT_muon")) {curstring >> varname >> value; PTCUT_muon = value;}
457 else if(strstr(temp_string.c_str(),"PTCUT_jet")) {curstring >> varname >> value; PTCUT_jet = value;}
458 else if(strstr(temp_string.c_str(),"PTCUT_gamma")) {curstring >> varname >> value; PTCUT_gamma = value;}
459 else if(strstr(temp_string.c_str(),"PTCUT_taujet")) {curstring >> varname >> value; PTCUT_taujet = value;}
460
461 else if(strstr(temp_string.c_str(),"JET_coneradius")) {curstring >> varname >> value; JET_coneradius = value;}
462 else if(strstr(temp_string.c_str(),"JET_jetalgo")) {curstring >> varname >> ivalue;JET_jetalgo = ivalue;}
463 else if(strstr(temp_string.c_str(),"JET_seed")) {curstring >> varname >> value; JET_seed = value;}
464
465 else if(strstr(temp_string.c_str(),"BTAG_b")) {curstring >> varname >> ivalue;BTAG_b = ivalue;}
466 else if(strstr(temp_string.c_str(),"BTAG_mistag_c")) {curstring >> varname >> ivalue;BTAG_mistag_c = ivalue;}
467 else if(strstr(temp_string.c_str(),"BTAG_mistag_l")) {curstring >> varname >> ivalue;BTAG_mistag_l = ivalue;}
468
469 else if(strstr(temp_string.c_str(),"FLAG_vfd")) {curstring >> varname >> ivalue; FLAG_vfd = ivalue;}
470 else if(strstr(temp_string.c_str(),"FLAG_trigger")) {curstring >> varname >> ivalue; FLAG_trigger = ivalue;}
471 else if(strstr(temp_string.c_str(),"FLAG_frog")) {curstring >> varname >> ivalue; FLAG_frog = ivalue;}
472 else if(strstr(temp_string.c_str(),"NEvents_Frog")) {curstring >> varname >> ivalue; NEvents_Frog = ivalue;}
473 }
474
475 //jet stuffs not defined in the input datacard
476 JET_overlap = 0.75;
477 // MidPoint algorithm definition
478 JET_M_coneareafraction = 0.25;
479 JET_M_maxpairsize = 2;
480 JET_M_maxiterations = 100;
481 // Define Cone algorithm.
482 JET_C_adjacencycut = 2;
483 JET_C_maxiterations = 100;
484 JET_C_iratch = 1;
485 //Define SISCone algorithm.
486 JET_S_npass = 0;
487 JET_S_protojet_ptmin= 0.0;
488
489 //For Tau-jet definition
490 TAU_energy_scone = 0.15; // radius R of the cone for tau definition, based on energy threshold
491 TAU_track_scone = 0.4; // radius R of the cone for tau definition, based on track number
492 TAU_track_pt = 2; // minimal pt [GeV] for tracks to be considered in tau definition
493 TAU_energy_frac = 0.95; // fraction of energy required in the central part of the cone, for tau jets
494
495}
496
497void RESOLution::Logfile(const string& LogName) {
498 //void RESOLution::Logfile(string outputfilename) {
499
500 ofstream f_out(LogName.c_str());
501
502 f_out<<"#*********************************************************************"<<"\n";
503 f_out<<"# *"<<"\n";
504 f_out<<"# ---- DELPHES release 1.0 ---- *"<<"\n";
505 f_out<<"# *"<<"\n";
506 f_out<<"# A Fast Simulator for general purpose LHC detector *"<<"\n";
507 f_out<<"# Written by S. Ovyn and X. Rouby *"<<"\n";
508 f_out<<"# severine.ovyn@uclouvain.be *"<<"\n";
509 f_out<<"# *"<<"\n";
510 f_out<<"# http: *"<<"\n";
511 f_out<<"# *"<<"\n";
512 f_out<<"# Center for Particle Physics and Phenomenology (CP3) *"<<"\n";
513 f_out<<"# Universite Catholique de Louvain (UCL) *"<<"\n";
514 f_out<<"# Louvain-la-Neuve, Belgium *"<<"\n";
515 f_out<<"# *"<<"\n";
516 f_out<<"#....................................................................*"<<"\n";
517 f_out<<"# *"<<"\n";
518 f_out<<"# This package uses: *"<<"\n";
519 f_out<<"# FastJet algorithm: Phys. Lett. B641 (2006) [hep-ph/0512210] *"<<"\n";
520 f_out<<"# Hector: JINST 2:P09005 (2007) [physics.acc-ph:0707.1198v2] *"<<"\n";
521 f_out<<"# ExRootAnalysis *"<<"\n";
522 f_out<<"# *"<<"\n";
523 f_out<<"#....................................................................*"<<"\n";
524 f_out<<"# *"<<"\n";
525 f_out<<"# This file contains all the running parameters (detector and cuts) *"<<"\n";
526 f_out<<"# necessary to reproduce the detector simulation *"<<"\n";
527 f_out<<"# *"<<"\n";
528 f_out<<"#....................................................................*"<<"\n";
529 f_out<<"#>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>"<<"\n";
530 f_out<<"* *"<<"\n";
531 f_out<<"#******************************** *"<<"\n";
532 f_out<<"# Central detector caracteristics *"<<"\n";
533 f_out<<"#******************************** *"<<"\n";
534 f_out<<"* *"<<"\n";
535 f_out << left << setw(30) <<"* Maximum tracking system: "<<""
536 << left << setw(10) <<CEN_max_tracker <<""<< right << setw(15)<<"*"<<"\n";
537 f_out << left << setw(30) <<"* Maximum central calorimeter: "<<""
538 << left << setw(10) <<CEN_max_calo_cen <<""<< right << setw(15)<<"*"<<"\n";
539 f_out << left << setw(30) <<"* Maximum forward calorimeter: "<<""
540 << left << setw(10) <<CEN_max_calo_fwd <<""<< right << setw(15)<<"*"<<"\n";
541 f_out << left << setw(30) <<"* Muon chambers coverage: "<<""
542 << left << setw(10) <<CEN_max_mu <<""<< right << setw(15)<<"*"<<"\n";
543 f_out<<"* *"<<"\n";
544 if(FLAG_vfd==1){
545 f_out<<"#********************************** *"<<"\n";
546 f_out<<"# Very forward detector switches on *"<<"\n";
547 f_out<<"#********************************** *"<<"\n";
548 f_out<<"* *"<<"\n";
549 f_out << left << setw(55) <<"* Minimum very forward calorimeter: "<<""
550 << left << setw(5) <<VFD_min_calo_vfd <<""<< right << setw(10)<<"*"<<"\n";
551 f_out << left << setw(55) <<"* Maximum very forward calorimeter: "<<""
552 << left << setw(5) <<VFD_max_calo_vfd <<""<< right << setw(10)<<"*"<<"\n";
553 f_out << left << setw(55) <<"* Minimum coverage zero_degree calorimeter "<<""
554 << left << setw(5) <<VFD_min_zdc <<""<< right << setw(10)<<"*"<<"\n";
555 f_out << left << setw(55) <<"* Distance of the ZDC to the IP, in meters: "<<""
556 << left << setw(5) <<VFD_s_zdc <<""<< right << setw(10)<<"*"<<"\n";
557 f_out << left << setw(55) <<"* Distance of the RP to the IP, in meters: "<<""
558 << left << setw(5) <<RP_220_s <<""<< right << setw(10)<<"*"<<"\n";
559 f_out << left << setw(55) <<"* Distance of the RP to the beam, in meters: "<<""
560 << left << setw(5) <<RP_220_x <<""<< right << setw(10)<<"*"<<"\n";
561 f_out << left << setw(55) <<"* Distance of the RP to the IP, in meters: "<<""
562 << left << setw(5) <<RP_420_s <<""<< right << setw(10)<<"*"<<"\n";
563 f_out << left << setw(55) <<"* Distance of the RP to the beam, in meters: "<<""
564 << left << setw(5) <<RP_420_x <<""<< right << setw(10)<<"*"<<"\n";
565 f_out << left << setw(35) <<"* Datacard for beam 1: "<<""
566 << left << setw(25) <<RP_beam1Card <<""<< right << setw(10)<<"*"<<"\n";
567 f_out << left << setw(35) <<"* Datacard for beam 2: "<<""
568 << left << setw(25) <<RP_beam2Card <<""<< right << setw(10)<<"*"<<"\n";
569 f_out << left << setw(40) <<"* Beam separation, in meters: "<<""
570 << left << setw(10) << RP_offsetEl_x <<""<< right << setw(20)<<"! not in datacard *"<<"\n";
571 f_out << left << setw(44) <<"* Distance from IP for Beam separation (m):"<<""
572 << left << setw(6) <<RP_offsetEl_s <<""<< right << setw(20)<<"! not in datacard *"<<"\n";
573
574
575
576 f_out<<"* *"<<"\n";
577 }
578 else {
579 f_out<<"#*********************************** *"<<"\n";
580 f_out<<"# Very forward detector switches off *"<<"\n";
581 f_out<<"#*********************************** *"<<"\n";
582 f_out<<"* *"<<"\n";
583 }
584 f_out<<"#************************************ *"<<"\n";
585 f_out<<"# Electromagnetic smearing parameters *"<<"\n";
586 f_out<<"#************************************ *"<<"\n";
587 f_out<<"* *"<<"\n";
588 //# \sigma/E = C + N/E + S/\sqrt{E}
589 f_out << left << setw(30) <<"* S term for central ECAL: "<<""
590 << left << setw(30) <<ELG_Scen <<""<< right << setw(10)<<"*"<<"\n";
591 f_out << left << setw(30) <<"* N term for central ECAL: "<<""
592 << left << setw(30) <<ELG_Ncen <<""<< right << setw(10)<<"*"<<"\n";
593 f_out << left << setw(30) <<"* C term for central ECAL: "<<""
594 << left << setw(30) <<ELG_Ccen <<""<< right << setw(10)<<"*"<<"\n";
595 f_out << left << setw(30) <<"* S term for forward ECAL: "<<""
596 << left << setw(30) <<ELG_Sfwd <<""<< right << setw(10)<<"*"<<"\n";
597 f_out << left << setw(30) <<"* N term for forward ECAL: "<<""
598 << left << setw(30) <<ELG_Nfwd <<""<< right << setw(10)<<"*"<<"\n";
599 f_out << left << setw(30) <<"* C term for forward ECAL: "<<""
600 << left << setw(30) <<ELG_Cfwd <<""<< right << setw(10)<<"*"<<"\n";
601 f_out<<"* *"<<"\n";
602 f_out<<"#***************************** *"<<"\n";
603 f_out<<"# Hadronic smearing parameters *"<<"\n";
604 f_out<<"#***************************** *"<<"\n";
605 f_out<<"* *"<<"\n";
606 f_out << left << setw(30) <<"* S term for central HCAL: "<<""
607 << left << setw(30) <<HAD_Shcal <<""<< right << setw(10)<<"*"<<"\n";
608 f_out << left << setw(30) <<"* N term for central HCAL: "<<""
609 << left << setw(30) <<HAD_Nhcal <<""<< right << setw(10)<<"*"<<"\n";
610 f_out << left << setw(30) <<"* C term for central HCAL: "<<""
611 << left << setw(30) <<HAD_Chcal <<""<< right << setw(10)<<"*"<<"\n";
612 f_out << left << setw(30) <<"* S term for forward HCAL: "<<""
613 << left << setw(30) <<HAD_Shf <<""<< right << setw(10)<<"*"<<"\n";
614 f_out << left << setw(30) <<"* N term for forward HCAL: "<<""
615 << left << setw(30) <<HAD_Nhf <<""<< right << setw(10)<<"*"<<"\n";
616 f_out << left << setw(30) <<"* C term for forward HCAL: "<<""
617 << left << setw(30) <<HAD_Chf <<""<< right << setw(10)<<"*"<<"\n";
618 f_out<<"* *"<<"\n";
619 f_out<<"#************************* *"<<"\n";
620 f_out<<"# Muon smearing parameters *"<<"\n";
621 f_out<<"#************************* *"<<"\n";
622 f_out<<"* *"<<"\n";
623 f_out << left << setw(55) <<"* PT resolution for muons : "<<""
624 << left << setw(5) <<MU_SmearPt <<""<< right << setw(10)<<"*"<<"\n";
625 f_out<<"* *"<<"\n";
626 if(FLAG_bfield==1){
627 f_out<<"#*************************** *"<<"\n";
628 f_out<<"# Magnetic field switches on *"<<"\n";
629 f_out<<"#*************************** *"<<"\n";
630 f_out<<"* *"<<"\n";
631 f_out << left << setw(55) <<"* Radius of the BField coverage: "<<""
632 << left << setw(5) <<TRACK_radius <<""<< right << setw(10)<<"*"<<"\n";
633 f_out << left << setw(55) <<"* Length of the BField coverage: "<<""
634 << left << setw(5) <<TRACK_length <<""<< right << setw(10)<<"*"<<"\n";
635 f_out << left << setw(55) <<"* BField X component: "<<""
636 << left << setw(5) <<TRACK_bfield_x <<""<< right << setw(10)<<"*"<<"\n";
637 f_out << left << setw(55) <<"* BField Y component: "<<""
638 << left << setw(5) <<TRACK_bfield_y <<""<< right << setw(10)<<"*"<<"\n";
639 f_out << left << setw(55) <<"* BField Z component: "<<""
640 << left << setw(5) <<TRACK_bfield_z <<""<< right << setw(10)<<"*"<<"\n";
641 f_out << left << setw(55) <<"* Minimal pT needed to reach the calorimeter [GeV]: "<<""
642 << left << setw(10) <<TRACK_ptmin <<""<< right << setw(5)<<"*"<<"\n";
643 f_out << left << setw(55) <<"* Efficiency associated to the tracking: "<<""
644 << left << setw(10) <<TRACK_eff <<""<< right << setw(5)<<"*"<<"\n";
645 f_out<<"* *"<<"\n";
646 }
647 else {
648 f_out<<"#**************************** *"<<"\n";
649 f_out<<"# Magnetic field switches off *"<<"\n";
650 f_out<<"#**************************** *"<<"\n";
651 f_out << left << setw(55) <<"* Minimal pT needed to reach the calorimeter [GeV]: "<<""
652 << left << setw(10) <<TRACK_ptmin <<""<< right << setw(5)<<"*"<<"\n";
653 f_out << left << setw(55) <<"* Efficiency associated to the tracking: "<<""
654 << left << setw(10) <<TRACK_eff <<""<< right << setw(5)<<"*"<<"\n";
655 f_out<<"* *"<<"\n";
656 }
657 f_out<<"#******************** *"<<"\n";
658 f_out<<"# Calorimetric Towers *"<<"\n";
659 f_out<<"#******************** *"<<"\n";
660 f_out << left << setw(55) <<"* Number of calorimetric towers in eta, for eta>0: "<<""
661 << left << setw(5) << TOWER_number <<""<< right << setw(10)<<"*"<<"\n";
662 f_out << left << setw(55) <<"* Tower edges in eta, for eta>0: "<<"" << right << setw(15)<<"*"<<"\n";
663 f_out << "* ";
664 for (unsigned int i=0; i<TOWER_number+1; i++) {
665 f_out << left << setw(7) << TOWER_eta_edges[i];
666 if(!( (i+1) %9 )) f_out << right << setw(3) << "*" << "\n" << "* ";
667 }
668 for (unsigned int i=(TOWER_number+1)%9; i<9; i++) f_out << left << setw(7) << "";
669 f_out << right << setw(3)<<"*"<<"\n";
670 f_out << left << setw(55) <<"* Tower sizes in phi, for eta>0 [degree]:"<<"" << right << setw(15)<<"*"<<"\n";
671 f_out << "* ";
672 for (unsigned int i=0; i<TOWER_number; i++) {
673 f_out << left << setw(7) << TOWER_dphi[i];
674 if(!( (i+1) %9 )) f_out << right << setw(3) << "*" << "\n" << "* ";
675 }
676 for (unsigned int i=(TOWER_number)%9; i<9; i++) f_out << left << setw(7) << "";
677 f_out << right << setw(3)<<"*"<<"\n";
678 f_out<<"* *"<<"\n";
679 f_out<<"#******************* *"<<"\n";
680 f_out<<"# Minimum pT's [GeV] *"<<"\n";
681 f_out<<"#******************* *"<<"\n";
682 f_out<<"* *"<<"\n";
683 f_out << left << setw(40) <<"* Minimum pT for electrons: "<<""
684 << left << setw(20) <<PTCUT_elec <<""<< right << setw(10)<<"*"<<"\n";
685 f_out << left << setw(40) <<"* Minimum pT for muons: "<<""
686 << left << setw(20) <<PTCUT_muon <<""<< right << setw(10)<<"*"<<"\n";
687 f_out << left << setw(40) <<"* Minimum pT for jets: "<<""
688 << left << setw(20) <<PTCUT_jet <<""<< right << setw(10)<<"*"<<"\n";
689 f_out << left << setw(40) <<"* Minimum pT for Tau-jets: "<<""
690 << left << setw(20) <<PTCUT_taujet <<""<< right << setw(10)<<"*"<<"\n";
691 f_out << left << setw(40) <<"* Minimum pT for photons: "<<""
692 << left << setw(20) <<PTCUT_gamma <<""<< right << setw(10)<<"*"<<"\n";
693 f_out<<"* *"<<"\n";
694 f_out<<"#*************** *"<<"\n";
695 f_out<<"# Jet definition *"<<"\n";
696 f_out<<"#*************** *"<<"\n";
697 f_out<<"* *"<<"\n";
698 f_out<<"* Six algorithms are currently available: *"<<"\n";
699 f_out<<"* - 1) CDF cone algorithm, *"<<"\n";
700 f_out<<"* - 2) CDF MidPoint algorithm, *"<<"\n";
701 f_out<<"* - 3) SIScone algorithm, *"<<"\n";
702 f_out<<"* - 4) kt algorithm, *"<<"\n";
703 f_out<<"* - 5) Cambrigde/Aachen algorithm, *"<<"\n";
704 f_out<<"* - 6) Anti-kt algorithm. *"<<"\n";
705 f_out<<"* *"<<"\n";
706 f_out<<"* You have chosen *"<<"\n";
707 switch(JET_jetalgo) {
708 default:
709 case 1: {
710 f_out<<"* CDF JetClu jet algorithm with parameters: *"<<"\n";
711 f_out << left << setw(40) <<"* - Seed threshold: "<<""
712 << left << setw(10) <<JET_seed <<""<< right << setw(20)<<"! not in datacard *"<<"\n";
713 f_out << left << setw(40) <<"* - Cone radius: "<<""
714 << left << setw(10) <<JET_coneradius <<""<< right << setw(20)<<"*"<<"\n";
715 f_out << left << setw(40) <<"* - Adjacency cut: "<<""
716 << left << setw(10) <<JET_C_adjacencycut <<""<< right << setw(20)<<"! not in datacard *"<<"\n";
717 f_out << left << setw(40) <<"* - Max iterations: "<<""
718 << left << setw(10) <<JET_C_maxiterations <<""<< right << setw(20)<<"! not in datacard *"<<"\n";
719 f_out << left << setw(40) <<"* - Iratch: "<<""
720 << left << setw(10) <<JET_C_iratch <<""<< right << setw(20)<<"! not in datacard *"<<"\n";
721 f_out << left << setw(40) <<"* - Overlap threshold: "<<""
722 << left << setw(10) <<JET_overlap <<""<< right << setw(20)<<"! not in datacard *"<<"\n";
723 }
724 break;
725 case 2: {
726 f_out<<"* CDF midpoint jet algorithm with parameters: *"<<"\n";
727 f_out << left << setw(40) <<"* - Seed threshold: "<<""
728 << left << setw(20) <<JET_seed <<""<< right << setw(10)<<"! not in datacard *"<<"\n";
729 f_out << left << setw(40) <<"* - Cone radius: "<<""
730 << left << setw(20) <<JET_coneradius <<""<< right << setw(10)<<"*"<<"\n";
731 f_out << left << setw(40) <<"* - Cone area fraction:"<<""
732 << left << setw(20) <<JET_M_coneareafraction <<""<< right << setw(10)<<"! not in datacard *"<<"\n";
733 f_out << left << setw(40) <<"* - Maximum pair size: "<<""
734 << left << setw(20) <<JET_M_maxpairsize <<""<< right << setw(10)<<"! not in datacard *"<<"\n";
735 f_out << left << setw(40) <<"* - Max iterations: "<<""
736 << left << setw(20) <<JET_M_maxiterations <<""<< right << setw(10)<<"! not in datacard *"<<"\n";
737 f_out << left << setw(40) <<"* - Overlap threshold: "<<""
738 << left << setw(20) <<JET_overlap <<""<< right << setw(10)<<"! not in datacard *"<<"\n";
739 }
740 break;
741 case 3: {
742 f_out <<"* SISCone jet algorithm with parameters: *"<<"\n";
743 f_out << left << setw(40) <<"* - Cone radius: "<<""
744 << left << setw(20) <<JET_coneradius <<""<< right << setw(10)<<"*"<<"\n";
745 f_out << left << setw(40) <<"* - Overlap threshold: "<<""
746 << left << setw(20) <<JET_overlap <<""<< right << setw(10)<<"! not in datacard *"<<"\n";
747 f_out << left << setw(40) <<"* - Number pass max: "<<""
748 << left << setw(20) <<JET_S_npass <<""<< right << setw(10)<<"! not in datacard *"<<"\n";
749 f_out << left << setw(40) <<"* - Minimum pT for protojet: "<<""
750 << left << setw(20) <<JET_S_protojet_ptmin <<""<< right << setw(10)<<"! not in datacard *"<<"\n";
751 }
752 break;
753 case 4: {
754 f_out <<"* KT jet algorithm with parameters: *"<<"\n";
755 f_out << left << setw(40) <<"* - Cone radius: "<<""
756 << left << setw(20) <<JET_coneradius <<""<< right << setw(10)<<"*"<<"\n";
757 }
758 break;
759 case 5: {
760 f_out <<"* Cambridge/Aachen jet algorithm with parameters: *"<<"\n";
761 f_out << left << setw(40) <<"* - Cone radius: "<<""
762 << left << setw(20) <<JET_coneradius <<""<< right << setw(10)<<"*"<<"\n";
763 }
764 break;
765 case 6: {
766 f_out <<"* Anti-kt jet algorithm with parameters: *"<<"\n";
767 f_out << left << setw(40) <<"* - Cone radius: "<<""
768 << left << setw(20) <<JET_coneradius <<""<< right << setw(10)<<"*"<<"\n";
769 }
770 break;
771 }
772 f_out<<"* *"<<"\n";
773 f_out<<"#****************************** *"<<"\n";
774 f_out<<"# Tau-jet definition parameters *"<<"\n";
775 f_out<<"#****************************** *"<<"\n";
776 f_out<<"* *"<<"\n";
777 f_out << left << setw(45) <<"* Cone radius for calorimeter tagging: "<<""
778 << left << setw(5) <<TAU_energy_scone <<""<< right << setw(20)<<"*"<<"\n";
779 f_out << left << setw(45) <<"* Fraction of energy in the small cone: "<<""
780 << left << setw(5) <<TAU_energy_frac*100 <<""<< right << setw(20)<<"! not in datacard *"<<"\n";
781 f_out << left << setw(45) <<"* Cone radius for tracking tagging: "<<""
782 << left << setw(5) <<TAU_track_scone <<""<< right << setw(20)<<"*"<<"\n";
783 f_out << left << setw(45) <<"* Minimum track pT [GeV]: "<<""
784 << left << setw(5) <<TAU_track_pt <<""<< right << setw(20)<<"*"<<"\n";
785 f_out<<"* *"<<"\n";
786 f_out<<"#*************************** *"<<"\n";
787 f_out<<"# B-tagging efficiencies [%] *"<<"\n";
788 f_out<<"#*************************** *"<<"\n";
789 f_out<<"* *"<<"\n";
790 f_out << left << setw(50) <<"* Efficiency to tag a \"b\" as a b-jet: "<<""
791 << left << setw(10) <<BTAG_b <<""<< right << setw(10)<<"*"<<"\n";
792 f_out << left << setw(50) <<"* Efficiency to mistag a c-jet as a b-jet: "<<""
793 << left << setw(10) <<BTAG_mistag_c <<""<< right << setw(10)<<"*"<<"\n";
794 f_out << left << setw(50) <<"* Efficiency to mistag a light jet as a b-jet: "<<""
795 << left << setw(10) <<BTAG_mistag_l <<""<< right << setw(10)<<"*"<<"\n";
796 f_out<<"* *"<<"\n";
797 f_out<<"* *"<<"\n";
798 f_out<<"#....................................................................*"<<"\n";
799 f_out<<"#>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>"<<"\n";
800
801}
802
803// **********Provides the smeared TLorentzVector for the electrons********
804// Smears the electron energy, and changes the 4-momentum accordingly
805// different smearing if the electron is central (eta < 2.5) or forward
806void RESOLution::SmearElectron(TLorentzVector &electron) {
807 // the 'electron' variable will be changed by the function
808 float energy = electron.E(); // before smearing
809 float energyS = 0.0; // after smearing // \sigma/E = C + N/E + S/\sqrt{E}
810
811 if(fabs(electron.Eta()) < CEN_max_tracker) { // if the electron is inside the tracker
812 energyS = gRandom->Gaus(energy, sqrt(
813 pow(ELG_Ncen,2) +
814 pow(ELG_Ccen*energy,2) +
815 pow(ELG_Scen*sqrt(energy),2) ));
816 }
817 if(fabs(electron.Eta()) > CEN_max_tracker && fabs(electron.Eta()) < CEN_max_calo_fwd){
818 energyS = gRandom->Gaus(energy, sqrt(
819 pow(ELG_Nfwd,2) +
820 pow(ELG_Cfwd*energy,2) +
821 pow(ELG_Sfwd*sqrt(energy),2) ) );
822 }
823 electron.SetPtEtaPhiE(energyS/cosh(electron.Eta()), electron.Eta(), electron.Phi(), energyS);
824 if(electron.E() < 0)electron.SetPxPyPzE(0,0,0,0); // no negative values after smearing !
825}
826
827
828// **********Provides the smeared TLorentzVector for the muons********
829// Smears the muon pT and changes the 4-momentum accordingly
830void RESOLution::SmearMu(TLorentzVector &muon) {
831 // the 'muon' variable will be changed by the function
832 float pt = muon.Pt(); // before smearing
833 float ptS=pt;
834
835 if(fabs(muon.Eta()) < CEN_max_mu )
836 {
837 ptS = gRandom->Gaus(pt, MU_SmearPt*pt ); // after smearing // \sigma/E = C + N/E + S/\sqrt{E}
838 }
839 muon.SetPtEtaPhiE(ptS, muon.Eta(), muon.Phi(), ptS*cosh(muon.Eta()));
840
841 if(muon.E() < 0)muon.SetPxPyPzE(0,0,0,0); // no negative values after smearing !
842}
843
844
845// **********Provides the smeared TLorentzVector for the hadrons********
846// Smears the hadron 4-momentum
847void RESOLution::SmearHadron(TLorentzVector &hadron, const float frac)
848 // the 'hadron' variable will be changed by the function
849 // the 'frac' variable describes the long-living particles. Should be 0.7 for K0S and Lambda, 1. otherwise
850{
851 float energy = hadron.E(); // before smearing
852 float energyS = 0.0; // after smearing // \sigma/E = C + N/E + S/\sqrt{E}
853 float energy_ecal = (1.0 - frac)*energy; // electromagnetic calorimeter
854 float energy_hcal = frac*energy; // hadronic calorimeter
855 // frac takes into account the decay of long-living particles, that decay in the calorimeters
856 // some of the particles decay mostly in the ecal, some mostly in the hcal
857
858 float energyS1,energyS2;
859 if(fabs(hadron.Eta()) < CEN_max_calo_cen) {
860 energyS1 = gRandom->Gaus(energy_hcal, sqrt(
861 pow(HAD_Nhcal,2) +
862 pow(HAD_Chcal*energy_hcal,2) +
863 pow(HAD_Shcal*sqrt(energy_hcal),2) )) ;
864
865
866 energyS2 = gRandom->Gaus(energy_ecal, sqrt(
867 pow(ELG_Ncen,2) +
868 pow(ELG_Ccen*energy_ecal,2) +
869 pow(ELG_Scen*sqrt(energy_ecal),2) ) );
870
871 energyS = ((energyS1>0)?energyS1:0) + ((energyS2>0)?energyS2:0);
872 }
873 if(fabs(hadron.Eta()) > CEN_max_calo_cen && fabs(hadron.Eta()) < CEN_max_calo_fwd){
874 energyS = gRandom->Gaus(energy, sqrt(
875 pow(HAD_Nhf,2) +
876 pow(HAD_Chf*energy,2) +
877 pow(HAD_Shf*sqrt(energy),2) ));
878}
879
880
881
882 hadron.SetPtEtaPhiE(energyS/cosh(hadron.Eta()),hadron.Eta(), hadron.Phi(), energyS);
883
884 if(hadron.E() < 0)hadron.SetPxPyPzE(0,0,0,0);
885}
886
887//******************************************************************************************
888
889void RESOLution::SortedVector(vector<ParticleUtil> &vect)
890{
891 int i,j = 0;
892 TLorentzVector tmp;
893 bool en_desordre = true;
894 int entries=vect.size();
895 for(i = 0 ; (i < entries) && en_desordre; i++)
896 {
897 en_desordre = false;
898 for(j = 1 ; j < entries - i ; j++)
899 {
900 if ( vect[j].Pt() > vect[j-1].Pt() )
901 {
902 ParticleUtil tmp = vect[j-1];
903 vect[j-1] = vect[j];
904 vect[j] = tmp;
905 en_desordre = true;
906 }
907 }
908 }
909}
910
911// **********Provides the energy in the cone of radius TAU_CONE_ENERGY for the tau identification********
912// to be taken into account, a calo tower should
913// 1) have a transverse energy \f$ E_T = \sqrt{E_X^2 + E_Y^2} \f$ above a given threshold
914// 2) be inside a cone with a radius R and the axis defined by (eta,phi)
915double RESOLution::EnergySmallCone(const vector<PhysicsTower> &towers, const float eta, const float phi) {
916 double Energie=0;
917 for(unsigned int i=0; i < towers.size(); i++) {
918 if(towers[i].fourVector.pt() < JET_seed) continue;
919 if((DeltaR(phi,eta,towers[i].fourVector.phi(),towers[i].fourVector.eta()) < TAU_energy_scone)) {
920 Energie += towers[i].fourVector.E;
921 }
922 }
923 return Energie;
924}
925
926
927// **********Provides the number of tracks in the cone of radius TAU_CONE_TRACKS for the tau identification********
928// to be taken into account, a track should
929// 1) avec a transverse momentum \$f p_T \$ above a given threshold
930// 2) be inside a cone with a radius R and the axis defined by (eta,phi)
931// IMPORTANT REMARK !!!!!
932// previously, the argument 'phi' was before the argument 'eta'
933// this has been changed for consistency with the other functions
934// double check your running code that uses NumTracks !
935unsigned int RESOLution::NumTracks(const vector<TLorentzVector> &tracks, const float pt_track, const float eta, const float phi) {
936 unsigned int numtrack=0;
937 for(unsigned int i=0; i < tracks.size(); i++) {
938 if((tracks[i].Pt() < pt_track )||
939 (DeltaR(phi,eta,tracks[i].Phi(),tracks[i].Eta()) > TAU_track_scone)
940 )continue;
941 numtrack++;
942 }
943 return numtrack;
944}
945
946
947//*** Returns the PID of the particle with the highest energy, in a cone with a radius CONERADIUS and an axis (eta,phi) *********
948//used by Btaggedjet
949///// Attention : bug removed => CONERADIUS/2 -> CONERADIUS !!
950int RESOLution::Bjets(const TSimpleArray<TRootGenParticle> &subarray, const float eta, const float phi) {
951 float emax=0;
952 int Ppid=0;
953 if(subarray.GetEntries()>0) {
954 for(int i=0; i < subarray.GetEntries();i++) { // should have pt>PT_JETMIN and a small cone radius (r<CONE_JET)
955 float genDeltaR = DeltaR(subarray[i]->Phi,subarray[i]->Eta,phi,eta);
956 if(genDeltaR < JET_coneradius && subarray[i]->E > emax) {
957 emax=subarray[i]->E;
958 Ppid=abs(subarray[i]->PID);
959 }
960 }
961 }
962 return Ppid;
963}
964
965
966//******************** Simulates the b-tagging efficiency for real bjet, or the misendentification for other jets****************
967bool RESOLution::Btaggedjet(const TLorentzVector &JET, const TSimpleArray<TRootGenParticle> &subarray) {
968 if( rand()%100 < (BTAG_b+1) && Bjets(subarray,JET.Eta(),JET.Phi())==pB ) return true; // b-tag of b-jets is 40%
969 else if( rand()%100 < (BTAG_mistag_c+1) && Bjets(subarray,JET.Eta(),JET.Phi())==pC ) return true; // b-tag of c-jets is 10%
970 else if( rand()%100 < (BTAG_mistag_l+1) && Bjets(subarray,JET.Eta(),JET.Phi())!=0) return true; // b-tag of light jets is 1%
971 return false;
972}
973
974//***********************Isolation criteria***********************
975//****************************************************************
976bool RESOLution::Isolation(const float phi, const float eta,const vector<TLorentzVector> &tracks, const float pt_second_track)
977{
978 bool isolated = false;
979 float deltar=5000.; // Initial value; should be high; no further repercussion
980 // loop on all final charged particles, with p_t >2, close enough from the electron
981 for(unsigned int i=0; i < tracks.size(); i++)
982 {
983 if(tracks[i].Pt() < pt_second_track)continue;
984 float genDeltaR = DeltaR(phi,eta,tracks[i].Phi(),tracks[i].Eta());
985 if(
986 (genDeltaR > deltar) ||
987 (genDeltaR==0)
988 ) continue ;
989 deltar=genDeltaR;
990 }
991 if(deltar > 0.5) isolated = true;
992 return isolated;
993}
994
995
996 //********** returns a segmented value for eta and phi, for calo towers *****
997void RESOLution::BinEtaPhi(const float phi, const float eta, float& iPhi, float& iEta){
998 iEta = -100;
999 int index=-100;
1000 for (unsigned int i=1; i< TOWER_number+1; i++) {
1001 if(fabs(eta)>TOWER_eta_edges[i-1] && fabs(eta)<TOWER_eta_edges[i]) {
1002 iEta = (eta>0) ? TOWER_eta_edges[i-1] : -TOWER_eta_edges[i];
1003 index = i-1;
1004 //cout << setw(15) << left << eta << "\t" << iEta << endl;
1005 break;
1006 }
1007 }
1008 if(index==-100) return;
1009 iPhi = -100;
1010 float dphi = TOWER_dphi[index]*pi/180.;
1011 for (unsigned int i=1; i < 360/TOWER_dphi[index]; i++ ) {
1012 float low = -pi+(i-1)*dphi;
1013 float high= low+dphi;
1014 if(phi > low && phi < high ){
1015 iPhi = low;
1016 break;
1017 }
1018 }
1019 if (phi > pi-dphi) iPhi = pi-dphi;
1020}
1021
1022//**************************** Returns the delta Phi ****************************
1023float DeltaPhi(const float phi1, const float phi2) {
1024 float deltaphi=phi1-phi2; // in here, -pi < phi < pi
1025 if(fabs(deltaphi) > pi) {
1026 deltaphi=2.*pi -fabs(deltaphi);// put deltaphi between 0 and pi
1027 }
1028 else deltaphi=fabs(deltaphi);
1029
1030 return deltaphi;
1031}
1032
1033//**************************** Returns the delta R****************************
1034float DeltaR(const float phi1, const float eta1, const float phi2, const float eta2) {
1035 return sqrt(pow(DeltaPhi(phi1,phi2),2) + pow(eta1-eta2,2));
1036}
1037
1038int sign(const int myint) {
1039 if (myint >0) return 1;
1040 else if (myint <0) return -1;
1041 else return 0;
1042}
1043
1044int sign(const float myfloat) {
1045 if (myfloat >0) return 1;
1046 else if (myfloat <0) return -1;
1047 else return 0;
1048}
1049
1050int Charge(const int pid)
1051{
1052 int charge;
1053 if(
1054 (pid == pGAMMA) ||
1055 (pid == pPI0) ||
1056 (pid == pK0L) ||
1057 (pid == pN) ||
1058 (pid == pSIGMA0) ||
1059 (pid == pDELTA0) ||
1060 (pid == pK0S) // not charged particles : invisible by tracker
1061 )
1062 charge = 0;
1063 else charge = (sign(pid));
1064 return charge;
1065
1066}
Note: See TracBrowser for help on using the repository browser.