Fork me on GitHub

source: git/classes/DelphesClasses.h @ 3b465ca

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

update ROOT tree classes

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