Fork me on GitHub

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

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

include statements have been cleaned ; copy-constructor ; operator= ; destructor

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