Fork me on GitHub

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

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

add shifts X,Y,angl

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