Fork me on GitHub

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

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

add parameters for RP

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