Fork me on GitHub

source: git/classes/DelphesClasses.h @ d07e957

ImprovedOutputFileTimingdual_readoutllp
Last change on this file since d07e957 was d07e957, checked in by pavel <pavel@…>, 7 years ago

add vertex block

  • Property mode set to 100644
File size: 10.5 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 Vertex: public TObject
147{
148public:
149  Float_t T; // vertex position (t component)
150  Float_t X; // vertex position (x component)
151  Float_t Y; // vertex position (y component)
152  Float_t Z; // vertex position (z component)
153
154  ClassDef(Vertex, 1)
155};
156
157//---------------------------------------------------------------------------
158
159class MissingET: public TObject
160{
161public:
162  Float_t MET; // mising transverse energy
163  Float_t Eta; // mising energy pseudorapidity
164  Float_t Phi; // mising energy azimuthal angle
165
166  TLorentzVector P4();
167
168  ClassDef(MissingET, 1)
169};
170
171//---------------------------------------------------------------------------
172
173class ScalarHT: public TObject
174{
175public:
176  Float_t HT; // scalar sum of transverse momenta
177
178  ClassDef(ScalarHT, 1)
179};
180
181//---------------------------------------------------------------------------
182
183class Rho: public TObject
184{
185public:
186  Float_t Rho; // rho energy density
187  Float_t Edges[2]; // pseudorapidity range edges
188
189  ClassDef(Rho, 1)
190};
191
192//---------------------------------------------------------------------------
193
194class Weight: public TObject
195{
196public:
197  Float_t Weight; // weight for the event
198
199  ClassDef(Weight, 1)
200};
201
202//---------------------------------------------------------------------------
203
204class Photon: public SortableObject
205{
206public:
207
208  Float_t PT; // photon transverse momentum
209  Float_t Eta; // photon pseudorapidity
210  Float_t Phi; // photon azimuthal angle
211
212  Float_t E; // photon energy
213
214  Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
215
216  TRefArray Particles; // references to generated particles
217
218  static CompBase *fgCompare; //!
219  const CompBase *GetCompare() const { return fgCompare; }
220
221  TLorentzVector P4();
222
223  ClassDef(Photon, 2)
224};
225
226//---------------------------------------------------------------------------
227
228class Electron: public SortableObject
229{
230public:
231
232  Float_t PT; // electron transverse momentum
233  Float_t Eta; // electron pseudorapidity
234  Float_t Phi; // electron azimuthal angle
235
236  Int_t Charge; // electron charge
237
238  Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
239
240  TRef Particle; // reference to generated particle
241
242  static CompBase *fgCompare; //!
243  const CompBase *GetCompare() const { return fgCompare; }
244
245  TLorentzVector P4();
246
247  ClassDef(Electron, 2)
248};
249
250//---------------------------------------------------------------------------
251
252class Muon: public SortableObject
253{
254public:
255
256  Float_t PT; // muon transverse momentum
257  Float_t Eta; // muon pseudorapidity
258  Float_t Phi; // muon azimuthal angle
259
260  Int_t Charge; // muon charge
261
262  TRef Particle; // reference to generated particle
263
264  static CompBase *fgCompare; //!
265  const CompBase *GetCompare() const { return fgCompare; }
266
267  TLorentzVector P4();
268
269  ClassDef(Muon, 2)
270};
271
272//---------------------------------------------------------------------------
273
274class Jet: public SortableObject
275{
276public:
277
278  Float_t PT; // jet transverse momentum
279  Float_t Eta; // jet pseudorapidity
280  Float_t Phi; // jet azimuthal angle
281
282  Float_t Mass; // jet invariant mass
283
284  Float_t DeltaEta;  // jet radius in pseudorapidity
285  Float_t DeltaPhi;  // jet radius in azimuthal angle
286
287  UInt_t BTag; // 0 or 1 for a jet that has been tagged as containing a heavy quark
288  UInt_t TauTag; // 0 or 1 for a jet that has been tagged as a tau
289
290  Int_t Charge; // tau charge
291
292  Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
293
294  TRefArray Constituents; // references to constituents
295  TRefArray Particles; // references to generated particles
296
297  static CompBase *fgCompare; //!
298  const CompBase *GetCompare() const { return fgCompare; }
299
300  TLorentzVector P4();
301
302  ClassDef(Jet, 2)
303};
304
305//---------------------------------------------------------------------------
306
307class Track: public SortableObject
308{
309public:
310  Int_t PID; // HEP ID number
311
312  Int_t Charge; // track charge
313
314  Float_t PT; // track transverse momentum
315
316  Float_t Eta; // track pseudorapidity
317  Float_t Phi; // track azimuthal angle
318
319  Float_t EtaOuter; // track pseudorapidity at the tracker edge
320  Float_t PhiOuter; // track azimuthal angle at the tracker edge
321
322  Float_t X; // track vertex position (x component)
323  Float_t Y; // track vertex position (y component)
324  Float_t Z; // track vertex position (z component)
325
326  Float_t XOuter; // track position (x component) at the tracker edge
327  Float_t YOuter; // track position (y component) at the tracker edge
328  Float_t ZOuter; // track position (z component) at the tracker edge
329
330  TRef Particle; // reference to generated particle
331
332  static CompBase *fgCompare; //!
333  const CompBase *GetCompare() const { return fgCompare; }
334
335  TLorentzVector P4();
336
337  ClassDef(Track, 1)
338};
339
340//---------------------------------------------------------------------------
341
342class Tower: public SortableObject
343{
344public:
345  Float_t ET; // calorimeter tower transverse energy
346  Float_t Eta; // calorimeter tower pseudorapidity
347  Float_t Phi; // calorimeter tower azimuthal angle
348
349  Float_t E; // calorimeter tower energy
350
351  Float_t Eem; // calorimeter tower electromagnetic energy
352  Float_t Ehad; // calorimeter tower hadronic energy
353
354  Float_t Edges[4]; // calorimeter tower edges
355
356  TRefArray Particles; // references to generated particles
357
358  static CompBase *fgCompare; //!
359  const CompBase *GetCompare() const { return fgCompare; }
360
361  TLorentzVector P4();
362
363  ClassDef(Tower, 1)
364};
365
366//---------------------------------------------------------------------------
367
368class Candidate: public SortableObject
369{
370  friend class DelphesFactory;
371
372public:
373  Candidate();
374
375  Int_t PID;
376
377  Int_t Status;
378  Int_t M1, M2, D1, D2;
379
380  Int_t Charge;
381
382  Float_t Mass;
383
384  Int_t IsPU;
385  Int_t IsConstituent;
386
387  UInt_t BTag;
388  UInt_t TauTag;
389
390  Float_t Eem;
391  Float_t Ehad;
392
393  Float_t Edges[4];
394  Float_t DeltaEta;
395  Float_t DeltaPhi;
396
397  TLorentzVector Momentum, Position, Area;
398
399  static CompBase *fgCompare; //!
400  const CompBase *GetCompare() const { return fgCompare; }
401
402  void AddCandidate(Candidate *object);
403  TObjArray *GetCandidates();
404
405  Bool_t Overlaps(const Candidate *object) const;
406
407  virtual void Copy(TObject &object) const;
408  virtual TObject *Clone(const char *newname = "") const;
409  virtual void Clear(Option_t* option = "");
410
411private:
412  DelphesFactory *fFactory; //!
413  TObjArray *fArray; //!
414
415  void SetFactory(DelphesFactory *factory) { fFactory = factory; }
416
417  ClassDef(Candidate, 1)
418};
419
420#endif // DelphesClasses_h
421
422
Note: See TracBrowser for help on using the repository browser.