Fork me on GitHub

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

Last change on this file since 261 was 261, checked in by severine ovyn, 15 years ago

small changes

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