Fork me on GitHub

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

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

small changes

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