Fork me on GitHub

source: svn/trunk/interface/SmearUtil.h@ 289

Last change on this file since 289 was 287, checked in by Xavier Rouby, 16 years ago

charge added in Taujets and tracks

File size: 10.4 KB
Line 
1#ifndef _SMEARUTIL_H_
2#define _SMEARUTIL_H_
3
4/***********************************************************************
5** **
6** /----------------------------------------------\ **
7** | Delphes, a framework for the fast simulation | **
8** | of a generic collider experiment | **
9** \----------------------------------------------/ **
10** **
11** **
12** This package uses: **
13** ------------------ **
14** FastJet algorithm: Phys. Lett. B641 (2006) [hep-ph/0512210] **
15** Hector: JINST 2:P09005 (2007) [physics.acc-ph:0707.1198v2] **
16** FROG: [hep-ex/0901.2718v1] **
17** **
18** ------------------------------------------------------------------ **
19** **
20** Main authors: **
21** ------------- **
22** **
23** Severine Ovyn Xavier Rouby **
24** severine.ovyn@uclouvain.be xavier.rouby@cern **
25** **
26** Center for Particle Physics and Phenomenology (CP3) **
27** Universite catholique de Louvain (UCL) **
28** Louvain-la-Neuve, Belgium **
29** **
30** Copyright (C) 2008-2009, **
31** All rights reserved. **
32** **
33***********************************************************************/
34
35/// \file SmearUtil.h
36/// \brief RESOLution class, and some generic definitions
37
38
39#include <vector>
40#include "TLorentzVector.h"
41
42#include "D_Constants.h"
43#include "CaloUtil.h"
44#include "BlockClasses.h"
45#include "TSimpleArray.h"
46#include "PhysicsTower.hh"
47
48using namespace std;
49
50class D_Particle {
51
52 public:
53 D_Particle(const TLorentzVector & p, const int pid, const float etacalo, const float phicalo) :
54 _fourmomentum(p), _pid(pid), _etaCalo(etacalo), _phiCalo(phicalo) {}
55 //D_Particle(const float e, const float eta, const float phi, const float pt, const int pid) :
56 // _pid(pid), _etaCalo(UNDEFINED), _phiCalo(UNDEFINED) { TLorentzVector p; p.SetPtEtaPhiE(pt,eta,phi,e); _fourmomentum = p; }
57 D_Particle(const float px, const float py, const float pz, const float e, const int pid) :
58 _fourmomentum(px,py,pz,e), _pid(pid), _etaCalo(UNDEFINED), _phiCalo(UNDEFINED) {}
59
60 const float E() const {return _fourmomentum.E();} // particle energy [GeV]
61 const float Px() const {return _fourmomentum.Px();} // horizontal coordinate of momentum [GeV]
62 const float Py() const {return _fourmomentum.Py();} // vertical coordinate of momentum [GeV]
63 const float Pz() const {return _fourmomentum.Pz();} // longitudinal coordinate of momentum [GeV]
64 const float Pt() const {return _fourmomentum.Pt();} // transverse momentum [GeV]
65 const float EtaCalo() const {return _etaCalo;} // pseudorapidity
66 const float Eta() const {return _fourmomentum.Eta();} // pseudorapidity
67 const float PhiCalo() const {return _phiCalo;} // azimuthal angle
68 const float Phi() const {return _fourmomentum.Phi();} // azimuthal angle
69 const int PID() const {return _pid;} // particle energy in [GeV]
70 const TLorentzVector& getFourMomentum() const {return _fourmomentum;}
71
72 private:
73 TLorentzVector _fourmomentum;
74 int _pid;
75 float _etaCalo, _phiCalo;
76};
77
78class RESOLution
79{
80 public:
81 /// Constructor
82 RESOLution();
83 RESOLution(const RESOLution & DET);
84 RESOLution& operator=(const RESOLution& DET);
85 ~RESOLution() { delete [] TOWER_eta_edges; delete [] TOWER_dphi;};
86
87 // Detector coverage
88 float CEN_max_tracker; // tracker pseudorapidity coverage
89 float CEN_max_calo_cen; // central calorimeter pseudorapidity coverage
90 float CEN_max_calo_fwd; // forward calorimeter pseudorapidity coverage
91 float CEN_max_mu; // muon chambers pseudorapidity coverage
92
93 float VFD_min_calo_vfd; // very forward calorimeter pseudorapidity coverage
94 float VFD_max_calo_vfd; // very forward calorimeter pseudorapidity coverage
95 float VFD_min_zdc; // coverage for Zero Degree Calorimeter, for photons and neutrons
96 float VFD_s_zdc; // distance of the Zero Degree Calorimeter, from the Interaction poin, in [m]
97
98 float RP_220_s; // distance of the RP to the IP, in meters
99 float RP_220_x; // distance of the RP to the beam, in meters
100 float RP_420_s; // distance of the RP to the IP, in meters
101 float RP_420_x; // distance of the RP to the beam, in meters
102 string RP_beam1Card; //
103 string RP_beam2Card; //
104 string RP_IP_name; //
105 float RP_offsetEl_s;
106 float RP_offsetEl_x;
107 float RP_cross_x;
108 float RP_cross_y;
109 float RP_cross_ang;
110
111
112 //energy resolution for electron/photon
113 // \sigma/E = C + N/E + S/\sqrt{E}
114 float ELG_Scen; // S term for central ECAL
115 float ELG_Ncen; // N term for central ECAL
116 float ELG_Ccen; // C term for central ECAL
117 float ELG_Sfwd; // S term for forward ECAL
118 float ELG_Cfwd; // C term for forward ECAL
119 float ELG_Nfwd; // N term for central ECAL
120
121 //energy resolution for hadrons in ecal/hcal/hf
122 // \sigma/E = C + N/E + S/\sqrt{E}
123 float HAD_Shcal; // S term for central HCAL // hadronic calorimeter
124 float HAD_Nhcal; // N term for central HCAL
125 float HAD_Chcal; // C term for central HCAL
126 float HAD_Shf; // S term for central HF // forward calorimeter
127 float HAD_Nhf; // N term for central HF
128 float HAD_Chf; // C term for central HF
129
130 // muon smearing
131 float MU_SmearPt;
132
133 //Magnetic Field information
134 int TRACK_radius; //radius of the BField coverage
135 int TRACK_length; //length of the BField coverage
136 float TRACK_bfield_x;
137 float TRACK_bfield_y;
138 float TRACK_bfield_z;
139 float TRACK_ptmin; // minimal pt needed to reach the calorimeter, in GeV
140 int TRACK_eff; // in percent, should be an integer
141
142
143 //Define Calorimetric towers
144 unsigned int TOWER_number;
145 float * TOWER_eta_edges;
146 float * TOWER_dphi;
147
148 //thresholds for reconstructed objetcs
149 float PTCUT_elec;
150 float PTCUT_muon;
151 float PTCUT_jet;
152 float PTCUT_gamma;
153 float PTCUT_taujet;
154
155 //General jet variable
156 double JET_coneradius;
157 int JET_jetalgo;
158 double JET_seed;
159 double JET_overlap;
160
161 // MidPoint algorithm definition
162 double JET_M_coneareafraction;
163 int JET_M_maxpairsize;
164 int JET_M_maxiterations;
165 // Define Cone algorithm.
166 int JET_C_adjacencycut;
167 int JET_C_maxiterations;
168 int JET_C_iratch;
169 //Define SISCone algorithm.
170 int JET_S_npass;
171 double JET_S_protojet_ptmin;
172
173 //For Tau-jet definition
174 // R = sqrt (phi^2 + eta^2)
175 float TAU_energy_scone; // radius R of the cone for tau definition, based on energy threshold
176 float TAU_track_scone; // radius R of the cone for tau definition, based on track number
177 float TAU_track_pt; // minimal pt [GeV] for tracks to be considered in tau definition
178 float TAU_energy_frac; // fraction of energy required in the central part of the cone, for tau jets
179
180 //tagging definition
181 int BTAG_b;
182 int BTAG_mistag_c;
183 int BTAG_mistag_l;
184
185
186 //trigger flag
187 int FLAG_trigger; //flag for trigger
188 int FLAG_frog; //flag for frog display
189 int FLAG_bfield; //flag for bfield propagation
190 int FLAG_vfd; //flag for very forward detector
191 int FLAG_zdc; //flag for very forward detector
192
193 int NEvents_Frog;
194 float PT_QUARKS_MIN; // minimal pt needed for quarks to reach the tracker, in GeV
195
196 // to sort a vector
197 //void SortedVector(vector<ParticleUtil> &vect);
198 void SortedVector(vector<D_Particle> &vect);
199
200 /// Reads the data card for the initialisation of the parameters
201 void ReadDataCard(const string datacard);
202
203 /// Create the output log file
204 void Logfile(const string& LogName);
205
206 /// Provides the smeared TLorentzVector for the electrons
207 void SmearElectron(TLorentzVector &electron);
208
209 /// Provides the smeared TLorentzVector for the muons
210 void SmearMu(TLorentzVector &muon);
211
212 /// Provides the smeared TLorentzVector for the hadrons
213 void SmearHadron(TLorentzVector &hadron, const float frac);
214
215 /// For electromagnetic collimation in tau jets
216 double EnergySmallCone(const vector<PhysicsTower> &towers, const float eta, const float phi);
217
218 /// Number of tracks in tau jet algo
219 //unsigned int NumTracks(float& charge, const vector<TLorentzVector> &tracks, const float pt_track, const float eta, const float phi);
220 unsigned int NumTracks(float& charge, const vector<TRootTracks> &tracks, const float pt_track, const float eta, const float phi);
221
222 /// b-jets
223 int Bjets(const TSimpleArray<TRootGenParticle> &subarray, const float eta, const float phi);
224
225 /// b-tag efficiency and misidentification
226 bool Btaggedjet(const TLorentzVector &JET, const TSimpleArray<TRootGenParticle> &subarray);
227
228 /// Lepton isolation
229 //bool Isolation(const float phi, const float eta,const vector<TLorentzVector> &tracks,float PT_TRACK2);
230 bool Isolation(const float phi, const float eta,const vector<TRootTracks> &tracks, const float PT_TRACK2);
231
232 //********************* returns a segmented value for eta and phi, for calo towers *****
233 void BinEtaPhi(const float phi, const float eta, float& iPhi, float& iEta);
234
235};
236
237// ** returns the sign (+1 or -1) or an integer
238int sign(const int myint);
239int sign(const float myfloat);
240
241// **************************** Return the Delta Phi****************************
242float DeltaPhi(const float phi1, const float phi2);
243
244// **************************** Returns the Delta R****************************
245float DeltaR(const float phi1, const float eta1, const float phi2, const float eta2);
246
247//************* Returns an array of the quarks sitting within the tracker acceptance ***************
248int ChargeVal(const int pid);
249
250#endif
Note: See TracBrowser for help on using the repository browser.