Fork me on GitHub

source: git/classes/DelphesClasses.h@ 769f65b

ImprovedOutputFile Timing llp
Last change on this file since 769f65b was 769f65b, checked in by Pavel Demin <pavel-demin@…>, 6 years ago

update version of LHEFEvent, HepMCEvent, Photon, Electron, Muon, Jet and Candidate

  • 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 "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, 3)
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, 3)
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, 4)
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 Float_t D0; // track transverse impact parameter
318 Float_t DZ; // track longitudinal impact parameter
319 Float_t ErrorD0; // track transverse impact parameter error
320 Float_t ErrorDZ; // track longitudinal impact parameter error
321
322 static CompBase *fgCompare; //!
323 const CompBase *GetCompare() const { return fgCompare; }
324
325 TLorentzVector P4() const;
326
327 ClassDef(Electron, 4)
328};
329
330//---------------------------------------------------------------------------
331
332class Muon: public SortableObject
333{
334public:
335
336 Float_t PT; // muon transverse momentum
337 Float_t Eta; // muon pseudorapidity
338 Float_t Phi; // muon azimuthal angle
339
340 Float_t T; // particle arrival time of flight
341
342 Int_t Charge; // muon charge
343
344 TRef Particle; // reference to generated particle
345
346 Float_t IsolationVar; // isolation variable
347 Float_t IsolationVarRhoCorr; // isolation variable
348 Float_t SumPtCharged; // isolation variable
349 Float_t SumPtNeutral; // isolation variable
350 Float_t SumPtChargedPU; // isolation variable
351 Float_t SumPt; // isolation variable
352
353 Float_t D0; // track transverse impact parameter
354 Float_t DZ; // track longitudinal impact parameter
355 Float_t ErrorD0; // track transverse impact parameter error
356 Float_t ErrorDZ; // track longitudinal impact parameter error
357
358 static CompBase *fgCompare; //!
359 const CompBase *GetCompare() const { return fgCompare; }
360
361 TLorentzVector P4() const;
362
363 ClassDef(Muon, 4)
364};
365
366//---------------------------------------------------------------------------
367
368class Jet: public SortableObject
369{
370public:
371
372 Float_t PT; // jet transverse momentum
373 Float_t Eta; // jet pseudorapidity
374 Float_t Phi; // jet azimuthal angle
375
376 Float_t T; //particle arrival time of flight
377
378 Float_t Mass; // jet invariant mass
379
380 Float_t DeltaEta; // jet radius in pseudorapidity
381 Float_t DeltaPhi; // jet radius in azimuthal angle
382
383 UInt_t Flavor; // jet flavor
384 UInt_t FlavorAlgo; // jet flavor
385 UInt_t FlavorPhys; // jet flavor
386
387 UInt_t BTag; // 0 or 1 for a jet that has been tagged as containing a heavy quark
388 UInt_t BTagAlgo; // 0 or 1 for a jet that has been tagged as containing a heavy quark
389 UInt_t BTagPhys; // 0 or 1 for a jet that has been tagged as containing a heavy quark
390
391 UInt_t TauTag; // 0 or 1 for a jet that has been tagged as a tau
392 Float_t TauWeight; // probability for jet to be identified as tau
393
394 Int_t Charge; // tau charge
395
396 Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
397
398 Int_t NCharged; // number of charged constituents
399 Int_t NNeutrals; // number of neutral constituents
400 Float_t Beta; // (sum pt of charged pile-up constituents)/(sum pt of charged constituents)
401 Float_t BetaStar; // (sum pt of charged constituents coming from hard interaction)/(sum pt of charged constituents)
402 Float_t MeanSqDeltaR; // average distance (squared) between constituent and jet weighted by pt (squared) of constituent
403 Float_t PTD; // average pt between constituent and jet weighted by pt of constituent
404 Float_t FracPt[5]; // (sum pt of constituents within a ring 0.1*i < DeltaR < 0.1*(i+1))/(sum pt of constituents)
405
406 Float_t Tau[5]; // N-subjettiness
407
408 TLorentzVector SoftDroppedJet;
409 TLorentzVector SoftDroppedSubJet1;
410 TLorentzVector SoftDroppedSubJet2;
411
412 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
413 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
414 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
415
416
417 Int_t NSubJetsTrimmed; // number of subjets trimmed
418 Int_t NSubJetsPruned; // number of subjets pruned
419 Int_t NSubJetsSoftDropped; // number of subjets soft-dropped
420
421 Double_t ExclYmerge23;
422 Double_t ExclYmerge34;
423 Double_t ExclYmerge45;
424 Double_t ExclYmerge56;
425
426 TRefArray Constituents; // references to constituents
427 TRefArray Particles; // references to generated particles
428
429 static CompBase *fgCompare; //!
430 const CompBase *GetCompare() const { return fgCompare; }
431
432 TLorentzVector P4() const;
433 TLorentzVector Area;
434
435 ClassDef(Jet, 4)
436};
437
438//---------------------------------------------------------------------------
439
440class Track: public SortableObject
441{
442public:
443 Int_t PID; // HEP ID number
444
445 Int_t Charge; // track charge
446
447 Float_t P; // track momentum
448 Float_t PT; // track transverse momentum
449 Float_t Eta; // track pseudorapidity
450 Float_t Phi; // track azimuthal angle
451 Float_t CtgTheta; // track cotangent of theta
452
453 Float_t EtaOuter; // track pseudorapidity at the tracker edge
454 Float_t PhiOuter; // track azimuthal angle at the tracker edge
455
456 Float_t T; // track vertex position (t component)
457 Float_t X; // track vertex position (x component)
458 Float_t Y; // track vertex position (y component)
459 Float_t Z; // track vertex position (z component)
460
461 Float_t TOuter; // track position (t component) at the tracker edge
462 Float_t XOuter; // track position (x component) at the tracker edge
463 Float_t YOuter; // track position (y component) at the tracker edge
464 Float_t ZOuter; // track position (z component) at the tracker edge
465
466 Float_t Xd; // X coordinate of point of closest approach to vertex
467 Float_t Yd; // Y coordinate of point of closest approach to vertex
468 Float_t Zd; // Z coordinate of point of closest approach to vertex
469
470 Float_t L; // track path length
471 Float_t D0; // track transverse impact parameter
472 Float_t DZ; // track longitudinal impact parameter
473
474 Float_t ErrorP; // track momentum error
475 Float_t ErrorPT; // track transverse momentum error
476 Float_t ErrorPhi; // track azimuthal angle error
477 Float_t ErrorCtgTheta; // track cotangent of theta error
478
479 Float_t ErrorT; // time measurement error
480 Float_t ErrorD0; // track transverse impact parameter error
481 Float_t ErrorDZ; // track longitudinal impact parameter error
482
483 TRef Particle; // reference to generated particle
484
485 Int_t VertexIndex; // reference to vertex
486
487 static CompBase *fgCompare; //!
488 const CompBase *GetCompare() const { return fgCompare; }
489
490 TLorentzVector P4() const;
491
492 ClassDef(Track, 3)
493};
494
495//---------------------------------------------------------------------------
496
497class Tower: public SortableObject
498{
499public:
500 Float_t ET; // calorimeter tower transverse energy
501 Float_t Eta; // calorimeter tower pseudorapidity
502 Float_t Phi; // calorimeter tower azimuthal angle
503
504 Float_t E; // calorimeter tower energy
505
506 Float_t T; // ecal deposit time, averaged by sqrt(EM energy) over all particles, not smeared
507 Int_t NTimeHits; // number of hits contributing to time measurement
508
509 Float_t Eem; // calorimeter tower electromagnetic energy
510 Float_t Ehad; // calorimeter tower hadronic energy
511
512 Float_t Edges[4]; // calorimeter tower edges
513
514 TRefArray Particles; // references to generated particles
515
516 static CompBase *fgCompare; //!
517 const CompBase *GetCompare() const { return fgCompare; }
518
519 TLorentzVector P4() const;
520
521 ClassDef(Tower, 2)
522};
523
524//---------------------------------------------------------------------------
525
526class HectorHit: public SortableObject
527{
528public:
529 Float_t E; // reconstructed energy [GeV]
530
531 Float_t Tx; // angle of the momentum in the horizontal (x,z) plane [urad]
532 Float_t Ty; // angle of the momentum in the verical (y,z) plane [urad]
533
534 Float_t T; // time of flight to the detector [s]
535
536 Float_t X; // horizontal distance to the beam [um]
537 Float_t Y; // vertical distance to the beam [um]
538 Float_t S; // distance to the interaction point [m]
539
540 TRef Particle; // reference to generated particle
541
542 static CompBase *fgCompare; //!
543 const CompBase *GetCompare() const { return fgCompare; }
544
545 ClassDef(HectorHit, 1)
546};
547
548//---------------------------------------------------------------------------
549
550class Candidate: public SortableObject
551{
552 friend class DelphesFactory;
553
554public:
555 Candidate();
556
557 Int_t PID;
558
559 Int_t Status;
560 Int_t M1, M2, D1, D2;
561
562 Int_t Charge;
563
564 Float_t Mass;
565
566 Int_t IsPU;
567 Int_t IsRecoPU;
568
569 Int_t IsConstituent;
570 Int_t IsFromConversion;
571
572 UInt_t Flavor;
573 UInt_t FlavorAlgo;
574 UInt_t FlavorPhys;
575
576 UInt_t BTag;
577 UInt_t BTagAlgo;
578 UInt_t BTagPhys;
579
580 UInt_t TauTag;
581 Float_t TauWeight;
582
583 Float_t Eem;
584 Float_t Ehad;
585
586 Float_t Edges[4];
587 Float_t DeltaEta;
588 Float_t DeltaPhi;
589
590 TLorentzVector Momentum, Position, InitialPosition, PositionError, Area;
591
592 Float_t L; // path length
593 Float_t ErrorT; // path length
594 Float_t D0;
595 Float_t ErrorD0;
596 Float_t DZ;
597 Float_t ErrorDZ;
598 Float_t P;
599 Float_t ErrorP;
600 Float_t PT;
601 Float_t ErrorPT;
602 Float_t CtgTheta;
603 Float_t ErrorCtgTheta;
604 Float_t Phi;
605 Float_t ErrorPhi;
606
607 Float_t Xd;
608 Float_t Yd;
609 Float_t Zd;
610
611 // tracking resolution
612
613 Float_t TrackResolution;
614
615 // PileUpJetID variables
616
617 Int_t NCharged;
618 Int_t NNeutrals;
619 Float_t Beta;
620 Float_t BetaStar;
621 Float_t MeanSqDeltaR;
622 Float_t PTD;
623 Float_t FracPt[5];
624
625 // Timing information
626
627 Int_t NTimeHits;
628 std::vector< std::pair< Float_t, Float_t > > ECalEnergyTimePairs;
629
630 // Isolation variables
631
632 Float_t IsolationVar;
633 Float_t IsolationVarRhoCorr;
634 Float_t SumPtCharged;
635 Float_t SumPtNeutral;
636 Float_t SumPtChargedPU;
637 Float_t SumPt;
638
639 // vertex variables
640
641 Int_t ClusterIndex;
642 Int_t ClusterNDF;
643 Double_t ClusterSigma;
644 Double_t SumPT2;
645 Double_t BTVSumPT2;
646 Double_t GenDeltaZ;
647 Double_t GenSumPT2;
648
649 // N-subjettiness variables
650
651 Float_t Tau[5];
652
653 // Other Substructure variables
654
655 TLorentzVector SoftDroppedJet;
656 TLorentzVector SoftDroppedSubJet1;
657 TLorentzVector SoftDroppedSubJet2;
658
659 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
660 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
661 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
662
663 Int_t NSubJetsTrimmed; // number of subjets trimmed
664 Int_t NSubJetsPruned; // number of subjets pruned
665 Int_t NSubJetsSoftDropped; // number of subjets soft-dropped
666
667 // Exclusive clustering variables
668 Double_t ExclYmerge23;
669 Double_t ExclYmerge34;
670 Double_t ExclYmerge45;
671 Double_t ExclYmerge56;
672
673 static CompBase *fgCompare; //!
674 const CompBase *GetCompare() const { return fgCompare; }
675
676 void AddCandidate(Candidate *object);
677 TObjArray *GetCandidates();
678
679 Bool_t Overlaps(const Candidate *object) const;
680
681 virtual void Copy(TObject &object) const;
682 virtual TObject *Clone(const char *newname = "") const;
683 virtual void Clear(Option_t* option = "");
684
685private:
686 DelphesFactory *fFactory; //!
687 TObjArray *fArray; //!
688
689 void SetFactory(DelphesFactory *factory) { fFactory = factory; }
690
691 ClassDef(Candidate, 6)
692};
693
694#endif // DelphesClasses_h
695
696
Note: See TracBrowser for help on using the repository browser.