Fork me on GitHub

source: git/classes/DelphesClasses.h @ 01f9722

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

read xsec from LHE files produced with Wizard

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