Fork me on GitHub

source: git/classes/DelphesClasses.h @ acd0621

ImprovedOutputFileTimingdual_readoutllp
Last change on this file since acd0621 was acd0621, checked in by Michele Selvaggi <michele.selvaggi@…>, 4 years ago

including tracking variables

  • Property mode set to 100644
File size: 16.9 KB
Line 
1/*
2 *  Delphes: a framework for fast simulation of a generic collider experiment
3 *  Copyright (C) 2012-2014  Universite catholique de Louvain (UCL), Belgium
4 *
5 *  This program is free software: you can redistribute it and/or modify
6 *  it under the terms of the GNU General Public License as published by
7 *  the Free Software Foundation, either version 3 of the License, or
8 *  (at your option) any later version.
9 *
10 *  This program is distributed in the hope that it will be useful,
11 *  but WITHOUT ANY WARRANTY; without even the implied warranty of
12 *  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
13 *  GNU General Public License for more details.
14 *
15 *  You should have received a copy of the GNU General Public License
16 *  along with this program.  If not, see <http://www.gnu.org/licenses/>.
17 */
18
19#ifndef DelphesClasses_h
20#define DelphesClasses_h
21
22/**
23 *
24 *  Definition of classes to be stored in the root tree.
25 *  Function CompareXYZ sorts objects by the variable XYZ that MUST be
26 *  present in the data members of the root tree class of the branch.
27 *
28 *  \author P. Demin - UCL, Louvain-la-Neuve
29 *
30 */
31
32// Dependencies (#includes)
33
34#include "TRef.h"
35#include "TObject.h"
36#include "TRefArray.h"
37#include "TLorentzVector.h"
38
39#include "classes/SortableObject.h"
40
41class DelphesFactory;
42
43//---------------------------------------------------------------------------
44
45class Event: public TObject
46{
47public:
48
49  Long64_t Number; // event number
50
51  Float_t ReadTime;
52  Float_t ProcTime;
53
54  ClassDef(Event, 1)
55};
56
57//---------------------------------------------------------------------------
58
59class LHCOEvent: public Event
60{
61public:
62
63  Int_t Trigger; // trigger word
64
65  ClassDef(LHCOEvent, 1)
66};
67
68//---------------------------------------------------------------------------
69
70class LHEFEvent: public Event
71{
72public:
73
74  Int_t ProcessID; // subprocess code for the event | hepup.IDPRUP
75
76  Float_t Weight; // weight for the event | hepup.XWGTUP
77  Float_t ScalePDF; // scale in GeV used in the calculation of the PDFs in the event | hepup.SCALUP
78  Float_t AlphaQED; // value of the QED coupling used in the event | hepup.AQEDUP
79  Float_t AlphaQCD; // value of the QCD coupling used in the event | hepup.AQCDUP
80
81  ClassDef(LHEFEvent, 2)
82};
83
84//---------------------------------------------------------------------------
85
86class LHEFWeight: public TObject
87{
88public:
89  Int_t ID; // weight ID
90  Float_t Weight; // weight value
91
92  ClassDef(LHEFWeight, 1)
93};
94
95//---------------------------------------------------------------------------
96
97class HepMCEvent: public Event
98{
99public:
100
101  Int_t ProcessID; // unique signal process id | signal_process_id()
102  Int_t MPI; // number of multi parton interactions | mpi ()
103
104  Float_t Weight; // weight for the event
105
106  Float_t Scale; // energy scale, see hep-ph/0109068 | event_scale()
107  Float_t AlphaQED; // QED coupling, see hep-ph/0109068 | alphaQED()
108  Float_t AlphaQCD; // QCD coupling, see hep-ph/0109068 | alphaQCD()
109
110  Int_t ID1; // flavour code of first parton | pdf_info()->id1()
111  Int_t ID2; // flavour code of second parton | pdf_info()->id2()
112
113  Float_t X1; // fraction of beam momentum carried by first parton ("beam side") | pdf_info()->x1()
114  Float_t X2; // fraction of beam momentum carried by second parton ("target side") | pdf_info()->x2()
115
116  Float_t ScalePDF; // Q-scale used in evaluation of PDF's (in GeV) | pdf_info()->scalePDF()
117
118  Float_t PDF1; // PDF (id1, x1, Q) | pdf_info()->pdf1()
119  Float_t PDF2; // PDF (id2, x2, Q) | pdf_info()->pdf2()
120
121  ClassDef(HepMCEvent, 2)
122};
123
124//---------------------------------------------------------------------------
125
126class GenParticle: public SortableObject
127{
128public:
129  Int_t PID; // particle HEP ID number | hepevt.idhep[number]
130
131  Int_t Status; // particle status | hepevt.isthep[number]
132  Int_t IsPU; // 0 or 1 for particles from pile-up interactions
133
134  Int_t M1; // particle 1st mother | hepevt.jmohep[number][0] - 1
135  Int_t M2; // particle 2nd mother | hepevt.jmohep[number][1] - 1
136
137  Int_t D1; // particle 1st daughter | hepevt.jdahep[number][0] - 1
138  Int_t D2; // particle last daughter | hepevt.jdahep[number][1] - 1
139
140  Int_t Charge; // particle charge
141
142  Float_t Mass; // particle mass
143
144  Float_t E; // particle energy | hepevt.phep[number][3]
145  Float_t Px; // particle momentum vector (x component) | hepevt.phep[number][0]
146  Float_t Py; // particle momentum vector (y component) | hepevt.phep[number][1]
147  Float_t Pz; // particle momentum vector (z component) | hepevt.phep[number][2]
148
149  Float_t D0;
150  Float_t DZ;
151  Float_t P;
152  Float_t PT;
153  Float_t CtgTheta;
154  Float_t Phi;
155  Float_t Eta; // particle pseudorapidity
156  Float_t Rapidity; // particle rapidity
157
158  Float_t T; // particle vertex position (t component) | hepevt.vhep[number][3]
159  Float_t X; // particle vertex position (x component) | hepevt.vhep[number][0]
160  Float_t Y; // particle vertex position (y component) | hepevt.vhep[number][1]
161  Float_t Z; // particle vertex position (z component) | hepevt.vhep[number][2]
162
163  static CompBase *fgCompare; //!
164  const CompBase *GetCompare() const { return fgCompare; }
165
166  TLorentzVector P4() const;
167
168  ClassDef(GenParticle, 1)
169};
170
171//---------------------------------------------------------------------------
172
173class Vertex: public TObject
174{
175public:
176  Float_t T; // vertex position (t component)
177  Float_t X; // vertex position (x component)
178  Float_t Y; // vertex position (y component)
179  Float_t Z; // vertex position (z component)
180
181  ClassDef(Vertex, 1)
182};
183
184//---------------------------------------------------------------------------
185
186class MissingET: public TObject
187{
188public:
189  Float_t MET; // mising transverse energy
190  Float_t Eta; // mising energy pseudorapidity
191  Float_t Phi; // mising energy azimuthal angle
192
193  TLorentzVector P4() const;
194
195  ClassDef(MissingET, 1)
196};
197
198//---------------------------------------------------------------------------
199
200class ScalarHT: public TObject
201{
202public:
203  Float_t HT; // scalar sum of transverse momenta
204
205  ClassDef(ScalarHT, 1)
206};
207
208//---------------------------------------------------------------------------
209
210class Rho: public TObject
211{
212public:
213  Float_t Rho; // rho energy density
214  Float_t Edges[2]; // pseudorapidity range edges
215
216  ClassDef(Rho, 1)
217};
218
219//---------------------------------------------------------------------------
220
221class Weight: public TObject
222{
223public:
224  Float_t Weight; // weight for the event
225
226  ClassDef(Weight, 1)
227};
228
229//---------------------------------------------------------------------------
230
231class Photon: public SortableObject
232{
233public:
234
235  Float_t PT; // photon transverse momentum
236  Float_t Eta; // photon pseudorapidity
237  Float_t Phi; // photon azimuthal angle
238
239  Float_t E; // photon energy
240
241  Float_t T; //particle arrival time of flight
242
243  Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
244
245  TRefArray Particles; // references to generated particles
246
247  // Isolation variables
248
249  Float_t IsolationVar;
250  Float_t IsolationVarRhoCorr;
251  Float_t SumPtCharged;
252  Float_t SumPtNeutral;
253  Float_t SumPtChargedPU;
254  Float_t SumPt;
255
256  static CompBase *fgCompare; //!
257  const CompBase *GetCompare() const { return fgCompare; }
258
259  TLorentzVector P4() const;
260
261  ClassDef(Photon, 3)
262};
263
264//---------------------------------------------------------------------------
265
266class Electron: public SortableObject
267{
268public:
269
270  Float_t PT; // electron transverse momentum
271  Float_t Eta; // electron pseudorapidity
272  Float_t Phi; // electron azimuthal angle
273
274  Float_t T; //particle arrival time of flight
275
276  Int_t Charge; // electron charge
277
278  Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
279
280  TRef Particle; // reference to generated particle
281
282  // Isolation variables
283
284  Float_t IsolationVar;
285  Float_t IsolationVarRhoCorr;
286  Float_t SumPtCharged;
287  Float_t SumPtNeutral;
288  Float_t SumPtChargedPU;
289  Float_t SumPt;
290
291  static CompBase *fgCompare; //!
292  const CompBase *GetCompare() const { return fgCompare; }
293
294  TLorentzVector P4() const;
295
296  ClassDef(Electron, 3)
297};
298
299//---------------------------------------------------------------------------
300
301class Muon: public SortableObject
302{
303public:
304
305  Float_t PT; // muon transverse momentum
306  Float_t Eta; // muon pseudorapidity
307  Float_t Phi; // muon azimuthal angle
308
309  Float_t T; //particle arrival time of flight
310
311  Int_t Charge; // muon charge
312
313  TRef Particle; // reference to generated particle
314
315   // Isolation variables
316
317  Float_t IsolationVar;
318  Float_t IsolationVarRhoCorr;
319  Float_t SumPtCharged;
320  Float_t SumPtNeutral;
321  Float_t SumPtChargedPU;
322  Float_t SumPt;
323
324  static CompBase *fgCompare; //!
325  const CompBase *GetCompare() const { return fgCompare; }
326
327  TLorentzVector P4() const;
328
329  ClassDef(Muon, 3)
330};
331
332//---------------------------------------------------------------------------
333
334class Jet: public SortableObject
335{
336public:
337
338  Float_t PT; // jet transverse momentum
339  Float_t Eta; // jet pseudorapidity
340  Float_t Phi; // jet azimuthal angle
341
342  Float_t T; //particle arrival time of flight
343
344  Float_t Mass; // jet invariant mass
345
346  Float_t DeltaEta;  // jet radius in pseudorapidity
347  Float_t DeltaPhi;  // jet radius in azimuthal angle
348
349  UInt_t Flavor;
350  UInt_t FlavorAlgo;
351  UInt_t FlavorPhys;
352
353  UInt_t BTag; // 0 or 1 for a jet that has been tagged as containing a heavy quark
354  UInt_t BTagAlgo;
355  UInt_t BTagPhys;
356
357  UInt_t TauTag; // 0 or 1 for a jet that has been tagged as a tau
358
359  Int_t Charge; // tau charge
360
361  Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
362
363  Int_t NCharged; // number of charged constituents
364  Int_t NNeutrals; // number of neutral constituents
365  Float_t Beta; // (sum pt of charged pile-up constituents)/(sum pt of charged constituents)
366  Float_t BetaStar; // (sum pt of charged constituents coming from hard interaction)/(sum pt of charged constituents)
367  Float_t MeanSqDeltaR; // average distance (squared) between constituent and jet weighted by pt (squared) of constituent
368  Float_t PTD; // average pt between constituent and jet weighted by pt of constituent
369  Float_t FracPt[5]; // (sum pt of constituents within a ring 0.1*i < DeltaR < 0.1*(i+1))/(sum pt of constituents)
370
371  Float_t Tau[5]; // N-subjettiness
372
373  TLorentzVector TrimmedP4[5]; // first entry (i = 0) is the total Trimmed Jet 4-momenta and from i = 1 to 4 are the trimmed subjets 4-momenta
374  TLorentzVector PrunedP4[5]; // first entry (i = 0) is the total Pruned Jet 4-momenta and from i = 1 to 4 are the pruned subjets 4-momenta
375  TLorentzVector SoftDroppedP4[5]; // first entry (i = 0) is the total SoftDropped Jet 4-momenta and from i = 1 to 4 are the pruned subjets 4-momenta
376
377  Int_t NSubJetsTrimmed; // number of subjets trimmed
378  Int_t NSubJetsPruned; // number of subjets pruned
379  Int_t NSubJetsSoftDropped; // number of subjets soft-dropped
380
381  TRefArray Constituents; // references to constituents
382  TRefArray Particles; // references to generated particles
383
384  static CompBase *fgCompare; //!
385  const CompBase *GetCompare() const { return fgCompare; }
386
387  TLorentzVector P4() const;
388  TLorentzVector Area;
389
390  ClassDef(Jet, 3)
391};
392
393//---------------------------------------------------------------------------
394
395class Track: public SortableObject
396{
397public:
398  Int_t PID; // HEP ID number
399
400  Int_t Charge; // track charge
401
402  Float_t Eta; // track pseudorapidity
403 
404  Float_t EtaOuter; // track pseudorapidity at the tracker edge
405  Float_t PhiOuter; // track azimuthal angle at the tracker edge
406
407  Float_t X; // track vertex position (x component)
408  Float_t Y; // track vertex position (y component)
409  Float_t Z; // track vertex position (z component)
410  Float_t T; // track vertex position (z component)
411
412  Float_t XOuter; // track position (x component) at the tracker edge
413  Float_t YOuter; // track position (y component) at the tracker edge
414  Float_t ZOuter; // track position (z component) at the tracker edge
415  Float_t TOuter; // track position (z component) at the tracker edge
416
417  Float_t L; // track path length
418   
419  Float_t D0;     // track signed transverse impact parameter
420  Float_t ErrorD0;    // signed error on the track signed transverse impact parameter
421 
422  Float_t DZ; // track transverse momentum
423  Float_t ErrorDZ; // track transverse momentum error
424 
425  Float_t P; // track transverse momentum
426  Float_t ErrorP; // track transverse momentum error
427 
428  Float_t PT; // track transverse momentum
429  Float_t ErrorPT; // track transverse momentum error
430 
431  Float_t CtgTheta; // track transverse momentum
432  Float_t ErrorCtgTheta; // track transverse momentum error
433 
434  Float_t Phi; // track azimuthal angle
435  Float_t ErrorPhi; // track azimuthal angle
436 
437  Float_t Xd;      // X coordinate of point of closest approach to vertex
438  Float_t Yd;      // Y coordinate of point of closest approach to vertex
439  Float_t Zd;      // Z coordinate of point of closest approach to vertex
440
441  TRef Particle; // reference to generated particle
442
443  static CompBase *fgCompare; //!
444  const CompBase *GetCompare() const { return fgCompare; }
445
446  TLorentzVector P4() const;
447
448  ClassDef(Track, 2)
449};
450
451//---------------------------------------------------------------------------
452
453class Tower: public SortableObject
454{
455public:
456  Float_t ET; // calorimeter tower transverse energy
457  Float_t Eta; // calorimeter tower pseudorapidity
458  Float_t Phi; // calorimeter tower azimuthal angle
459
460  Float_t E; // calorimeter tower energy
461
462  Float_t T; // ecal deposit time, averaged by sqrt(EM energy) over all particles, not smeared
463  Int_t NTimeHits; // number of hits contributing to time measurement
464
465  Float_t Eem; // calorimeter tower electromagnetic energy
466  Float_t Ehad; // calorimeter tower hadronic energy
467
468  Float_t Edges[4]; // calorimeter tower edges
469
470  TRefArray Particles; // references to generated particles
471
472  static CompBase *fgCompare; //!
473  const CompBase *GetCompare() const { return fgCompare; }
474
475  TLorentzVector P4() const;
476
477  ClassDef(Tower, 2)
478};
479
480//---------------------------------------------------------------------------
481
482class HectorHit: public SortableObject
483{
484public:
485  Float_t E; // reconstructed energy [GeV]
486
487  Float_t Tx; // angle of the momentum in the horizontal (x,z) plane [urad]
488  Float_t Ty; // angle of the momentum in the verical (y,z) plane [urad]
489
490  Float_t T; // time of flight to the detector [s]
491
492  Float_t X; // horizontal distance to the beam [um]
493  Float_t Y; // vertical distance to the beam [um]
494  Float_t S; // distance to the interaction point [m]
495
496  TRef Particle; // reference to generated particle
497
498  static CompBase *fgCompare; //!
499  const CompBase *GetCompare() const { return fgCompare; }
500
501  ClassDef(HectorHit, 1)
502};
503
504//---------------------------------------------------------------------------
505
506class Candidate: public SortableObject
507{
508  friend class DelphesFactory;
509
510public:
511  Candidate();
512
513  Int_t PID;
514
515  Int_t Status;
516  Int_t M1, M2, D1, D2;
517
518  Int_t Charge;
519
520  Float_t Mass;
521
522  Int_t IsPU;
523  Int_t IsRecoPU;
524
525  Int_t IsConstituent;
526
527  Int_t IsFromConversion;
528
529  UInt_t Flavor;
530  UInt_t FlavorAlgo;
531  UInt_t FlavorPhys;
532
533  UInt_t BTag;
534  UInt_t BTagAlgo;
535  UInt_t BTagPhys;
536
537  UInt_t TauTag;
538
539  Float_t Eem;
540  Float_t Ehad;
541
542  Float_t Edges[4];
543  Float_t DeltaEta;
544  Float_t DeltaPhi;
545
546  TLorentzVector Momentum, Position, InitialPosition, Area;
547
548  Float_t L; // path length
549  Float_t D0;
550  Float_t ErrorD0;
551  Float_t DZ;
552  Float_t ErrorDZ;
553  Float_t P;
554  Float_t ErrorP;
555  Float_t PT;
556  Float_t ErrorPT;
557  Float_t CtgTheta;
558  Float_t ErrorCtgTheta;
559  Float_t Phi;
560  Float_t ErrorPhi;
561
562  Float_t Xd;
563  Float_t Yd;
564  Float_t Zd;
565
566  // tracking resolution
567 
568  Float_t TrackResolution;
569
570  // PileUpJetID variables
571
572  Int_t NCharged;
573  Int_t NNeutrals;
574  Float_t Beta;
575  Float_t BetaStar;
576  Float_t MeanSqDeltaR;
577  Float_t PTD;
578  Float_t FracPt[5];
579
580  // Timing information
581
582  Int_t NTimeHits;
583  std::vector< std::pair< Float_t, Float_t > > ECalEnergyTimePairs;
584
585  // Isolation variables
586
587  Float_t IsolationVar;
588  Float_t IsolationVarRhoCorr;
589  Float_t SumPtCharged;
590  Float_t SumPtNeutral;
591  Float_t SumPtChargedPU;
592  Float_t SumPt;
593
594  // N-subjettiness variables
595
596  Float_t Tau[5];
597
598  // Other Substructure variables
599
600  TLorentzVector TrimmedP4[5]; // first entry (i = 0) is the total Trimmed Jet 4-momenta and from i = 1 to 4 are the trimmed subjets 4-momenta
601  TLorentzVector PrunedP4[5]; // first entry (i = 0) is the total Pruned Jet 4-momenta and from i = 1 to 4 are the pruned subjets 4-momenta
602  TLorentzVector SoftDroppedP4[5]; // first entry (i = 0) is the total SoftDropped Jet 4-momenta and from i = 1 to 4 are the pruned subjets 4-momenta
603
604  Int_t NSubJetsTrimmed; // number of subjets trimmed
605  Int_t NSubJetsPruned; // number of subjets pruned
606  Int_t NSubJetsSoftDropped; // number of subjets soft-dropped
607
608
609  static CompBase *fgCompare; //!
610  const CompBase *GetCompare() const { return fgCompare; }
611
612  void AddCandidate(Candidate *object);
613  TObjArray *GetCandidates();
614
615  Bool_t Overlaps(const Candidate *object) const;
616
617  virtual void Copy(TObject &object) const;
618  virtual TObject *Clone(const char *newname = "") const;
619  virtual void Clear(Option_t* option = "");
620
621private:
622  DelphesFactory *fFactory; //!
623  TObjArray *fArray; //!
624
625  void SetFactory(DelphesFactory *factory) { fFactory = factory; }
626
627  ClassDef(Candidate, 4)
628};
629
630#endif // DelphesClasses_h
631
632
Note: See TracBrowser for help on using the repository browser.