Fork me on GitHub

source: git/classes/DelphesClasses.h@ cdbca0d

ImprovedOutputFile Timing dual_readout llp
Last change on this file since cdbca0d was edeb0f0, checked in by Michele Selvaggi <michele.selvaggi@…>, 6 years ago

added cross section and error in hepmc event class and reader

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