Fork me on GitHub

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

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

OK flags very forward

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