Fork me on GitHub

source: git/classes/DelphesClasses.h@ e57c062

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

Merge branch 'master' into photonId

Conflicts:

cards/CMS_PhaseII/CMS_PhaseII_200PU_v03.tcl

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