Fork me on GitHub

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

ImprovedOutputFile Timing dual_readout llp 3.0.9
Last change on this file since 3ccc2b6e was 71648c2, checked in by pavel <pavel@…>, 11 years ago

add Rho branch

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