Fork me on GitHub

source: svn/trunk/classes/DelphesClasses.h@ 1318

Last change on this file since 1318 was 1318, checked in by Pavel Demin, 11 years ago

add Eta to MissingET

File size: 10.2 KB
Line 
1#ifndef DelphesClasses_h
2#define DelphesClasses_h
3
4/**
5 *
6 * Definition of classes to be stored in the root tree.
7 * Function CompareXYZ sorts objects by the variable XYZ that MUST be
8 * present in the data members of the root tree class of the branch.
9 *
10 * $Date: 2008-06-04 13:57:24 $
11 * $Revision: 1.1 $
12 *
13 *
14 * \author P. Demin - UCL, Louvain-la-Neuve
15 *
16 */
17
18// Dependencies (#includes)
19
20#include "TRef.h"
21#include "TObject.h"
22#include "TRefArray.h"
23#include "TLorentzVector.h"
24
25#include "classes/SortableObject.h"
26
27class DelphesFactory;
28
29//---------------------------------------------------------------------------
30
31class Event: public TObject
32{
33public:
34
35 Long64_t Number; // event number
36
37 Float_t ReadTime;
38 Float_t ProcTime;
39
40 ClassDef(Event, 1)
41};
42
43//---------------------------------------------------------------------------
44
45class LHCOEvent: public Event
46{
47public:
48
49 Int_t Trigger; // trigger word
50
51 ClassDef(LHCOEvent, 1)
52};
53
54//---------------------------------------------------------------------------
55
56class LHEFEvent: public Event
57{
58public:
59
60 Int_t ProcessID; // subprocess code for the event | hepup.IDPRUP
61
62 Float_t Weight; // weight for the event | hepup.XWGTUP
63 Float_t ScalePDF; // scale in GeV used in the calculation of the PDFs in the event | hepup.SCALUP
64 Float_t AlphaQED; // value of the QED coupling used in the event | hepup.AQEDUP
65 Float_t AlphaQCD; // value of the QCD coupling used in the event | hepup.AQCDUP
66
67 ClassDef(LHEFEvent, 2)
68};
69
70//---------------------------------------------------------------------------
71
72class HepMCEvent: public Event
73{
74public:
75
76 Int_t ProcessID; // unique signal process id | signal_process_id()
77 Int_t MPI; // number of multi parton interactions | mpi ()
78
79 Float_t Weight; // weight for the event
80
81 Float_t Scale; // energy scale, see hep-ph/0109068 | event_scale()
82 Float_t AlphaQED; // QED coupling, see hep-ph/0109068 | alphaQED()
83 Float_t AlphaQCD; // QCD coupling, see hep-ph/0109068 | alphaQCD()
84
85 Int_t ID1; // flavour code of first parton | pdf_info()->id1()
86 Int_t ID2; // flavour code of second parton | pdf_info()->id2()
87
88 Float_t X1; // fraction of beam momentum carried by first parton ("beam side") | pdf_info()->x1()
89 Float_t X2; // fraction of beam momentum carried by second parton ("target side") | pdf_info()->x2()
90
91 Float_t ScalePDF; // Q-scale used in evaluation of PDF's (in GeV) | pdf_info()->scalePDF()
92
93 Float_t PDF1; // PDF (id1, x1, Q) | pdf_info()->pdf1()
94 Float_t PDF2; // PDF (id2, x2, Q) | pdf_info()->pdf2()
95
96 ClassDef(HepMCEvent, 2)
97};
98
99//---------------------------------------------------------------------------
100
101class GenParticle: public SortableObject
102{
103public:
104 Int_t PID; // particle HEP ID number | hepevt.idhep[number]
105
106 Int_t Status; // particle status | hepevt.isthep[number]
107 Int_t IsPU; // 0 or 1 for particles from pile-up interactions
108
109
110 Int_t M1; // particle 1st mother | hepevt.jmohep[number][0] - 1
111 Int_t M2; // particle 2nd mother | hepevt.jmohep[number][1] - 1
112
113 Int_t D1; // particle 1st daughter | hepevt.jdahep[number][0] - 1
114 Int_t D2; // particle last daughter | hepevt.jdahep[number][1] - 1
115
116 Int_t Charge; // particle charge
117
118 Float_t Mass; // particle mass
119
120 Float_t E; // particle energy | hepevt.phep[number][3]
121 Float_t Px; // particle momentum vector (x component) | hepevt.phep[number][0]
122 Float_t Py; // particle momentum vector (y component) | hepevt.phep[number][1]
123 Float_t Pz; // particle momentum vector (z component) | hepevt.phep[number][2]
124
125 Float_t PT; // particle transverse momentum
126 Float_t Eta; // particle pseudorapidity
127 Float_t Phi; // particle azimuthal angle
128
129 Float_t Rapidity; // particle rapidity
130
131 Float_t T; // particle vertex position (t component) | hepevt.vhep[number][3]
132 Float_t X; // particle vertex position (x component) | hepevt.vhep[number][0]
133 Float_t Y; // particle vertex position (y component) | hepevt.vhep[number][1]
134 Float_t Z; // particle vertex position (z component) | hepevt.vhep[number][2]
135
136 static CompBase *fgCompare; //!
137 const CompBase *GetCompare() const { return fgCompare; }
138
139 TLorentzVector P4();
140
141 ClassDef(GenParticle, 1)
142};
143
144//---------------------------------------------------------------------------
145
146class MissingET: public TObject
147{
148public:
149 Float_t MET; // mising transverse energy
150 Float_t Eta; // mising energy pseudorapidity
151 Float_t Phi; // mising energy azimuthal angle
152
153 TLorentzVector P4();
154
155 ClassDef(MissingET, 1)
156};
157
158//---------------------------------------------------------------------------
159
160class ScalarHT: public TObject
161{
162public:
163 Float_t HT; // scalar sum of transverse momenta
164
165 ClassDef(ScalarHT, 1)
166};
167
168//---------------------------------------------------------------------------
169
170class Rho: public TObject
171{
172public:
173 Float_t Rho; // rho energy density
174
175 ClassDef(Rho, 1)
176};
177
178//---------------------------------------------------------------------------
179
180class Weight: public TObject
181{
182public:
183 Float_t Weight; // weight for the event
184
185 ClassDef(Weight, 1)
186};
187
188//---------------------------------------------------------------------------
189
190class Photon: public SortableObject
191{
192public:
193
194 Float_t PT; // photon transverse momentum
195 Float_t Eta; // photon pseudorapidity
196 Float_t Phi; // photon azimuthal angle
197
198 Float_t E; // photon energy
199
200 Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
201
202 TRefArray Particles; // references to generated particles
203
204 static CompBase *fgCompare; //!
205 const CompBase *GetCompare() const { return fgCompare; }
206
207 TLorentzVector P4();
208
209 ClassDef(Photon, 2)
210};
211
212//---------------------------------------------------------------------------
213
214class Electron: public SortableObject
215{
216public:
217
218 Float_t PT; // electron transverse momentum
219 Float_t Eta; // electron pseudorapidity
220 Float_t Phi; // electron azimuthal angle
221
222 Int_t Charge; // electron charge
223
224 Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
225
226 TRef Particle; // reference to generated particle
227
228 static CompBase *fgCompare; //!
229 const CompBase *GetCompare() const { return fgCompare; }
230
231 TLorentzVector P4();
232
233 ClassDef(Electron, 2)
234};
235
236//---------------------------------------------------------------------------
237
238class Muon: public SortableObject
239{
240public:
241
242 Float_t PT; // muon transverse momentum
243 Float_t Eta; // muon pseudorapidity
244 Float_t Phi; // muon azimuthal angle
245
246 Int_t Charge; // muon charge
247
248 TRef Particle; // reference to generated particle
249
250 static CompBase *fgCompare; //!
251 const CompBase *GetCompare() const { return fgCompare; }
252
253 TLorentzVector P4();
254
255 ClassDef(Muon, 2)
256};
257
258//---------------------------------------------------------------------------
259
260class Jet: public SortableObject
261{
262public:
263
264 Float_t PT; // jet transverse momentum
265 Float_t Eta; // jet pseudorapidity
266 Float_t Phi; // jet azimuthal angle
267
268 Float_t Mass; // jet invariant mass
269
270 Float_t DeltaEta; // jet radius in pseudorapidity
271 Float_t DeltaPhi; // jet radius in azimuthal angle
272
273 UInt_t BTag; // 0 or 1 for a jet that has been tagged as containing a heavy quark
274 UInt_t TauTag; // 0 or 1 for a jet that has been tagged as a tau
275
276 Int_t Charge; // tau charge
277
278 Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
279
280 TRefArray Constituents; // references to constituents
281 TRefArray Particles; // references to generated particles
282
283 static CompBase *fgCompare; //!
284 const CompBase *GetCompare() const { return fgCompare; }
285
286 TLorentzVector P4();
287
288 ClassDef(Jet, 2)
289};
290
291//---------------------------------------------------------------------------
292
293class Track: public SortableObject
294{
295public:
296 Int_t PID; // HEP ID number
297
298 Int_t Charge; // track charge
299
300 Float_t PT; // track transverse momentum
301
302 Float_t Eta; // track pseudorapidity
303 Float_t Phi; // track azimuthal angle
304
305 Float_t EtaOuter; // track pseudorapidity at the tracker edge
306 Float_t PhiOuter; // track azimuthal angle at the tracker edge
307
308 Float_t X; // track vertex position (x component)
309 Float_t Y; // track vertex position (y component)
310 Float_t Z; // track vertex position (z component)
311
312 Float_t XOuter; // track position (x component) at the tracker edge
313 Float_t YOuter; // track position (y component) at the tracker edge
314 Float_t ZOuter; // track position (z component) at the tracker edge
315
316 TRef Particle; // reference to generated particle
317
318 static CompBase *fgCompare; //!
319 const CompBase *GetCompare() const { return fgCompare; }
320
321 TLorentzVector P4();
322
323 ClassDef(Track, 1)
324};
325
326//---------------------------------------------------------------------------
327
328class Tower: public SortableObject
329{
330public:
331 Float_t ET; // calorimeter tower transverse energy
332 Float_t Eta; // calorimeter tower pseudorapidity
333 Float_t Phi; // calorimeter tower azimuthal angle
334
335 Float_t E; // calorimeter tower energy
336
337 Float_t Eem; // calorimeter tower electromagnetic energy
338 Float_t Ehad; // calorimeter tower hadronic energy
339
340 Float_t Edges[4]; // calorimeter tower edges
341
342 TRefArray Particles; // references to generated particles
343
344 static CompBase *fgCompare; //!
345 const CompBase *GetCompare() const { return fgCompare; }
346
347 TLorentzVector P4();
348
349 ClassDef(Tower, 1)
350};
351
352//---------------------------------------------------------------------------
353
354class Candidate: public SortableObject
355{
356 friend class DelphesFactory;
357
358public:
359 Candidate();
360
361 Int_t PID;
362
363 Int_t Status;
364 Int_t M1, M2, D1, D2;
365
366 Int_t Charge;
367
368 Float_t Mass;
369
370 Int_t IsPU;
371 Int_t IsConstituent;
372
373 UInt_t BTag;
374 UInt_t TauTag;
375
376 Float_t Eem;
377 Float_t Ehad;
378
379 Float_t Edges[4];
380 Float_t DeltaEta;
381 Float_t DeltaPhi;
382
383 TLorentzVector Momentum, Position, Area;
384
385 static CompBase *fgCompare; //!
386 const CompBase *GetCompare() const { return fgCompare; }
387
388 void AddCandidate(Candidate *object);
389 TObjArray *GetCandidates();
390
391 Bool_t Overlaps(const Candidate *object) const;
392
393 virtual void Copy(TObject &object) const;
394 virtual TObject *Clone(const char *newname = "") const;
395 virtual void Clear(Option_t* option = "");
396
397private:
398 DelphesFactory *fFactory; //!
399 TObjArray *fArray; //!
400
401 void SetFactory(DelphesFactory *factory) { fFactory = factory; }
402
403 ClassDef(Candidate, 1)
404};
405
406#endif // DelphesClasses_h
407
408
Note: See TracBrowser for help on using the repository browser.