Fork me on GitHub

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

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

moving from PI to pi

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