Fork me on GitHub

source: git/classes/DelphesClasses.h @ 341014c

ImprovedOutputFileTimingllp
Last change on this file since 341014c was 341014c, checked in by Pavel Demin <pavel-demin@…>, 18 months ago

apply .clang-format to all .h, .cc and .cpp files

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