Fork me on GitHub

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

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

add shifts X,Y,angl

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