Fork me on GitHub

source: svn/trunk/Utilities/ExRootAnalysis/interface/BlockClasses.h@ 268

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

nettoyage

File size: 13.5 KB
Line 
1#ifndef BLOCKCLASSES_H
2#define BLOCKCLASSES_H
3
4/***********************************************************************
5** **
6** /----------------------------------------------\ **
7** | Delphes, a framework for the fast simulation | **
8** | of a generic collider experiment | **
9** \----------------------------------------------/ **
10** **
11** **
12** This package uses: **
13** ------------------ **
14** FastJet algorithm: Phys. Lett. B641 (2006) [hep-ph/0512210] **
15** Hector: JINST 2:P09005 (2007) [physics.acc-ph:0707.1198v2] **
16** FROG: [hep-ex/0901.2718v1] **
17** **
18** ------------------------------------------------------------------ **
19** **
20** Main authors: **
21** ------------- **
22** **
23** Severine Ovyn Xavier Rouby **
24** severine.ovyn@uclouvain.be xavier.rouby@cern **
25** **
26** Center for Particle Physics and Phenomenology (CP3) **
27** Universite catholique de Louvain (UCL) **
28** Louvain-la-Neuve, Belgium **
29** **
30** Copyright (C) 2008-2009, **
31** All rights reserved. **
32** **
33***********************************************************************/
34
35#include "TLorentzVector.h"
36#include "TObject.h"
37#include "BlockCompare.h"
38#include "interface/D_Constants.h"
39#include "interface/CaloUtil.h"
40
41class TSortableObject: public TObject
42{
43public:
44 TSortableObject() {};
45 Bool_t IsSortable() const { return GetCompare() ? GetCompare()->IsSortable(this) : kFALSE; }
46 Int_t Compare(const TObject *obj) const { return GetCompare()->Compare(this, obj); }
47
48 virtual const TCompare *GetCompare() const = 0;
49
50 ClassDef(TSortableObject, 1)
51};
52
53//---------------------------------------------------------------------------
54//
55class TRootLHEFEvent: public TObject
56{
57public:
58 TRootLHEFEvent() {};
59
60 Long64_t Number; // event number
61
62 int Nparticles; // number of particles in the event | hepup.NUP
63 int ProcessID; // subprocess code for the event | hepup.IDPRUP
64
65 Double_t Weight; // weight for the event | hepup.XWGTUP
66 Double_t ScalePDF; // scale in GeV used in the calculation of the PDFs in the event | hepup.SCALUP
67 Double_t CouplingQED; // value of the QED coupling used in the event | hepup.AQEDUP
68 Double_t CouplingQCD; // value of the QCD coupling used in the event | hepup.AQCDUP
69
70 ClassDef(TRootLHEFEvent, 2)
71};
72
73//---------------------------------------------------------------------------
74
75class TRootLHEFParticle: public TSortableObject
76{
77public:
78 TRootLHEFParticle() {};
79 int PID; // particle HEP ID number | hepup.IDUP[number]
80 int Status; // particle status code | hepup.ISTUP[number]
81 int Mother1; // index for the particle first mother | hepup.MOTHUP[number][0]
82 int Mother2; // index for the particle last mother | hepup.MOTHUP[number][1]
83 int ColorLine1; // index for the particle color-line | hepup.ICOLUP[number][0]
84 int ColorLine2; // index for the particle anti-color-line | hepup.ICOLUP[number][1]
85
86 double Px; // particle momentum vector (x component) | hepup.PUP[number][0]
87 double Py; // particle momentum vector (y component) | hepup.PUP[number][1]
88 double Pz; // particle momentum vector (z component) | hepup.PUP[number][2]
89 double E; // particle energy | hepup.PUP[number][3]
90 double M; // particle mass | hepup.PUP[number][4]
91
92 double PT; // particle transverse momentum
93 double Eta; // particle pseudorapidity
94 double Phi; // particle azimuthal angle
95
96 double Rapidity; // particle rapidity
97
98 double LifeTime; // particle invariant lifetime
99 // (c*tau, distance from production to decay in mm)
100 // | hepup.VTIMUP[number]
101
102 double Spin; // cosine of the angle between the particle spin vector
103 // and the decaying particle 3-momentum,
104 // specified in the lab frame. | hepup.SPINUP[number]
105
106 static TCompare *fgCompare; //!
107 const TCompare *GetCompare() const { return fgCompare; }
108 ClassDef(TRootLHEFParticle, 2)
109
110};
111
112//---------------------------------------------------------------------------
113
114class TRootSelectorInfo: public TObject
115{
116public:
117 TRootSelectorInfo() {};
118 int Processed; // current number of processed events
119 int Accepted; // current number of accepted events
120
121 ClassDef(TRootSelectorInfo, 1)
122};
123
124
125class TRootGenEvent: public TObject
126{
127public:
128 TRootGenEvent() {};
129 Long64_t Number; // event number | hepevt.nevhep
130
131 static TCompare *fgCompare; //!
132 const TCompare *GetCompare() const { return fgCompare; }
133
134 ClassDef(TRootGenEvent, 1)
135};
136
137
138class TRootEvent: public TObject {
139
140public:
141 TRootEvent() {};
142 int Run; // run number [G3EventProxy::simSignal().id().runNumber()]
143 int Event; // event number [G3EventProxy::simSignal().id().eventInRun()]
144
145// Short_t L1Decision; // L1 trigger global decision [L1Trigger::decision()]
146// Short_t HLTDecision; // HLT trigger global decision [HighLevelTriggerResult::getGlobalDecision()]
147
148 ClassDef(TRootEvent, 1)
149};
150
151//---------------------------------------------------------------------------
152
153class TRootParticle: public TSortableObject {
154
155public:
156
157 TRootParticle() {};
158 float E; // particle energy in GeV
159 float Px; // particle momentum vector (x component) in GeV
160 float Py; // particle momentum vector (y component) in GeV
161 float Pz; // particle momentum vector (z component) in GeV
162
163 float Eta; // particle pseudorapidity
164 float Phi; // particle azimuthal angle in rad
165
166 float EtaCalo; // particle pseudorapidity when entering the calo,
167 float PhiCalo; // particle azimuthal angle in rad when entering the calo
168
169 void Set(const TLorentzVector& momentum);
170 void Set(const float px, const float py, const float pz, const float e);
171 void SetEtaPhi(const float eta, const float phi) {Eta=eta; Phi=phi;};
172 void SetEtaPhiCalo(const float eta, const float phi) {EtaCalo=eta; PhiCalo=phi;};
173 void SetEtaPhiEET(const float eta, const float phi, const float e, const float et);
174 static TCompare *fgCompare; //!
175 const TCompare *GetCompare() const { return fgCompare; }
176 float PT; // particle transverse momentum in GeV
177
178 ClassDef(TRootParticle, 1)
179};
180
181
182//---------------------------------------------------------------------------
183
184class TRootGenParticle: public TRootParticle {
185
186public:
187 TRootGenParticle() {_initialised=false;}
188 int PID; // particle HEP ID number [RawHepEventParticle::pid()]
189 int Status; // particle status [RawHepEventParticle::status()]
190 int M1; // particle 1st mother [RawHepEventParticle::mother1() - 1]
191 int M2; // particle 2nd mother [RawHepEventParticle::mother2() - 1]
192 int D1; // particle 1st daughter [RawHepEventParticle::daughter1() - 1]
193 int D2; // particle 2nd daughter [RawHepEventParticle::daughter2() - 1]
194
195 float T; // particle vertex position (t component) [RawHepEventParticle::t()]
196 float X; // particle vertex position (x component) [RawHepEventParticle::x()]
197 float Y; // particle vertex position (y component) [RawHepEventParticle::y()]
198 float Z; // particle vertex position (z component) [RawHepEventParticle::z()]
199 float M;
200 void setFractions();
201 const float getCharge() const {return Charge;};
202 const float getFem() {if(!_initialised) setFractions(); return _Fem;}
203 const float getFhad() {if(!_initialised) setFractions(); return _Fhad;}
204
205 float Charge; // electrical charge
206
207 static TCompare *fgCompare; //!
208 protected:
209 float _Fem, _Fhad; // fractions of energy deposit
210 bool _initialised;
211 ClassDef(TRootGenParticle, 1)
212};
213
214
215//------------------------------------------------------------------------------
216
217class TRootElectron: public TRootParticle
218{
219public:
220 TRootElectron() {};
221 int Charge; // particle Charge [RawHepEventParticle::pid()]
222 static TCompare *fgCompare; //!
223
224 bool IsolFlag;
225
226 ClassDef(TRootElectron, 1)
227};
228
229//------------------------------------------------------------------------------
230
231class TRootPhoton: public TRootParticle
232{
233public:
234 TRootPhoton() {};
235 static TCompare *fgCompare; //!
236
237 ClassDef(TRootPhoton, 1)
238};
239
240
241//------------------------------------------------------------------------------
242
243class TRootMuon: public TRootParticle
244{
245public:
246 TRootMuon() {};
247 int Charge; // particle Charge [RawHepEventParticle::pid()]
248 bool IsolFlag;
249 static TCompare *fgCompare; //!
250
251 ClassDef(TRootMuon, 1)
252};
253
254//---------------------------------------------------------------------------
255
256class TRootTracks : public TSortableObject {
257 public:
258 TRootTracks(); // needed for storage in ExRootAnalysis
259 TRootTracks(const TRootTracks& track);
260 TRootTracks(const float inEta, const float inPhi, const float outEta, const float outPhi, const float pt);
261 TRootTracks& operator=(const TRootTracks& track);
262 void Set(const float inEta, const float inPhi, const float outEta, const float outPhi, const float pt);
263 const TLorentzVector GetFourVector() const;
264 const float getEta() const {return Eta;}
265 const float getPhi() const {return Phi;}
266 const float getEtaOuter() const {return EtaOuter;}
267 const float getPhiOuter() const {return PhiOuter;}
268
269 static TCompare *fgCompare; //!
270 const TCompare *GetCompare() const { return fgCompare; }
271
272 float Eta, Phi; // (eta,phi) at the beginning of the track
273 float EtaOuter, PhiOuter; // (eta,phi) at the end of the track
274 float PT, E, Px, Py, Pz; // transverse momentum
275 ClassDef(TRootTracks, 1)
276};
277
278//---------------------------------------------------------------------------
279
280class TRootCalo: public TSortableObject
281{
282//class TRootCalo: public TRootParticle {
283 public:
284 float Eta;
285 float Phi;
286 float E;
287 TRootCalo() ;
288 TRootCalo(const TRootCalo& cal);
289 TRootCalo& operator=(const TRootCalo& cal);
290 void set(const D_CaloTower& cal);
291 static TCompare *fgCompare; //!
292 const TCompare *GetCompare() const { return fgCompare; }
293 const float getET() const {return ET;}
294
295 protected:
296 float E_em, E_had; // electromagnetic and hadronic components of the tower energy
297 float ET; // total energy and transverse energy
298 ClassDef(TRootCalo, 1)
299};
300
301//---------------------------------------------------------------------------
302class TRootZdcHits: public TRootParticle
303{
304public:
305 TRootZdcHits() {};
306 float T; // time of flight [s]
307 int side; // -1 or +1
308 static TCompare *fgCompare; //!
309
310 ClassDef(TRootZdcHits, 1)
311};
312
313//---------------------------------------------------------------------------
314
315class TRootTauJet: public TSortableObject
316{
317public:
318 TRootTauJet() {};
319
320 float E; // particle energy in GeV
321 float Px; // particle momentum vector (x component) in GeV
322 float Py; // particle momentum vector (y component) in GeV
323 float Pz; // particle momentum vector (z component) in GeV
324
325 float Eta; // particle pseudorapidity
326 float Phi; // particle azimuthal angle in rad
327
328 void Set(const TLorentzVector& momentum);// { return TRootParticle::Set(momentum); }
329
330 static TCompare *fgCompare; //!
331 const TCompare *GetCompare() const { return fgCompare; }
332
333 float PT; // particle transverse momentum in GeV
334
335 ClassDef(TRootTauJet, 1)
336};
337
338//---------------------------------------------------------------------------
339
340class TRootJet: public TRootTauJet
341{
342public:
343 TRootJet() {};
344
345 static TCompare *fgCompare; //!
346
347 bool Btag;
348
349 ClassDef(TRootJet, 1)
350};
351
352//------------------------------------------------------------------------------
353
354class TRootTrigger: public TSortableObject
355{
356public:
357 TRootTrigger() {};
358
359 int Accepted;
360
361 static TCompare *fgCompare; //!
362 const TCompare *GetCompare() const { return fgCompare; }
363
364 ClassDef(TRootTrigger, 1)
365};
366//---------------------------------------------------------------------------
367
368class TRootETmis: public TSortableObject
369{
370public:
371 TRootETmis() {};
372 float ET; // jet energy [RecJet::getEnergy()]
373 float Phi; // jet azimuthal angle [RecJet::getPhi()]
374 float Px;
375 float Py;
376
377 static TCompare *fgCompare; //!
378 const TCompare *GetCompare() const { return fgCompare; }
379
380 ClassDef(TRootETmis, 1)
381};
382
383//---------------------------------------------------------------------------
384
385class TRootRomanPotHits: public TSortableObject
386{
387public:
388 TRootRomanPotHits() {};
389 float T; // time of flight to the detector [s]
390 float S; // distance to the IP [m]
391 float E; // reconstructed energy [GeV]
392 float q2; // reconstructed squared momentum transfer [GeV^2]
393
394 float X; // horizontal distance to the beam [um]
395 float Y; // vertical distance to the beam [um]
396
397 float Tx; // angle of the momentum in the horizontal (x,z) plane [urad]
398 float Ty; // angle of the momentum in the verical (y,z) plane [urad]
399
400 int side; // -1 or 1
401
402 static TCompare *fgCompare; //!
403 const TCompare *GetCompare() const { return fgCompare; }
404
405 ClassDef(TRootRomanPotHits, 1)
406};
407
408#endif // BLOCKCLASSES_H
409
Note: See TracBrowser for help on using the repository browser.