Fork me on GitHub

source: git/classes/DelphesClasses.h@ 0ac1afd

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 0ac1afd was 6d8a29a, checked in by Pavel Demin <pavel.demin@…>, 8 years ago

update comments and rearrange some variables in DelphesClasses.h

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