Fork me on GitHub

source: git/classes/DelphesClasses.h@ 282f591

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 282f591 was e4c3fef, checked in by pavel <pavel@…>, 11 years ago

several minor corrections

  • Property mode set to 100644
File size: 12.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 Int_t M1; // particle 1st mother | hepevt.jmohep[number][0] - 1
110 Int_t M2; // particle 2nd mother | hepevt.jmohep[number][1] - 1
111
112 Int_t D1; // particle 1st daughter | hepevt.jdahep[number][0] - 1
113 Int_t D2; // particle last daughter | hepevt.jdahep[number][1] - 1
114
115 Int_t Charge; // particle charge
116
117 Float_t Mass; // particle mass
118
119 Float_t E; // particle energy | hepevt.phep[number][3]
120 Float_t Px; // particle momentum vector (x component) | hepevt.phep[number][0]
121 Float_t Py; // particle momentum vector (y component) | hepevt.phep[number][1]
122 Float_t Pz; // particle momentum vector (z component) | hepevt.phep[number][2]
123
124 Float_t PT; // particle transverse momentum
125 Float_t Eta; // particle pseudorapidity
126 Float_t Phi; // particle azimuthal angle
127
128 Float_t Rapidity; // particle rapidity
129
130 Float_t T; // particle vertex position (t component) | hepevt.vhep[number][3]
131 Float_t X; // particle vertex position (x component) | hepevt.vhep[number][0]
132 Float_t Y; // particle vertex position (y component) | hepevt.vhep[number][1]
133 Float_t Z; // particle vertex position (z component) | hepevt.vhep[number][2]
134
135 static CompBase *fgCompare; //!
136 const CompBase *GetCompare() const { return fgCompare; }
137
138 TLorentzVector P4();
139
140 ClassDef(GenParticle, 1)
141};
142
143//---------------------------------------------------------------------------
144
145class Vertex: public TObject
146{
147public:
148 Float_t T; // vertex position (t component)
149 Float_t X; // vertex position (x component)
150 Float_t Y; // vertex position (y component)
151 Float_t Z; // vertex position (z component)
152
153 ClassDef(Vertex, 1)
154};
155
156//---------------------------------------------------------------------------
157
158class MissingET: public TObject
159{
160public:
161 Float_t MET; // mising transverse energy
162 Float_t Eta; // mising energy pseudorapidity
163 Float_t Phi; // mising energy azimuthal angle
164
165 TLorentzVector P4();
166
167 ClassDef(MissingET, 1)
168};
169
170//---------------------------------------------------------------------------
171
172class ScalarHT: public TObject
173{
174public:
175 Float_t HT; // scalar sum of transverse momenta
176
177 ClassDef(ScalarHT, 1)
178};
179
180//---------------------------------------------------------------------------
181
182class Rho: public TObject
183{
184public:
185 Float_t Rho; // rho energy density
186 Float_t Edges[2]; // pseudorapidity range edges
187
188 ClassDef(Rho, 1)
189};
190
191//---------------------------------------------------------------------------
192
193class Weight: public TObject
194{
195public:
196 Float_t Weight; // weight for the event
197
198 ClassDef(Weight, 1)
199};
200
201//---------------------------------------------------------------------------
202
203class Photon: public SortableObject
204{
205public:
206
207 Float_t PT; // photon transverse momentum
208 Float_t Eta; // photon pseudorapidity
209 Float_t Phi; // photon azimuthal angle
210
211 Float_t E; // photon energy
212
213 Float_t T; //particle arrival time of flight
214
215 Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
216
217 TRefArray Particles; // references to generated particles
218
219 static CompBase *fgCompare; //!
220 const CompBase *GetCompare() const { return fgCompare; }
221
222 TLorentzVector P4();
223
224 ClassDef(Photon, 2)
225};
226
227//---------------------------------------------------------------------------
228
229class Electron: public SortableObject
230{
231public:
232
233 Float_t PT; // electron transverse momentum
234 Float_t Eta; // electron pseudorapidity
235 Float_t Phi; // electron azimuthal angle
236
237 Float_t T; //particle arrival time of flight
238
239 Int_t Charge; // electron charge
240
241 Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
242
243 TRef Particle; // reference to generated particle
244
245 static CompBase *fgCompare; //!
246 const CompBase *GetCompare() const { return fgCompare; }
247
248 TLorentzVector P4();
249
250 ClassDef(Electron, 2)
251};
252
253//---------------------------------------------------------------------------
254
255class Muon: public SortableObject
256{
257public:
258
259 Float_t PT; // muon transverse momentum
260 Float_t Eta; // muon pseudorapidity
261 Float_t Phi; // muon azimuthal angle
262
263 Float_t T; //particle arrival time of flight
264
265 Int_t Charge; // muon charge
266
267 TRef Particle; // reference to generated particle
268
269 static CompBase *fgCompare; //!
270 const CompBase *GetCompare() const { return fgCompare; }
271
272 TLorentzVector P4();
273
274 ClassDef(Muon, 2)
275};
276
277//---------------------------------------------------------------------------
278
279class Jet: public SortableObject
280{
281public:
282
283 Float_t PT; // jet transverse momentum
284 Float_t Eta; // jet pseudorapidity
285 Float_t Phi; // jet azimuthal angle
286
287 Float_t T; //particle arrival time of flight
288
289 Float_t Mass; // jet invariant mass
290
291 Float_t DeltaEta; // jet radius in pseudorapidity
292 Float_t DeltaPhi; // jet radius in azimuthal angle
293
294 UInt_t BTag; // 0 or 1 for a jet that has been tagged as containing a heavy quark
295 UInt_t TauTag; // 0 or 1 for a jet that has been tagged as a tau
296
297 Int_t Charge; // tau charge
298
299 Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
300
301 // -- PileUpJetID variables ---
302
303 Int_t NCharged;
304 Int_t NNeutrals;
305 Float_t Beta;
306 Float_t BetaStar;
307 Float_t MeanSqDeltaR;
308 Float_t PTD;
309 Float_t FracPt[5];
310
311 // -- Nsubjettiness variables ---
312
313 Float_t Tau1;
314 Float_t Tau2;
315 Float_t Tau3;
316 Float_t Tau4;
317 Float_t Tau5;
318
319 TRefArray Constituents; // references to constituents
320 TRefArray Particles; // references to generated particles
321
322 static CompBase *fgCompare; //!
323 const CompBase *GetCompare() const { return fgCompare; }
324
325 TLorentzVector P4();
326
327 ClassDef(Jet, 2)
328};
329
330//---------------------------------------------------------------------------
331
332class Track: public SortableObject
333{
334public:
335 Int_t PID; // HEP ID number
336
337 Int_t Charge; // track charge
338
339 Float_t PT; // track transverse momentum
340
341 Float_t Eta; // track pseudorapidity
342 Float_t Phi; // track azimuthal angle
343
344 Float_t EtaOuter; // track pseudorapidity at the tracker edge
345 Float_t PhiOuter; // track azimuthal angle at the tracker edge
346
347 Float_t X; // track vertex position (x component)
348 Float_t Y; // track vertex position (y component)
349 Float_t Z; // track vertex position (z component)
350 Float_t T; // track vertex position (z component)
351
352 Float_t XOuter; // track position (x component) at the tracker edge
353 Float_t YOuter; // track position (y component) at the tracker edge
354 Float_t ZOuter; // track position (z component) at the tracker edge
355 Float_t TOuter; // track position (z component) at the tracker edge
356
357 Float_t Dxy; // track signed transverse impact parameter
358 Float_t SDxy; // signed error on the track signed transverse impact parameter
359 Float_t Xd; // X coordinate of point of closest approach to vertex
360 Float_t Yd; // Y coordinate of point of closest approach to vertex
361 Float_t Zd; // Z coordinate of point of closest approach to vertex
362
363 TRef Particle; // reference to generated particle
364
365 static CompBase *fgCompare; //!
366 const CompBase *GetCompare() const { return fgCompare; }
367
368 TLorentzVector P4();
369
370 ClassDef(Track, 2)
371};
372
373//---------------------------------------------------------------------------
374
375class Tower: public SortableObject
376{
377public:
378 Float_t ET; // calorimeter tower transverse energy
379 Float_t Eta; // calorimeter tower pseudorapidity
380 Float_t Phi; // calorimeter tower azimuthal angle
381
382 Float_t E; // calorimeter tower energy
383
384 Float_t T; //particle arrival time of flight
385
386 Float_t Eem; // calorimeter tower electromagnetic energy
387 Float_t Ehad; // calorimeter tower hadronic energy
388
389 Float_t Edges[4]; // calorimeter tower edges
390
391 TRefArray Particles; // references to generated particles
392
393 static CompBase *fgCompare; //!
394 const CompBase *GetCompare() const { return fgCompare; }
395
396 TLorentzVector P4();
397
398 ClassDef(Tower, 1)
399};
400
401//---------------------------------------------------------------------------
402
403class HectorHit: public SortableObject
404{
405public:
406 Float_t E; // reconstructed energy [GeV]
407
408 Float_t Tx; // angle of the momentum in the horizontal (x,z) plane [urad]
409 Float_t Ty; // angle of the momentum in the verical (y,z) plane [urad]
410
411 Float_t T; // time of flight to the detector [s]
412
413 Float_t X; // horizontal distance to the beam [um]
414 Float_t Y; // vertical distance to the beam [um]
415 Float_t S; // distance to the interaction point [m]
416
417 TRef Particle; // reference to generated particle
418
419 static CompBase *fgCompare; //!
420 const CompBase *GetCompare() const { return fgCompare; }
421
422 ClassDef(HectorHit, 1)
423};
424
425//---------------------------------------------------------------------------
426
427class Candidate: public SortableObject
428{
429 friend class DelphesFactory;
430
431public:
432 Candidate();
433
434 Int_t PID;
435
436 Int_t Status;
437 Int_t M1, M2, D1, D2;
438
439 Int_t Charge;
440
441 Float_t Mass;
442
443 Int_t IsPU;
444 Int_t IsConstituent;
445
446 UInt_t BTag;
447 UInt_t TauTag;
448
449 Float_t Eem;
450 Float_t Ehad;
451
452 Float_t Edges[4];
453 Float_t DeltaEta;
454 Float_t DeltaPhi;
455
456 TLorentzVector Momentum, Position, Area;
457
458 Float_t Dxy;
459 Float_t SDxy;
460 Float_t Xd;
461 Float_t Yd;
462 Float_t Zd;
463
464 // PileUpJetID variables
465
466 Int_t NCharged;
467 Int_t NNeutrals;
468 Float_t Beta;
469 Float_t BetaStar;
470 Float_t MeanSqDeltaR;
471 Float_t PTD;
472 Float_t FracPt[5];
473
474 // -- Nsubjettiness variables ---
475
476 Float_t Tau[5];
477
478 static CompBase *fgCompare; //!
479 const CompBase *GetCompare() const { return fgCompare; }
480
481 void AddCandidate(Candidate *object);
482 TObjArray *GetCandidates();
483
484 Bool_t Overlaps(const Candidate *object) const;
485
486 virtual void Copy(TObject &object) const;
487 virtual TObject *Clone(const char *newname = "") const;
488 virtual void Clear(Option_t* option = "");
489
490private:
491 DelphesFactory *fFactory; //!
492 TObjArray *fArray; //!
493
494 void SetFactory(DelphesFactory *factory) { fFactory = factory; }
495
496 ClassDef(Candidate, 2)
497};
498
499#endif // DelphesClasses_h
500
501
Note: See TracBrowser for help on using the repository browser.